xref: /petsc/src/sys/objects/init.c (revision af359df37090132bfea54ad0f69df140649f0fd8)
1e5c89e4eSSatish Balay /*
2e5c89e4eSSatish Balay 
3e5c89e4eSSatish Balay    This file defines part of the initialization of PETSc
4e5c89e4eSSatish Balay 
5e5c89e4eSSatish Balay   This file uses regular malloc and free because it cannot know
6e5c89e4eSSatish Balay   what malloc is being used until it has already processed the input.
7e5c89e4eSSatish Balay */
8e5c89e4eSSatish Balay 
951d315f7SKerry Stevens #define _GNU_SOURCE
1051d315f7SKerry Stevens #include <sched.h>
11c6db04a5SJed Brown #include <petscsys.h>        /*I  "petscsys.h"   I*/
12ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD)
1351dcc849SKerry Stevens #include <pthread.h>
14ba61063dSBarry Smith #endif
15ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
1651d315f7SKerry Stevens #include <sys/sysinfo.h>
17ba61063dSBarry Smith #endif
1851d315f7SKerry Stevens #include <unistd.h>
19e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
20e5c89e4eSSatish Balay #include <stdlib.h>
21e5c89e4eSSatish Balay #endif
22e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
23e5c89e4eSSatish Balay #include <malloc.h>
24e5c89e4eSSatish Balay #endif
25555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
26555d055bSBarry Smith #include <valgrind/valgrind.h>
27555d055bSBarry Smith #endif
28555d055bSBarry Smith 
29e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
30e5c89e4eSSatish Balay /*
31e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
32e5c89e4eSSatish Balay    already started before PETSc was initialized.
33e5c89e4eSSatish Balay */
347087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
357087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
367087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
3751dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
3851dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
397087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
407087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
41ba61063dSBarry Smith 
42ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
4351dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
4451dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
45*af359df3SBarry Smith #define PETSC_HAVE_PTHREAD_BARRIER
46ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
47ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
48ba61063dSBarry Smith #endif
4951d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
50f09cb4aaSKerry Stevens int*         pVal;
5151dcc849SKerry Stevens 
52ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
5351d315f7SKerry Stevens int* ThreadCoreAffinity;
5451d315f7SKerry Stevens 
55ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
5651d315f7SKerry Stevens 
5751d315f7SKerry Stevens typedef struct {
5851d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
5951d315f7SKerry Stevens   pthread_cond_t**  cond1array;
6051d315f7SKerry Stevens   pthread_cond_t** cond2array;
6151d315f7SKerry Stevens   void* (*pfunc)(void*);
6251d315f7SKerry Stevens   void** pdata;
6351d315f7SKerry Stevens   PetscBool startJob;
6451d315f7SKerry Stevens   estat eJobStat;
6551d315f7SKerry Stevens   PetscBool** arrThreadStarted;
6651d315f7SKerry Stevens   PetscBool** arrThreadReady;
6751d315f7SKerry Stevens } sjob_tree;
6851d315f7SKerry Stevens sjob_tree job_tree;
6951d315f7SKerry Stevens typedef struct {
7051d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
7151d315f7SKerry Stevens   pthread_cond_t**  cond1array;
7251d315f7SKerry Stevens   pthread_cond_t** cond2array;
7351d315f7SKerry Stevens   void* (*pfunc)(void*);
7451d315f7SKerry Stevens   void** pdata;
7551d315f7SKerry Stevens   PetscBool** arrThreadReady;
7651d315f7SKerry Stevens } sjob_main;
7751d315f7SKerry Stevens sjob_main job_main;
7851d315f7SKerry Stevens typedef struct {
7951d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8051d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8151d315f7SKerry Stevens   pthread_cond_t** cond2array;
8251d315f7SKerry Stevens   void* (*pfunc)(void*);
8351d315f7SKerry Stevens   void** pdata;
8451d315f7SKerry Stevens   PetscBool startJob;
8551d315f7SKerry Stevens   estat eJobStat;
8651d315f7SKerry Stevens   PetscBool** arrThreadStarted;
8751d315f7SKerry Stevens   PetscBool** arrThreadReady;
8851d315f7SKerry Stevens } sjob_chain;
8951d315f7SKerry Stevens sjob_chain job_chain;
90ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
9151dcc849SKerry Stevens typedef struct {
9251dcc849SKerry Stevens   pthread_mutex_t mutex;
9351dcc849SKerry Stevens   pthread_cond_t cond;
9451dcc849SKerry Stevens   void* (*pfunc)(void*);
9551dcc849SKerry Stevens   void** pdata;
9651dcc849SKerry Stevens   pthread_barrier_t* pbarr;
9751dcc849SKerry Stevens   int iNumJobThreads;
9851dcc849SKerry Stevens   int iNumReadyThreads;
9951dcc849SKerry Stevens   PetscBool startJob;
10051d315f7SKerry Stevens } sjob_true;
10151d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
102ba61063dSBarry Smith #endif
10351dcc849SKerry Stevens 
104ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
105ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
106ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
107ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
108ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
109ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
11051dcc849SKerry Stevens 
11151d315f7SKerry Stevens /* Function Pointers */
11251d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
11351d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
11451d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
11551d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
11651d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
11751d315f7SKerry Stevens /**** Tree Functions ****/
11851d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
11951d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
12051d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
12151d315f7SKerry Stevens void           MainWait_Tree(void);
12251d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
12351d315f7SKerry Stevens /**** Main Functions ****/
12451d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
12551d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
12651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
12751d315f7SKerry Stevens void           MainWait_Main(void);
12851d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
12951d315f7SKerry Stevens /**** Chain Functions ****/
13051d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
13151d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
13251d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
13351d315f7SKerry Stevens void           MainWait_Chain(void);
13451d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
13551d315f7SKerry Stevens /**** True Functions ****/
13651d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
13751d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
13851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
13951d315f7SKerry Stevens void           MainWait_True(void);
14051d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
14151d315f7SKerry Stevens /****  ****/
14251d315f7SKerry Stevens 
14351dcc849SKerry Stevens void* FuncFinish(void*);
1440ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1450ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
146ba61063dSBarry Smith #endif
147e5c89e4eSSatish Balay 
148e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
149e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
150e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
151e5c89e4eSSatish Balay #endif
1522c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1537087cfbeSBarry Smith MPI_Datatype   MPI_C_DOUBLE_COMPLEX;
1547087cfbeSBarry Smith MPI_Datatype   MPI_C_COMPLEX;
1552c876bd9SBarry Smith #endif
1567087cfbeSBarry Smith PetscScalar    PETSC_i;
157e5c89e4eSSatish Balay #else
1587087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
159e5c89e4eSSatish Balay #endif
160ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
161c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
162c90a1750SBarry Smith #endif
1637087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1647087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
16575567043SBarry Smith 
166e5c89e4eSSatish Balay /*
167e5c89e4eSSatish Balay      These are needed by petscbt.h
168e5c89e4eSSatish Balay */
169c6db04a5SJed Brown #include <petscbt.h>
1707087cfbeSBarry Smith char      _BT_mask = ' ';
1717087cfbeSBarry Smith char      _BT_c = ' ';
1727087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
173e5c89e4eSSatish Balay 
174e5c89e4eSSatish Balay /*
175e5c89e4eSSatish Balay        Function that is called to display all error messages
176e5c89e4eSSatish Balay */
1777087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1787087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
179238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1807087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
181238ccf28SShri Abhyankar #else
1827087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
183238ccf28SShri Abhyankar #endif
184bab1f7e6SVictor Minden /*
1858154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1868154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
187bab1f7e6SVictor Minden 
188e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
189e5c89e4eSSatish Balay /*
190e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
191e5c89e4eSSatish Balay */
192e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
193e5c89e4eSSatish Balay 
194e5c89e4eSSatish Balay #undef __FUNCT__
195f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
1967087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
197e5c89e4eSSatish Balay {
198e5c89e4eSSatish Balay   PetscErrorCode ierr;
199e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
200e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
201e5c89e4eSSatish Balay   char           version[256];
202e5c89e4eSSatish Balay 
203e5c89e4eSSatish Balay   PetscFunctionBegin;
204e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
205e5c89e4eSSatish Balay   if (!rank) {
206e5c89e4eSSatish Balay     char        arch[10];
207f56c2debSBarry Smith     int         err;
20888c29154SBarry Smith     PetscViewer viewer;
209f56c2debSBarry Smith 
210e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
211e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
212a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
213e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
214e5c89e4eSSatish Balay     if (filename) {
215e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
216e5c89e4eSSatish Balay     } else {
217e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
218e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
219e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
220e5c89e4eSSatish Balay     }
221e5c89e4eSSatish Balay 
222e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
223e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
224e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
225e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
226e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
22788c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
22888c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2296bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
230e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
231f56c2debSBarry Smith     err = fflush(*fd);
232e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
233e5c89e4eSSatish Balay   }
234e5c89e4eSSatish Balay   PetscFunctionReturn(0);
235e5c89e4eSSatish Balay }
236e5c89e4eSSatish Balay 
237e5c89e4eSSatish Balay #undef __FUNCT__
238f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2397087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
240e5c89e4eSSatish Balay {
241e5c89e4eSSatish Balay   PetscErrorCode ierr;
242e5c89e4eSSatish Balay   PetscMPIInt    rank;
243e5c89e4eSSatish Balay   char           date[64];
244f56c2debSBarry Smith   int            err;
245e5c89e4eSSatish Balay 
246e5c89e4eSSatish Balay   PetscFunctionBegin;
247e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
248e5c89e4eSSatish Balay   if (!rank) {
249e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
250e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
251e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
252e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
253f56c2debSBarry Smith     err = fflush(*fd);
254e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
255f56c2debSBarry Smith     err = fclose(*fd);
256e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
257e5c89e4eSSatish Balay   }
258e5c89e4eSSatish Balay   PetscFunctionReturn(0);
259e5c89e4eSSatish Balay }
260e5c89e4eSSatish Balay 
261e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
262e5c89e4eSSatish Balay 
263e5c89e4eSSatish Balay /*
264e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
265e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
266e5c89e4eSSatish Balay 
267e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2683c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2693c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
270e5c89e4eSSatish Balay */
271e5c89e4eSSatish Balay 
272e5c89e4eSSatish Balay #undef __FUNCT__
273e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
274e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
275e5c89e4eSSatish Balay {
276e5c89e4eSSatish Balay   PetscFunctionBegin;
2773c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
278e5c89e4eSSatish Balay   abort();
279e5c89e4eSSatish Balay }
280e5c89e4eSSatish Balay 
281e5c89e4eSSatish Balay #undef __FUNCT__
282e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
283e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
284e5c89e4eSSatish Balay {
285e5c89e4eSSatish Balay   PetscErrorCode ierr;
286e5c89e4eSSatish Balay 
287e5c89e4eSSatish Balay   PetscFunctionBegin;
2883c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
289e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
290e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
2913c311c98SBarry Smith     MPI_Abort(*comm,*flag);
292e5c89e4eSSatish Balay   }
293e5c89e4eSSatish Balay }
294e5c89e4eSSatish Balay 
295e5c89e4eSSatish Balay #undef __FUNCT__
296e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
297e5c89e4eSSatish Balay /*@C
298e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
299e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
300e5c89e4eSSatish Balay 
301e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
302e5c89e4eSSatish Balay 
303e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
304e5c89e4eSSatish Balay 
305e5c89e4eSSatish Balay    Level: advanced
306e5c89e4eSSatish Balay 
307e5c89e4eSSatish Balay    Note:
308e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
309e5c89e4eSSatish Balay 
31088c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
311e5c89e4eSSatish Balay @*/
3127087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
313e5c89e4eSSatish Balay {
314e5c89e4eSSatish Balay   PetscFunctionBegin;
315e5c89e4eSSatish Balay   PetscFinalize();
316e5c89e4eSSatish Balay   exit(0);
317e5c89e4eSSatish Balay   return 0;
318e5c89e4eSSatish Balay }
319e5c89e4eSSatish Balay 
320ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
32109573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
322ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
323e5c89e4eSSatish Balay static char       emacsmachinename[256];
324e5c89e4eSSatish Balay 
325e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
326e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
327e5c89e4eSSatish Balay 
328e5c89e4eSSatish Balay #undef __FUNCT__
329e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
330e5c89e4eSSatish Balay /*@C
331e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
332e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
333e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
334e5c89e4eSSatish Balay 
335e5c89e4eSSatish Balay    Input Parameter:
336e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
337da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
338e5c89e4eSSatish Balay 
339e5c89e4eSSatish Balay    Level: developer
340e5c89e4eSSatish Balay 
341e5c89e4eSSatish Balay    Concepts: package help message
342e5c89e4eSSatish Balay 
343e5c89e4eSSatish Balay @*/
3447087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
345e5c89e4eSSatish Balay {
346e5c89e4eSSatish Balay   PetscFunctionBegin;
347e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
348e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
349e5c89e4eSSatish Balay   PetscFunctionReturn(0);
350e5c89e4eSSatish Balay }
351e5c89e4eSSatish Balay 
352e5c89e4eSSatish Balay #undef __FUNCT__
353e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3547087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
355e5c89e4eSSatish Balay {
356e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
357e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
358ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
359e5c89e4eSSatish Balay   PetscErrorCode ierr;
360a6d0e24fSJed Brown   PetscReal      si;
361e5c89e4eSSatish Balay   int            i;
362e5c89e4eSSatish Balay   PetscMPIInt    rank;
363e5c89e4eSSatish Balay   char           version[256];
364e5c89e4eSSatish Balay 
365e5c89e4eSSatish Balay   PetscFunctionBegin;
366e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
367e5c89e4eSSatish Balay 
368e5c89e4eSSatish Balay   /*
369e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
370e5c89e4eSSatish Balay   */
3718bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
37281b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
373acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
374e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
375555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
376555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
377555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
378555d055bSBarry Smith #endif
379e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
380555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
381555d055bSBarry Smith     }
382555d055bSBarry Smith #endif
383e5c89e4eSSatish Balay   }
384e5c89e4eSSatish Balay #else
385acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
386acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
387e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
388e5c89e4eSSatish Balay #endif
389e5c89e4eSSatish Balay   if (flg3) {
390e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
391e5c89e4eSSatish Balay   }
39290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
393acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
394e5c89e4eSSatish Balay   if (flg1) {
395e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
396e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
397e5c89e4eSSatish Balay   }
398e5c89e4eSSatish Balay 
39990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
400acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4017783f70dSSatish Balay   if (!flg1) {
40290d69ab7SBarry Smith     flg1 = PETSC_FALSE;
403acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4047783f70dSSatish Balay   }
405e5c89e4eSSatish Balay   if (flg1) {
406e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
407e5c89e4eSSatish Balay   }
408e5c89e4eSSatish Balay 
409e5c89e4eSSatish Balay   /*
410e5c89e4eSSatish Balay       Set the display variable for graphics
411e5c89e4eSSatish Balay   */
412e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
413e5c89e4eSSatish Balay 
414fc087633SBarry Smith 
41551dcc849SKerry Stevens   /*
416e5c89e4eSSatish Balay       Print the PETSc version information
417e5c89e4eSSatish Balay   */
418e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
419e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
420e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
421e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
422e5c89e4eSSatish Balay 
423e5c89e4eSSatish Balay     /*
424e5c89e4eSSatish Balay        Print "higher-level" package version message
425e5c89e4eSSatish Balay     */
426e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
427e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
428e5c89e4eSSatish Balay     }
429e5c89e4eSSatish Balay 
430a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
431e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
432e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
433e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
434e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
435e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
43684e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
437e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
438e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
439e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
440e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
441e5c89e4eSSatish Balay   }
442e5c89e4eSSatish Balay 
443e5c89e4eSSatish Balay   /*
444e5c89e4eSSatish Balay        Print "higher-level" package help message
445e5c89e4eSSatish Balay   */
446e5c89e4eSSatish Balay   if (flg3){
447e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
448e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
449e5c89e4eSSatish Balay     }
450e5c89e4eSSatish Balay   }
451e5c89e4eSSatish Balay 
452e5c89e4eSSatish Balay   /*
453e5c89e4eSSatish Balay       Setup the error handling
454e5c89e4eSSatish Balay   */
45590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
456acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
457cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
45890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
459acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
460cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
46190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
462acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
463e5c89e4eSSatish Balay   if (flg1) {
464e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
465e5c89e4eSSatish Balay   }
46690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
467acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
468cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
46996cc47afSJed Brown   flg1 = PETSC_FALSE;
470acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
47196cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
472e5c89e4eSSatish Balay 
473e5c89e4eSSatish Balay   /*
474e5c89e4eSSatish Balay       Setup debugger information
475e5c89e4eSSatish Balay   */
476e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
477e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
478e5c89e4eSSatish Balay   if (flg1) {
479e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
480e5c89e4eSSatish Balay 
481e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
482e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
483e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
484e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
485e5c89e4eSSatish Balay   }
4865e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
4875e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
488e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
489e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
490e5c89e4eSSatish Balay   if (flg1 || flg2) {
491e5c89e4eSSatish Balay     PetscMPIInt    size;
492e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
493e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
494e5c89e4eSSatish Balay     /*
495e5c89e4eSSatish Balay        we have to make sure that all processors have opened
496e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
497e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
498e5c89e4eSSatish Balay        and kill the program.
499e5c89e4eSSatish Balay     */
500e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
501e5c89e4eSSatish Balay     if (size > 2) {
502533163c2SBarry Smith       PetscMPIInt dummy = 0;
503e5c89e4eSSatish Balay       MPI_Status  status;
504e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
505e5c89e4eSSatish Balay         if (rank != i) {
506e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
507e5c89e4eSSatish Balay         }
508e5c89e4eSSatish Balay       }
509e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
510e5c89e4eSSatish Balay         if (rank != i) {
511e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
512e5c89e4eSSatish Balay         }
513e5c89e4eSSatish Balay       }
514e5c89e4eSSatish Balay     }
515e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
516e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
517e5c89e4eSSatish Balay     lsize = size;
518e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
519e5c89e4eSSatish Balay     if (flag) {
520e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
521e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
522e5c89e4eSSatish Balay       }
523e5c89e4eSSatish Balay     }
524e5c89e4eSSatish Balay     if (!flag) {
525e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
526e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
527e5c89e4eSSatish Balay       if (flg1) {
528e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
529e5c89e4eSSatish Balay       } else {
530e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
531e5c89e4eSSatish Balay       }
532e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
533e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
534e5c89e4eSSatish Balay     }
535e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
536e5c89e4eSSatish Balay   }
537e5c89e4eSSatish Balay 
538e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
539cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
540e5c89e4eSSatish Balay 
54193ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
54222b84c2fSbcordonn   /*
54322b84c2fSbcordonn     Activates new sockets for zope if needed
54422b84c2fSbcordonn   */
54584ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
546d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5476dc8fec2Sbcordonn   if (flgz){
54822b84c2fSbcordonn     int  sockfd;
549f1384234SBarry Smith     char hostname[256];
55022b84c2fSbcordonn     char username[256];
5516dc8fec2Sbcordonn     int  remoteport = 9999;
5529c4c166aSBarry Smith 
55384ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
55484ab5442Sbcordonn     if (!hostname[0]){
5559c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5569c4c166aSBarry Smith     }
55722b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5589c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
55922b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
56022b84c2fSbcordonn     if (flgzout){
56122b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
562606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5636dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5649c4c166aSBarry Smith     } else {
565d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
566d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5679c4c166aSBarry Smith     }
5689c4c166aSBarry Smith   }
56993ba235fSBarry Smith #endif
570ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
571ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
572ffc871a5SBarry Smith   if (flgz){
573ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
574ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
575ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
576ffc871a5SBarry Smith   }
577ffc871a5SBarry Smith #endif
5786dc8fec2Sbcordonn 
579e5c89e4eSSatish Balay   /*
580e5c89e4eSSatish Balay         Setup profiling and logging
581e5c89e4eSSatish Balay   */
5826cf91177SBarry Smith #if defined (PETSC_USE_INFO)
5838bb29257SSatish Balay   {
584e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
5856cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
5868bb29257SSatish Balay     if (flg1 && logname[0]) {
587fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
5888bb29257SSatish Balay     } else if (flg1) {
589fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
590e5c89e4eSSatish Balay     }
591e5c89e4eSSatish Balay   }
592865f6aa8SSatish Balay #endif
593865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
594865f6aa8SSatish Balay   mname[0] = 0;
595f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
596865f6aa8SSatish Balay   if (flg1) {
597865f6aa8SSatish Balay     if (mname[0]) {
598f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
599865f6aa8SSatish Balay     } else {
600f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
601865f6aa8SSatish Balay     }
602865f6aa8SSatish Balay   }
603e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
60490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
605fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
606e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
607e5c89e4eSSatish Balay #endif
60890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
60990d69ab7SBarry Smith   flg2 = PETSC_FALSE;
61090d69ab7SBarry Smith   flg3 = PETSC_FALSE;
611acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
612acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
613d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6149f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
615e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6169f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
617e5c89e4eSSatish Balay 
618e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
619e5c89e4eSSatish Balay   if (flg1) {
620e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
621e5c89e4eSSatish Balay     FILE *file;
622e5c89e4eSSatish Balay     if (mname[0]) {
623e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
624e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
625e5c89e4eSSatish Balay       file = fopen(fname,"w");
626f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
627e5c89e4eSSatish Balay     } else {
628da9f1d6bSBarry Smith       file = PETSC_STDOUT;
629e5c89e4eSSatish Balay     }
630e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
631e5c89e4eSSatish Balay   }
632e5c89e4eSSatish Balay #endif
633e5c89e4eSSatish Balay 
634e5c89e4eSSatish Balay   /*
635e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
636e5c89e4eSSatish Balay   */
63763d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
638e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
639e5c89e4eSSatish Balay #endif
640e5c89e4eSSatish Balay 
641acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
642e5c89e4eSSatish Balay 
643*af359df3SBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
644*af359df3SBarry Smith   /*
645*af359df3SBarry Smith       Determine whether user specified maximum number of threads
646*af359df3SBarry Smith    */
647*af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
648*af359df3SBarry Smith 
649*af359df3SBarry Smith   /*
650*af359df3SBarry Smith       Determine whether to use thread pool
651*af359df3SBarry Smith    */
652*af359df3SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
653*af359df3SBarry Smith   if (flg1) {
654*af359df3SBarry Smith     PetscUseThreadPool = PETSC_TRUE;
655*af359df3SBarry Smith     PetscInt N_CORES = get_nprocs();
656*af359df3SBarry Smith     ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
657*af359df3SBarry Smith     char tstr[9];
658*af359df3SBarry Smith     char tbuf[2];
659*af359df3SBarry Smith     strcpy(tstr,"-thread");
660*af359df3SBarry Smith     for(i=0;i<PetscMaxThreads;i++) {
661*af359df3SBarry Smith       ThreadCoreAffinity[i] = i;
662*af359df3SBarry Smith       sprintf(tbuf,"%d",i);
663*af359df3SBarry Smith       strcat(tstr,tbuf);
664*af359df3SBarry Smith       ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
665*af359df3SBarry Smith       if(flg1) {
666*af359df3SBarry Smith         ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
667*af359df3SBarry Smith         ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
668*af359df3SBarry Smith       }
669*af359df3SBarry Smith       tstr[7] = '\0';
670*af359df3SBarry Smith     }
671*af359df3SBarry Smith     /* get the thread pool type */
672*af359df3SBarry Smith     PetscInt ipool = 0;
673*af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
674*af359df3SBarry Smith 
675*af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
676*af359df3SBarry Smith     switch(ipool) {
677*af359df3SBarry Smith     case 1:
678*af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
679*af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
680*af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
681*af359df3SBarry Smith       MainWait              = &MainWait_Tree;
682*af359df3SBarry Smith       MainJob               = &MainJob_Tree;
683*af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
684*af359df3SBarry Smith       break;
685*af359df3SBarry Smith     case 2:
686*af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
687*af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
688*af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
689*af359df3SBarry Smith       MainWait              = &MainWait_Main;
690*af359df3SBarry Smith       MainJob               = &MainJob_Main;
691*af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
692*af359df3SBarry Smith       break;
693*af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
694*af359df3SBarry Smith     case 3:
695*af359df3SBarry Smith #else
696*af359df3SBarry Smith     default:
697*af359df3SBarry Smith #endif
698*af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
699*af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
700*af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
701*af359df3SBarry Smith       MainWait              = &MainWait_Chain;
702*af359df3SBarry Smith       MainJob               = &MainJob_Chain;
703*af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
704*af359df3SBarry Smith       break;
705*af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
706*af359df3SBarry Smith     default:
707*af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
708*af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
709*af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
710*af359df3SBarry Smith       MainWait              = &MainWait_True;
711*af359df3SBarry Smith       MainJob               = &MainJob_True;
712*af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
713*af359df3SBarry Smith       break;
714*af359df3SBarry Smith #endif
715*af359df3SBarry Smith     }
716*af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
717*af359df3SBarry Smith   }
718*af359df3SBarry Smith #endif
719e5c89e4eSSatish Balay   /*
720e5c89e4eSSatish Balay        Print basic help message
721e5c89e4eSSatish Balay   */
722e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
723e5c89e4eSSatish Balay   if (flg1) {
724e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
725301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
726301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
727301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
728e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
729e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
730e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
731e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
732e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
733e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
734e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
735e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
736e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
737e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
738e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
739e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
740e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
741e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
742e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
743e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
744e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
745e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
746e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7474161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7484161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
749e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
750e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
751e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
752e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
753e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
754e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
755a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
756e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
757e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
758e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
759e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary]: logging objects and events\n");CHKERRQ(ierr);
760e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
761e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
762e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
763e5c89e4eSSatish Balay #endif
7646cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
765e5c89e4eSSatish Balay #endif
766e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
767e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
768e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
769e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
770e5c89e4eSSatish Balay   }
771e5c89e4eSSatish Balay 
772a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
773e5c89e4eSSatish Balay   if (flg1) {
774e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
775e5c89e4eSSatish Balay   }
776e5c89e4eSSatish Balay 
7776cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
778e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
779e5c89e4eSSatish Balay   if (f) {
7806cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
781e5c89e4eSSatish Balay   }
782827f890bSBarry Smith 
7838154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
784c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
78573113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
78673113deaSBarry Smith   else flg1 = PETSC_FALSE;
7878154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
7888154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
789bab1f7e6SVictor Minden #endif
790192daf7cSBarry Smith 
791e5c89e4eSSatish Balay   PetscFunctionReturn(0);
792e5c89e4eSSatish Balay }
793df413903SBarry Smith 
794ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
795ba61063dSBarry Smith 
79651d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
79751d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
79851d315f7SKerry Stevens   PetscErrorCode iterr;
79951d315f7SKerry Stevens   int icorr,ierr;
80051d315f7SKerry Stevens   int* pId = (int*)arg;
80151d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
80251d315f7SKerry Stevens   PetscBool PeeOn;
80351d315f7SKerry Stevens   cpu_set_t mset;
80451d315f7SKerry Stevens 
80551d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
80651d315f7SKerry Stevens   CPU_ZERO(&mset);
80751d315f7SKerry Stevens   CPU_SET(icorr,&mset);
80851d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
80951d315f7SKerry Stevens 
81051d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
81151d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
81251d315f7SKerry Stevens   }
81351d315f7SKerry Stevens   else {
81451d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
81551d315f7SKerry Stevens   }
81651d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
817ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
81851d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
81951d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
82051d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
82151d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
82251d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
823ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
824ba61063dSBarry Smith            upon return, has the lock */
82551d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
82651d315f7SKerry Stevens         }
82751d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
82851d315f7SKerry Stevens       }
82951d315f7SKerry Stevens     }
830ba61063dSBarry Smith     /* your subordinates are now ready */
83151d315f7SKerry Stevens   }
83251d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
833ba61063dSBarry Smith   /* update your ready status */
83451d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
83551d315f7SKerry Stevens   if(ThreadId==0) {
83651d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
837ba61063dSBarry Smith     /* ignal main */
83851d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
83951d315f7SKerry Stevens   }
84051d315f7SKerry Stevens   else {
841ba61063dSBarry Smith     /* tell your boss that you're ready to work */
84251d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
84351d315f7SKerry Stevens   }
844ba61063dSBarry Smith   /* the while loop needs to have an exit
845ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
846ba61063dSBarry Smith    and calling FuncFinish */
84751d315f7SKerry Stevens   while(PetscThreadGo) {
848ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
849ba61063dSBarry Smith       waiting when you don't have to causes problems
850ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
85151d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
852ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
853ba61063dSBarry Smith        upon return, has the lock */
85451d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
85551d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
85651d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
85751d315f7SKerry Stevens     }
85851d315f7SKerry Stevens     if(ThreadId==0) {
85951d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
86051d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
86151d315f7SKerry Stevens     }
86251d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
86351d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
864ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
86551d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
86651d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
86751d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
86851d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
86951d315f7SKerry Stevens         }
87051d315f7SKerry Stevens       }
87151d315f7SKerry Stevens     }
872ba61063dSBarry Smith     /* do your job */
87351d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
87451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
87551d315f7SKerry Stevens     }
87651d315f7SKerry Stevens     else {
87751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
87851d315f7SKerry Stevens     }
87951d315f7SKerry Stevens     if(iterr!=0) {
88051d315f7SKerry Stevens       ithreaderr = 1;
88151d315f7SKerry Stevens     }
88251d315f7SKerry Stevens     if(PetscThreadGo) {
883ba61063dSBarry Smith       /* reset job, get ready for more */
88451d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
885ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
886ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
88751d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
88851d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
88951d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
89051d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
89151d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
892ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
893ba61063dSBarry Smith                upon return, has the lock */
89451d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
89551d315f7SKerry Stevens             }
89651d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
89751d315f7SKerry Stevens           }
89851d315f7SKerry Stevens 	}
899ba61063dSBarry Smith         /* your subordinates are now ready */
90051d315f7SKerry Stevens       }
90151d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
90251d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
90351d315f7SKerry Stevens       if(ThreadId==0) {
904ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
905ba61063dSBarry Smith         /* root thread signals 'main' */
90651d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
90751d315f7SKerry Stevens       }
90851d315f7SKerry Stevens       else {
909ba61063dSBarry Smith         /* signal your boss before you go to sleep */
91051d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
91151d315f7SKerry Stevens       }
91251d315f7SKerry Stevens     }
91351d315f7SKerry Stevens   }
91451d315f7SKerry Stevens   return NULL;
91551d315f7SKerry Stevens }
91651d315f7SKerry Stevens 
91751d315f7SKerry Stevens #undef __FUNCT__
91851d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
91951d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
92051d315f7SKerry Stevens   PetscInt i,ierr;
92151d315f7SKerry Stevens   int status;
92251d315f7SKerry Stevens 
92351d315f7SKerry Stevens   if(PetscUseThreadPool) {
92451d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
92551d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
92651d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
92751d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
92851d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
92951d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
93051d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
93151d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
93251d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
93351d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
93451d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
93551d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
936ba61063dSBarry Smith     /* initialize job structure */
93751d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
93851d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
93951d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
94051d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
94151d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
94251d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
94351d315f7SKerry Stevens     }
94451d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
94551d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
94651d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
94751d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
94851d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
94951d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
95051d315f7SKerry Stevens     }
95151d315f7SKerry Stevens     job_tree.pfunc = NULL;
95251d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
95351d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
95451d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
95551d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
956ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
95751d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
958ba61063dSBarry Smith     /* create threads */
95951d315f7SKerry Stevens     for(i=0; i<N; i++) {
96051d315f7SKerry Stevens       pVal[i] = i;
96151d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
962ba61063dSBarry Smith       /* should check status */
96351d315f7SKerry Stevens     }
96451d315f7SKerry Stevens   }
96551d315f7SKerry Stevens   return NULL;
96651d315f7SKerry Stevens }
96751d315f7SKerry Stevens 
96851d315f7SKerry Stevens #undef __FUNCT__
96951d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
97051d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
97151d315f7SKerry Stevens   int i,ierr;
97251d315f7SKerry Stevens   void* jstatus;
97351d315f7SKerry Stevens 
97451d315f7SKerry Stevens   PetscFunctionBegin;
97551d315f7SKerry Stevens 
97651d315f7SKerry Stevens   if(PetscUseThreadPool) {
977ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
978ba61063dSBarry Smith     /* join the threads */
97951d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
98051d315f7SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
981ba61063dSBarry Smith       /* do error checking*/
98251d315f7SKerry Stevens     }
98351d315f7SKerry Stevens     free(PetscThreadPoint);
98451d315f7SKerry Stevens     free(arrmutex);
98551d315f7SKerry Stevens     free(arrcond1);
98651d315f7SKerry Stevens     free(arrcond2);
98751d315f7SKerry Stevens     free(arrstart);
98851d315f7SKerry Stevens     free(arrready);
98951d315f7SKerry Stevens     free(job_tree.pdata);
99051d315f7SKerry Stevens     free(pVal);
99151d315f7SKerry Stevens   }
99251d315f7SKerry Stevens   else {
99351d315f7SKerry Stevens   }
99451d315f7SKerry Stevens   PetscFunctionReturn(0);
99551d315f7SKerry Stevens }
99651d315f7SKerry Stevens 
99751d315f7SKerry Stevens #undef __FUNCT__
99851d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
99951d315f7SKerry Stevens void MainWait_Tree() {
100051d315f7SKerry Stevens   int ierr;
100151d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
100251d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
100351d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
100451d315f7SKerry Stevens   }
100551d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
100651d315f7SKerry Stevens }
100751d315f7SKerry Stevens 
100851d315f7SKerry Stevens #undef __FUNCT__
100951d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
101051d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
101151d315f7SKerry Stevens   int i,ierr;
101251d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
101351d315f7SKerry Stevens   if(PetscUseThreadPool) {
101451d315f7SKerry Stevens     MainWait();
101551d315f7SKerry Stevens     job_tree.pfunc = pFunc;
101651d315f7SKerry Stevens     job_tree.pdata = data;
101751d315f7SKerry Stevens     job_tree.startJob = PETSC_TRUE;
101851d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
101951d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
102051d315f7SKerry Stevens     }
102151d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
102251d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond2array[0]);
102351d315f7SKerry Stevens     if(pFunc!=FuncFinish) {
1024ba61063dSBarry Smith       MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
102551d315f7SKerry Stevens     }
102651d315f7SKerry Stevens   }
102751d315f7SKerry Stevens   else {
102851d315f7SKerry Stevens     pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
102951d315f7SKerry Stevens     PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1030ba61063dSBarry Smith     PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
103151d315f7SKerry Stevens     free(apThread);
103251d315f7SKerry Stevens   }
103351d315f7SKerry Stevens   if(ithreaderr) {
103451d315f7SKerry Stevens     ijoberr = ithreaderr;
103551d315f7SKerry Stevens   }
103651d315f7SKerry Stevens   return ijoberr;
103751d315f7SKerry Stevens }
103851d315f7SKerry Stevens /****  ****/
103951d315f7SKerry Stevens 
104051d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
104151d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
104251d315f7SKerry Stevens   PetscErrorCode iterr;
104351d315f7SKerry Stevens   int icorr,ierr;
104451d315f7SKerry Stevens   int* pId = (int*)arg;
104551d315f7SKerry Stevens   int ThreadId = *pId;
104651d315f7SKerry Stevens   cpu_set_t mset;
104751d315f7SKerry Stevens 
104851d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
104951d315f7SKerry Stevens   CPU_ZERO(&mset);
105051d315f7SKerry Stevens   CPU_SET(icorr,&mset);
105151d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
105251d315f7SKerry Stevens 
105351d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1054ba61063dSBarry Smith   /* update your ready status */
105551d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1056ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
105751d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
105851d315f7SKerry Stevens 
1059ba61063dSBarry Smith   /* the while loop needs to have an exit
1060ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1061ba61063dSBarry Smith      and calling FuncFinish */
106251d315f7SKerry Stevens   while(PetscThreadGo) {
1063ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1064ba61063dSBarry Smith        waiting when you don't have to causes problems
1065ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
106651d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1067ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1068ba61063dSBarry Smith        upon return, has the lock */
106951d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1070ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
107151d315f7SKerry Stevens     }
107251d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
107351d315f7SKerry Stevens     if(job_main.pdata==NULL) {
107451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
107551d315f7SKerry Stevens     }
107651d315f7SKerry Stevens     else {
107751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
107851d315f7SKerry Stevens     }
107951d315f7SKerry Stevens     if(iterr!=0) {
108051d315f7SKerry Stevens       ithreaderr = 1;
108151d315f7SKerry Stevens     }
108251d315f7SKerry Stevens     if(PetscThreadGo) {
1083ba61063dSBarry Smith       /* reset job, get ready for more */
108451d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
108551d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1086ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
108751d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
108851d315f7SKerry Stevens     }
108951d315f7SKerry Stevens   }
109051d315f7SKerry Stevens   return NULL;
109151d315f7SKerry Stevens }
109251d315f7SKerry Stevens 
109351d315f7SKerry Stevens #undef __FUNCT__
109451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
109551d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
109651d315f7SKerry Stevens   PetscInt i,ierr;
109751d315f7SKerry Stevens   int status;
109851d315f7SKerry Stevens 
109951d315f7SKerry Stevens   if(PetscUseThreadPool) {
110051d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
110151d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
110251d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
110351d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
110451d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
110551d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
110651d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
110751d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
110851d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
110951d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
111051d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1111ba61063dSBarry Smith     /* initialize job structure */
111251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
111351d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
111451d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
111551d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
111651d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
111751d315f7SKerry Stevens     }
111851d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
111951d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
112051d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
112151d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
112251d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
112351d315f7SKerry Stevens     }
112451d315f7SKerry Stevens     job_main.pfunc = NULL;
112551d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
112651d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1127ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
112851d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1129ba61063dSBarry Smith     /* create threads */
113051d315f7SKerry Stevens     for(i=0; i<N; i++) {
113151d315f7SKerry Stevens       pVal[i] = i;
113251d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1133ba61063dSBarry Smith       /* error check */
113451d315f7SKerry Stevens     }
113551d315f7SKerry Stevens   }
113651d315f7SKerry Stevens   else {
113751d315f7SKerry Stevens   }
113851d315f7SKerry Stevens   return NULL;
113951d315f7SKerry Stevens }
114051d315f7SKerry Stevens 
114151d315f7SKerry Stevens #undef __FUNCT__
114251d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
114351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
114451d315f7SKerry Stevens   int i,ierr;
114551d315f7SKerry Stevens   void* jstatus;
114651d315f7SKerry Stevens 
114751d315f7SKerry Stevens   PetscFunctionBegin;
114851d315f7SKerry Stevens 
114951d315f7SKerry Stevens   if(PetscUseThreadPool) {
1150ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1151ba61063dSBarry Smith     /* join the threads */
115251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
1153ba61063dSBarry Smith       ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
115451d315f7SKerry Stevens     }
115551d315f7SKerry Stevens     free(PetscThreadPoint);
115651d315f7SKerry Stevens     free(arrmutex);
115751d315f7SKerry Stevens     free(arrcond1);
115851d315f7SKerry Stevens     free(arrcond2);
115951d315f7SKerry Stevens     free(arrstart);
116051d315f7SKerry Stevens     free(arrready);
116151d315f7SKerry Stevens     free(job_main.pdata);
116251d315f7SKerry Stevens     free(pVal);
116351d315f7SKerry Stevens   }
116451d315f7SKerry Stevens   PetscFunctionReturn(0);
116551d315f7SKerry Stevens }
116651d315f7SKerry Stevens 
116751d315f7SKerry Stevens #undef __FUNCT__
116851d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
116951d315f7SKerry Stevens void MainWait_Main() {
117051d315f7SKerry Stevens   int i,ierr;
117151d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
117251d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
117351d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
117451d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
117551d315f7SKerry Stevens     }
117651d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
117751d315f7SKerry Stevens   }
117851d315f7SKerry Stevens }
117951d315f7SKerry Stevens 
118051d315f7SKerry Stevens #undef __FUNCT__
118151d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
118251d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
118351d315f7SKerry Stevens   int i,ierr;
118451d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
118551d315f7SKerry Stevens   if(PetscUseThreadPool) {
1186ba61063dSBarry Smith     MainWait(); /* you know everyone is waiting to be signalled! */
118751d315f7SKerry Stevens     job_main.pfunc = pFunc;
118851d315f7SKerry Stevens     job_main.pdata = data;
118951d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
1190ba61063dSBarry Smith       *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
119151d315f7SKerry Stevens     }
1192ba61063dSBarry Smith     /* tell the threads to go to work */
119351d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
119451d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond2array[i]);
119551d315f7SKerry Stevens     }
119651d315f7SKerry Stevens     if(pFunc!=FuncFinish) {
1197ba61063dSBarry Smith       MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
119851d315f7SKerry Stevens     }
119951d315f7SKerry Stevens   }
120051d315f7SKerry Stevens   else {
120151d315f7SKerry Stevens     pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
120251d315f7SKerry Stevens     PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1203ba61063dSBarry Smith     PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
120451d315f7SKerry Stevens     free(apThread);
120551d315f7SKerry Stevens   }
120651d315f7SKerry Stevens   if(ithreaderr) {
120751d315f7SKerry Stevens     ijoberr = ithreaderr;
120851d315f7SKerry Stevens   }
120951d315f7SKerry Stevens   return ijoberr;
121051d315f7SKerry Stevens }
121151d315f7SKerry Stevens /****  ****/
121251d315f7SKerry Stevens 
121351d315f7SKerry Stevens /**** Chain Thread Functions ****/
121451d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
121551d315f7SKerry Stevens   PetscErrorCode iterr;
121651d315f7SKerry Stevens   int icorr,ierr;
121751d315f7SKerry Stevens   int* pId = (int*)arg;
121851d315f7SKerry Stevens   int ThreadId = *pId;
121951d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
122051d315f7SKerry Stevens   PetscBool PeeOn;
122151d315f7SKerry Stevens   cpu_set_t mset;
122251d315f7SKerry Stevens 
122351d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
122451d315f7SKerry Stevens   CPU_ZERO(&mset);
122551d315f7SKerry Stevens   CPU_SET(icorr,&mset);
122651d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
122751d315f7SKerry Stevens 
122851d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
122951d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
123051d315f7SKerry Stevens   }
123151d315f7SKerry Stevens   else {
123251d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
123351d315f7SKerry Stevens   }
123451d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1235ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
123651d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
123751d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1238ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1239ba61063dSBarry Smith        upon return, has the lock */
124051d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
124151d315f7SKerry Stevens     }
124251d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1243ba61063dSBarry Smith     /* your subordinate is now ready*/
124451d315f7SKerry Stevens   }
124551d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1246ba61063dSBarry Smith   /* update your ready status */
124751d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
124851d315f7SKerry Stevens   if(ThreadId==0) {
124951d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1250ba61063dSBarry Smith     /* signal main */
125151d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
125251d315f7SKerry Stevens   }
125351d315f7SKerry Stevens   else {
1254ba61063dSBarry Smith     /* tell your boss that you're ready to work */
125551d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
125651d315f7SKerry Stevens   }
1257ba61063dSBarry Smith   /*  the while loop needs to have an exit
1258ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1259ba61063dSBarry Smith    and calling FuncFinish */
126051d315f7SKerry Stevens   while(PetscThreadGo) {
1261ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1262ba61063dSBarry Smith        waiting when you don't have to causes problems
1263ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
126451d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1265ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1266ba61063dSBarry Smith        upon return, has the lock */
126751d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
126851d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
126951d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
127051d315f7SKerry Stevens     }
127151d315f7SKerry Stevens     if(ThreadId==0) {
127251d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
127351d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
127451d315f7SKerry Stevens     }
127551d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
127651d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1277ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
127851d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
127951d315f7SKerry Stevens     }
1280ba61063dSBarry Smith     /* do your job */
128151d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
128251d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
128351d315f7SKerry Stevens     }
128451d315f7SKerry Stevens     else {
128551d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
128651d315f7SKerry Stevens     }
128751d315f7SKerry Stevens     if(iterr!=0) {
128851d315f7SKerry Stevens       ithreaderr = 1;
128951d315f7SKerry Stevens     }
129051d315f7SKerry Stevens     if(PetscThreadGo) {
1291ba61063dSBarry Smith       /* reset job, get ready for more */
129251d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1293ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1294ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
129551d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
129651d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1297ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1298ba61063dSBarry Smith            upon return, has the lock */
129951d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
130051d315f7SKerry Stevens         }
130151d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1302ba61063dSBarry Smith         /* your subordinate is now ready */
130351d315f7SKerry Stevens       }
130451d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
130551d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
130651d315f7SKerry Stevens       if(ThreadId==0) {
1307ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1308ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
130951d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
131051d315f7SKerry Stevens       }
131151d315f7SKerry Stevens       else {
1312ba61063dSBarry Smith         /* signal your boss before you go to sleep */
131351d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
131451d315f7SKerry Stevens       }
131551d315f7SKerry Stevens     }
131651d315f7SKerry Stevens   }
131751d315f7SKerry Stevens   return NULL;
131851d315f7SKerry Stevens }
131951d315f7SKerry Stevens 
132051d315f7SKerry Stevens #undef __FUNCT__
132151d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
132251d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
132351d315f7SKerry Stevens   PetscInt i,ierr;
132451d315f7SKerry Stevens   int status;
132551d315f7SKerry Stevens 
132651d315f7SKerry Stevens   if(PetscUseThreadPool) {
132751d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
132851d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
132951d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
133051d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
133151d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
133251d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
133351d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
133451d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
133551d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
133651d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
133751d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
133851d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1339ba61063dSBarry Smith     /* initialize job structure */
134051d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
134151d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
134251d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
134351d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
134451d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
134551d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
134651d315f7SKerry Stevens     }
134751d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
134851d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
134951d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
135051d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
135151d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
135251d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
135351d315f7SKerry Stevens     }
135451d315f7SKerry Stevens     job_chain.pfunc = NULL;
135551d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
135651d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
135751d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
135851d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1359ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
136051d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1361ba61063dSBarry Smith     /* create threads */
136251d315f7SKerry Stevens     for(i=0; i<N; i++) {
136351d315f7SKerry Stevens       pVal[i] = i;
136451d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1365ba61063dSBarry Smith       /* should check error */
136651d315f7SKerry Stevens     }
136751d315f7SKerry Stevens   }
136851d315f7SKerry Stevens   else {
136951d315f7SKerry Stevens   }
137051d315f7SKerry Stevens   return NULL;
137151d315f7SKerry Stevens }
137251d315f7SKerry Stevens 
137351d315f7SKerry Stevens 
137451d315f7SKerry Stevens #undef __FUNCT__
137551d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
137651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
137751d315f7SKerry Stevens   int i,ierr;
137851d315f7SKerry Stevens   void* jstatus;
137951d315f7SKerry Stevens 
138051d315f7SKerry Stevens   PetscFunctionBegin;
138151d315f7SKerry Stevens 
138251d315f7SKerry Stevens   if(PetscUseThreadPool) {
1383ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1384ba61063dSBarry Smith     /* join the threads */
138551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
138651d315f7SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1387ba61063dSBarry Smith       /* should check error */
138851d315f7SKerry Stevens     }
138951d315f7SKerry Stevens     free(PetscThreadPoint);
139051d315f7SKerry Stevens     free(arrmutex);
139151d315f7SKerry Stevens     free(arrcond1);
139251d315f7SKerry Stevens     free(arrcond2);
139351d315f7SKerry Stevens     free(arrstart);
139451d315f7SKerry Stevens     free(arrready);
139551d315f7SKerry Stevens     free(job_chain.pdata);
139651d315f7SKerry Stevens     free(pVal);
139751d315f7SKerry Stevens   }
139851d315f7SKerry Stevens   else {
139951d315f7SKerry Stevens   }
140051d315f7SKerry Stevens   PetscFunctionReturn(0);
140151d315f7SKerry Stevens }
140251d315f7SKerry Stevens 
140351d315f7SKerry Stevens #undef __FUNCT__
140451d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
140551d315f7SKerry Stevens void MainWait_Chain() {
140651d315f7SKerry Stevens   int ierr;
140751d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
140851d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
140951d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
141051d315f7SKerry Stevens   }
141151d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
141251d315f7SKerry Stevens }
141351d315f7SKerry Stevens 
141451d315f7SKerry Stevens #undef __FUNCT__
141551d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
141651d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
141751d315f7SKerry Stevens   int i,ierr;
141851d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
141951d315f7SKerry Stevens   if(PetscUseThreadPool) {
142051d315f7SKerry Stevens     MainWait();
142151d315f7SKerry Stevens     job_chain.pfunc = pFunc;
142251d315f7SKerry Stevens     job_chain.pdata = data;
142351d315f7SKerry Stevens     job_chain.startJob = PETSC_TRUE;
142451d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
142551d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
142651d315f7SKerry Stevens     }
142751d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
142851d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond2array[0]);
142951d315f7SKerry Stevens     if(pFunc!=FuncFinish) {
1430ba61063dSBarry Smith       MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
143151d315f7SKerry Stevens     }
143251d315f7SKerry Stevens   }
143351d315f7SKerry Stevens   else {
143451d315f7SKerry Stevens     pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
143551d315f7SKerry Stevens     PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1436ba61063dSBarry Smith     PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
143751d315f7SKerry Stevens     free(apThread);
143851d315f7SKerry Stevens   }
143951d315f7SKerry Stevens   if(ithreaderr) {
144051d315f7SKerry Stevens     ijoberr = ithreaderr;
144151d315f7SKerry Stevens   }
144251d315f7SKerry Stevens   return ijoberr;
144351d315f7SKerry Stevens }
144451d315f7SKerry Stevens /****  ****/
144551d315f7SKerry Stevens 
1446ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
144751d315f7SKerry Stevens /**** True Thread Functions ****/
144851d315f7SKerry Stevens void* PetscThreadFunc_True(void* arg) {
144951d315f7SKerry Stevens   int icorr,ierr,iVal;
145051dcc849SKerry Stevens   int* pId = (int*)arg;
145151dcc849SKerry Stevens   int ThreadId = *pId;
14520ca81413SKerry Stevens   PetscErrorCode iterr;
145351d315f7SKerry Stevens   cpu_set_t mset;
145451dcc849SKerry Stevens 
145551d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
145651d315f7SKerry Stevens   CPU_ZERO(&mset);
145751d315f7SKerry Stevens   CPU_SET(icorr,&mset);
145851d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
145951d315f7SKerry Stevens 
146051d315f7SKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
146151d315f7SKerry Stevens   job_true.iNumReadyThreads++;
146251d315f7SKerry Stevens   if(job_true.iNumReadyThreads==PetscMaxThreads) {
146351dcc849SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
146451dcc849SKerry Stevens   }
1465ba61063dSBarry Smith   /*the while loop needs to have an exit
1466ba61063dSBarry Smith     the 'main' thread can terminate all the threads by performing a broadcast
1467ba61063dSBarry Smith    and calling FuncFinish */
146851dcc849SKerry Stevens   while(PetscThreadGo) {
1469ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
1470ba61063dSBarry Smith       waiting when you don't have to causes problems
1471ba61063dSBarry Smith      also need to wait if another thread sneaks in and messes with the predicate */
147251d315f7SKerry Stevens     while(job_true.startJob==PETSC_FALSE&&job_true.iNumJobThreads==0) {
1473ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1474ba61063dSBarry Smith        upon return, has the lock */
147551d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
147651dcc849SKerry Stevens     }
147751d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
147851d315f7SKerry Stevens     job_true.iNumJobThreads--;
147951d315f7SKerry Stevens     job_true.iNumReadyThreads--;
148051d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
148151d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
148251d315f7SKerry Stevens     if(job_true.pdata==NULL) {
148351d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
148451dcc849SKerry Stevens     }
148551dcc849SKerry Stevens     else {
148651d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
148751dcc849SKerry Stevens     }
14880ca81413SKerry Stevens     if(iterr!=0) {
14890ca81413SKerry Stevens       ithreaderr = 1;
14900ca81413SKerry Stevens     }
1491ba61063dSBarry Smith     /* the barrier is necessary BECAUSE: look at job_true.iNumReadyThreads
1492ba61063dSBarry Smith       what happens if a thread finishes before they all start? BAD!
1493ba61063dSBarry Smith      what happens if a thread finishes before any else start? BAD! */
1494ba61063dSBarry Smith     pthread_barrier_wait(job_true.pbarr); /* ensures all threads are finished */
1495ba61063dSBarry Smith     /* reset job */
149651dcc849SKerry Stevens     if(PetscThreadGo) {
149751d315f7SKerry Stevens       pthread_mutex_lock(&job_true.mutex);
149851d315f7SKerry Stevens       job_true.iNumReadyThreads++;
149951d315f7SKerry Stevens       if(job_true.iNumReadyThreads==PetscMaxThreads) {
1500ba61063dSBarry Smith 	/* signal the 'main' thread that the job is done! (only done once) */
150151dcc849SKerry Stevens 	ierr = pthread_cond_signal(&main_cond);
150251dcc849SKerry Stevens       }
150351dcc849SKerry Stevens     }
150451dcc849SKerry Stevens   }
150551dcc849SKerry Stevens   return NULL;
150651dcc849SKerry Stevens }
150751dcc849SKerry Stevens 
1508f09cb4aaSKerry Stevens #undef __FUNCT__
150951d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_True"
151051d315f7SKerry Stevens void* PetscThreadInitialize_True(PetscInt N) {
151151dcc849SKerry Stevens   PetscInt i;
151251dcc849SKerry Stevens   int status;
15130ca81413SKerry Stevens 
15140ca81413SKerry Stevens   if(PetscUseThreadPool) {
1515f09cb4aaSKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1516ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
151751dcc849SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1518ba61063dSBarry Smith     BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
151951d315f7SKerry Stevens     job_true.pdata = (void**)malloc(N*sizeof(void*));
152051dcc849SKerry Stevens     for(i=0; i<N; i++) {
1521f09cb4aaSKerry Stevens       pVal[i] = i;
1522f09cb4aaSKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1523ba61063dSBarry Smith       /* error check to ensure proper thread creation */
152451dcc849SKerry Stevens       status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1525ba61063dSBarry Smith       /* should check error */
152651dcc849SKerry Stevens     }
15270ca81413SKerry Stevens   }
15280ca81413SKerry Stevens   else {
15290ca81413SKerry Stevens   }
153051dcc849SKerry Stevens   return NULL;
153151dcc849SKerry Stevens }
153251dcc849SKerry Stevens 
1533f09cb4aaSKerry Stevens 
1534f09cb4aaSKerry Stevens #undef __FUNCT__
153551d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
153651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
153751dcc849SKerry Stevens   int i,ierr;
153851dcc849SKerry Stevens   void* jstatus;
153951dcc849SKerry Stevens 
154051dcc849SKerry Stevens   PetscFunctionBegin;
15410ca81413SKerry Stevens 
15420ca81413SKerry Stevens   if(PetscUseThreadPool) {
1543ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1544ba61063dSBarry Smith     /* join the threads */
154551dcc849SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
154651dcc849SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1547ba61063dSBarry Smith       /* should check error */
154851dcc849SKerry Stevens     }
154951dcc849SKerry Stevens     free(BarrPoint);
155051dcc849SKerry Stevens     free(PetscThreadPoint);
15510ca81413SKerry Stevens   }
15520ca81413SKerry Stevens   else {
15530ca81413SKerry Stevens   }
155451dcc849SKerry Stevens   PetscFunctionReturn(0);
155551dcc849SKerry Stevens }
155651dcc849SKerry Stevens 
1557f09cb4aaSKerry Stevens #undef __FUNCT__
155851d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
155951d315f7SKerry Stevens void MainWait_True() {
156051dcc849SKerry Stevens   int ierr;
156151d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
156251d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
156351dcc849SKerry Stevens   }
156451d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
156551dcc849SKerry Stevens }
156651dcc849SKerry Stevens 
1567f09cb4aaSKerry Stevens #undef __FUNCT__
156851d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
156951d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
157051dcc849SKerry Stevens   int ierr;
15710ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
15720ca81413SKerry Stevens   if(PetscUseThreadPool) {
15730ca81413SKerry Stevens     MainWait();
157451d315f7SKerry Stevens     job_true.pfunc = pFunc;
157551d315f7SKerry Stevens     job_true.pdata = data;
157651d315f7SKerry Stevens     job_true.pbarr = &BarrPoint[n];
157751d315f7SKerry Stevens     job_true.iNumJobThreads = n;
157851d315f7SKerry Stevens     job_true.startJob = PETSC_TRUE;
157951d315f7SKerry Stevens     ierr = pthread_cond_broadcast(&job_true.cond);
15800ca81413SKerry Stevens     if(pFunc!=FuncFinish) {
1581ba61063dSBarry Smith       MainWait(); /* why wait after? guarantees that job gets done */
15820ca81413SKerry Stevens     }
15830ca81413SKerry Stevens   }
15840ca81413SKerry Stevens   else {
15850ca81413SKerry Stevens     pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
15860ca81413SKerry Stevens     PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1587ba61063dSBarry Smith     PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
15880ca81413SKerry Stevens     free(apThread);
15890ca81413SKerry Stevens   }
15900ca81413SKerry Stevens   if(ithreaderr) {
15910ca81413SKerry Stevens     ijoberr = ithreaderr;
15920ca81413SKerry Stevens   }
15930ca81413SKerry Stevens   return ijoberr;
159451dcc849SKerry Stevens }
159551d315f7SKerry Stevens /****  ****/
1596ba61063dSBarry Smith #endif
159751dcc849SKerry Stevens 
159851dcc849SKerry Stevens void* FuncFinish(void* arg) {
159951dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16000ca81413SKerry Stevens   return(0);
160151dcc849SKerry Stevens }
1602ba61063dSBarry Smith 
1603ba61063dSBarry Smith #endif
1604