xref: /petsc/src/sys/objects/init.c (revision bda8bf91218abc8014db6f91a57cb182fae6533f)
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 */
8ef386f4bSSatish Balay 
9ef386f4bSSatish Balay #include <petscsys.h>        /*I  "petscsys.h"   I*/
10ef386f4bSSatish Balay 
118f54b378SBarry Smith #ifndef _GNU_SOURCE
128f54b378SBarry Smith #define _GNU_SOURCE
138f54b378SBarry Smith #endif
14cd723cd1SBarry Smith #if defined(PETSC_HAVE_SCHED_H)
158f54b378SBarry Smith #include <sched.h>
16cd723cd1SBarry Smith #endif
17*bda8bf91SBarry Smith #if defined(PETSC_HAVE_PTHREAD_H)
1851dcc849SKerry Stevens #include <pthread.h>
19ba61063dSBarry Smith #endif
20ef386f4bSSatish Balay 
21ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
2251d315f7SKerry Stevens #include <sys/sysinfo.h>
23ba61063dSBarry Smith #endif
24121deb67SSatish Balay #if defined(PETSC_HAVE_UNISTD_H)
2551d315f7SKerry Stevens #include <unistd.h>
26121deb67SSatish Balay #endif
27e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
28e5c89e4eSSatish Balay #include <stdlib.h>
29e5c89e4eSSatish Balay #endif
30e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
31e5c89e4eSSatish Balay #include <malloc.h>
32e5c89e4eSSatish Balay #endif
33555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
34555d055bSBarry Smith #include <valgrind/valgrind.h>
35555d055bSBarry Smith #endif
36555d055bSBarry Smith 
37e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
38e5c89e4eSSatish Balay /*
39e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
40e5c89e4eSSatish Balay    already started before PETSc was initialized.
41e5c89e4eSSatish Balay */
427087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
437087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
447087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
4551dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
4651dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
477087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
487087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
49ba61063dSBarry Smith 
50ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
5151dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
5251dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
53af359df3SBarry Smith #define PETSC_HAVE_PTHREAD_BARRIER
54ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
55ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
56ba61063dSBarry Smith #endif
5751d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
58f09cb4aaSKerry Stevens int*         pVal;
5951dcc849SKerry Stevens 
60ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
6151d315f7SKerry Stevens int* ThreadCoreAffinity;
6251d315f7SKerry Stevens 
63ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
6451d315f7SKerry Stevens 
6551d315f7SKerry Stevens typedef struct {
6651d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
6751d315f7SKerry Stevens   pthread_cond_t**  cond1array;
6851d315f7SKerry Stevens   pthread_cond_t** cond2array;
6951d315f7SKerry Stevens   void* (*pfunc)(void*);
7051d315f7SKerry Stevens   void** pdata;
7151d315f7SKerry Stevens   PetscBool startJob;
7251d315f7SKerry Stevens   estat eJobStat;
7351d315f7SKerry Stevens   PetscBool** arrThreadStarted;
7451d315f7SKerry Stevens   PetscBool** arrThreadReady;
7551d315f7SKerry Stevens } sjob_tree;
7651d315f7SKerry Stevens sjob_tree job_tree;
7751d315f7SKerry Stevens typedef struct {
7851d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
7951d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8051d315f7SKerry Stevens   pthread_cond_t** cond2array;
8151d315f7SKerry Stevens   void* (*pfunc)(void*);
8251d315f7SKerry Stevens   void** pdata;
8351d315f7SKerry Stevens   PetscBool** arrThreadReady;
8451d315f7SKerry Stevens } sjob_main;
8551d315f7SKerry Stevens sjob_main job_main;
8651d315f7SKerry Stevens typedef struct {
8751d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8851d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8951d315f7SKerry Stevens   pthread_cond_t** cond2array;
9051d315f7SKerry Stevens   void* (*pfunc)(void*);
9151d315f7SKerry Stevens   void** pdata;
9251d315f7SKerry Stevens   PetscBool startJob;
9351d315f7SKerry Stevens   estat eJobStat;
9451d315f7SKerry Stevens   PetscBool** arrThreadStarted;
9551d315f7SKerry Stevens   PetscBool** arrThreadReady;
9651d315f7SKerry Stevens } sjob_chain;
9751d315f7SKerry Stevens sjob_chain job_chain;
98ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
9951dcc849SKerry Stevens typedef struct {
10051dcc849SKerry Stevens   pthread_mutex_t mutex;
10151dcc849SKerry Stevens   pthread_cond_t cond;
10251dcc849SKerry Stevens   void* (*pfunc)(void*);
10351dcc849SKerry Stevens   void** pdata;
10451dcc849SKerry Stevens   pthread_barrier_t* pbarr;
10551dcc849SKerry Stevens   int iNumJobThreads;
10651dcc849SKerry Stevens   int iNumReadyThreads;
10751dcc849SKerry Stevens   PetscBool startJob;
10851d315f7SKerry Stevens } sjob_true;
10951d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
110ba61063dSBarry Smith #endif
11151dcc849SKerry Stevens 
112ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
113ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
114ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
115ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
116ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
117ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
11851dcc849SKerry Stevens 
11951d315f7SKerry Stevens /* Function Pointers */
12051d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
12151d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
12251d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
12351d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
12451d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
12536d20dc5SKerry Stevens /**** Tree Thread Pool Functions ****/
12651d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
12751d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
12851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
12951d315f7SKerry Stevens void           MainWait_Tree(void);
13051d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
13136d20dc5SKerry Stevens /**** Main Thread Pool Functions ****/
13251d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
13351d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
13451d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
13551d315f7SKerry Stevens void           MainWait_Main(void);
13651d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
13736d20dc5SKerry Stevens /**** Chain Thread Pool Functions ****/
13851d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
13951d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
14051d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
14151d315f7SKerry Stevens void           MainWait_Chain(void);
14251d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
14336d20dc5SKerry Stevens /**** True Thread Pool Functions ****/
14451d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
14551d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
14651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
14751d315f7SKerry Stevens void           MainWait_True(void);
14851d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
14936d20dc5SKerry Stevens /**** NO Thread Pool Function  ****/
150683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void**,PetscInt);
15136d20dc5SKerry Stevens /****  ****/
15251dcc849SKerry Stevens void* FuncFinish(void*);
1530ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1540ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
155ba61063dSBarry Smith #endif
156e5c89e4eSSatish Balay 
157e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
158e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
159e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
160e5c89e4eSSatish Balay #endif
1612c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1627087cfbeSBarry Smith MPI_Datatype   MPI_C_DOUBLE_COMPLEX;
1637087cfbeSBarry Smith MPI_Datatype   MPI_C_COMPLEX;
1642c876bd9SBarry Smith #endif
1657087cfbeSBarry Smith PetscScalar    PETSC_i;
166e5c89e4eSSatish Balay #else
1677087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
168e5c89e4eSSatish Balay #endif
169ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
170c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
171c90a1750SBarry Smith #endif
1727087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1737087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
17475567043SBarry Smith 
175e5c89e4eSSatish Balay /*
176e5c89e4eSSatish Balay      These are needed by petscbt.h
177e5c89e4eSSatish Balay */
178c6db04a5SJed Brown #include <petscbt.h>
1797087cfbeSBarry Smith char      _BT_mask = ' ';
1807087cfbeSBarry Smith char      _BT_c = ' ';
1817087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
182e5c89e4eSSatish Balay 
183e5c89e4eSSatish Balay /*
184e5c89e4eSSatish Balay        Function that is called to display all error messages
185e5c89e4eSSatish Balay */
1867087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1877087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
188238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1897087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
190238ccf28SShri Abhyankar #else
1917087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
192238ccf28SShri Abhyankar #endif
193bab1f7e6SVictor Minden /*
1948154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1958154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
196bab1f7e6SVictor Minden 
197e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
198e5c89e4eSSatish Balay /*
199e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
200e5c89e4eSSatish Balay */
201e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
202e5c89e4eSSatish Balay 
203e5c89e4eSSatish Balay #undef __FUNCT__
204f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
2057087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
206e5c89e4eSSatish Balay {
207e5c89e4eSSatish Balay   PetscErrorCode ierr;
208e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
209e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
210e5c89e4eSSatish Balay   char           version[256];
211e5c89e4eSSatish Balay 
212e5c89e4eSSatish Balay   PetscFunctionBegin;
213e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
214e5c89e4eSSatish Balay   if (!rank) {
215e5c89e4eSSatish Balay     char        arch[10];
216f56c2debSBarry Smith     int         err;
21788c29154SBarry Smith     PetscViewer viewer;
218f56c2debSBarry Smith 
219e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
220e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
221a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
222e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
223e5c89e4eSSatish Balay     if (filename) {
224e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
225e5c89e4eSSatish Balay     } else {
226e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
227e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
228e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
229e5c89e4eSSatish Balay     }
230e5c89e4eSSatish Balay 
231e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
232e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
233e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
234e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
235e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
23688c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
23788c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2386bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
239e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
240f56c2debSBarry Smith     err = fflush(*fd);
241e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
242e5c89e4eSSatish Balay   }
243e5c89e4eSSatish Balay   PetscFunctionReturn(0);
244e5c89e4eSSatish Balay }
245e5c89e4eSSatish Balay 
246e5c89e4eSSatish Balay #undef __FUNCT__
247f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2487087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
249e5c89e4eSSatish Balay {
250e5c89e4eSSatish Balay   PetscErrorCode ierr;
251e5c89e4eSSatish Balay   PetscMPIInt    rank;
252e5c89e4eSSatish Balay   char           date[64];
253f56c2debSBarry Smith   int            err;
254e5c89e4eSSatish Balay 
255e5c89e4eSSatish Balay   PetscFunctionBegin;
256e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
257e5c89e4eSSatish Balay   if (!rank) {
258e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
259e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
260e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
261e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
262f56c2debSBarry Smith     err = fflush(*fd);
263e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
264f56c2debSBarry Smith     err = fclose(*fd);
265e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
266e5c89e4eSSatish Balay   }
267e5c89e4eSSatish Balay   PetscFunctionReturn(0);
268e5c89e4eSSatish Balay }
269e5c89e4eSSatish Balay 
270e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
271e5c89e4eSSatish Balay 
272e5c89e4eSSatish Balay /*
273e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
274e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
275e5c89e4eSSatish Balay 
276e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2773c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2783c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
279e5c89e4eSSatish Balay */
280e5c89e4eSSatish Balay 
281e5c89e4eSSatish Balay #undef __FUNCT__
282e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
283e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
284e5c89e4eSSatish Balay {
285e5c89e4eSSatish Balay   PetscFunctionBegin;
2863c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
287e5c89e4eSSatish Balay   abort();
288e5c89e4eSSatish Balay }
289e5c89e4eSSatish Balay 
290e5c89e4eSSatish Balay #undef __FUNCT__
291e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
292e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
293e5c89e4eSSatish Balay {
294e5c89e4eSSatish Balay   PetscErrorCode ierr;
295e5c89e4eSSatish Balay 
296e5c89e4eSSatish Balay   PetscFunctionBegin;
2973c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
298e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
299e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
3003c311c98SBarry Smith     MPI_Abort(*comm,*flag);
301e5c89e4eSSatish Balay   }
302e5c89e4eSSatish Balay }
303e5c89e4eSSatish Balay 
304e5c89e4eSSatish Balay #undef __FUNCT__
305e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
306e5c89e4eSSatish Balay /*@C
307e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
308e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
309e5c89e4eSSatish Balay 
310e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
311e5c89e4eSSatish Balay 
312e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
313e5c89e4eSSatish Balay 
314e5c89e4eSSatish Balay    Level: advanced
315e5c89e4eSSatish Balay 
316e5c89e4eSSatish Balay    Note:
317e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
318e5c89e4eSSatish Balay 
31988c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
320e5c89e4eSSatish Balay @*/
3217087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
322e5c89e4eSSatish Balay {
323e5c89e4eSSatish Balay   PetscFunctionBegin;
324e5c89e4eSSatish Balay   PetscFinalize();
325e5c89e4eSSatish Balay   exit(0);
326e5c89e4eSSatish Balay   return 0;
327e5c89e4eSSatish Balay }
328e5c89e4eSSatish Balay 
329ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
33009573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
331ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
332e5c89e4eSSatish Balay static char       emacsmachinename[256];
333e5c89e4eSSatish Balay 
334e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
335e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
336e5c89e4eSSatish Balay 
337e5c89e4eSSatish Balay #undef __FUNCT__
338e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
339e5c89e4eSSatish Balay /*@C
340e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
341e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
342e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
343e5c89e4eSSatish Balay 
344e5c89e4eSSatish Balay    Input Parameter:
345e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
346da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
347e5c89e4eSSatish Balay 
348e5c89e4eSSatish Balay    Level: developer
349e5c89e4eSSatish Balay 
350e5c89e4eSSatish Balay    Concepts: package help message
351e5c89e4eSSatish Balay 
352e5c89e4eSSatish Balay @*/
3537087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
354e5c89e4eSSatish Balay {
355e5c89e4eSSatish Balay   PetscFunctionBegin;
356e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
357e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
358e5c89e4eSSatish Balay   PetscFunctionReturn(0);
359e5c89e4eSSatish Balay }
360e5c89e4eSSatish Balay 
361e5c89e4eSSatish Balay #undef __FUNCT__
362e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3637087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
364e5c89e4eSSatish Balay {
365e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
366e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
367ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
368e5c89e4eSSatish Balay   PetscErrorCode ierr;
369a6d0e24fSJed Brown   PetscReal      si;
370e5c89e4eSSatish Balay   int            i;
371e5c89e4eSSatish Balay   PetscMPIInt    rank;
372e5c89e4eSSatish Balay   char           version[256];
373e5c89e4eSSatish Balay 
374e5c89e4eSSatish Balay   PetscFunctionBegin;
375e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
376e5c89e4eSSatish Balay 
377e5c89e4eSSatish Balay   /*
378e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
379e5c89e4eSSatish Balay   */
3808bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
38181b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
382acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
383e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
384555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
385555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
386555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
387555d055bSBarry Smith #endif
388e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
389555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
390555d055bSBarry Smith     }
391555d055bSBarry Smith #endif
392e5c89e4eSSatish Balay   }
393e5c89e4eSSatish Balay #else
394acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
395acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
396e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
397e5c89e4eSSatish Balay #endif
398e5c89e4eSSatish Balay   if (flg3) {
399e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
400e5c89e4eSSatish Balay   }
40190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
402acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
403e5c89e4eSSatish Balay   if (flg1) {
404e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
405e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
406e5c89e4eSSatish Balay   }
407e5c89e4eSSatish Balay 
40890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
409acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4107783f70dSSatish Balay   if (!flg1) {
41190d69ab7SBarry Smith     flg1 = PETSC_FALSE;
412acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4137783f70dSSatish Balay   }
414e5c89e4eSSatish Balay   if (flg1) {
415e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
416e5c89e4eSSatish Balay   }
417e5c89e4eSSatish Balay 
418e5c89e4eSSatish Balay   /*
419e5c89e4eSSatish Balay       Set the display variable for graphics
420e5c89e4eSSatish Balay   */
421e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
422e5c89e4eSSatish Balay 
42351dcc849SKerry Stevens 
42451dcc849SKerry Stevens   /*
425e5c89e4eSSatish Balay       Print the PETSc version information
426e5c89e4eSSatish Balay   */
427e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
428e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
429e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
430e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
431e5c89e4eSSatish Balay 
432e5c89e4eSSatish Balay     /*
433e5c89e4eSSatish Balay        Print "higher-level" package version message
434e5c89e4eSSatish Balay     */
435e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
436e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
437e5c89e4eSSatish Balay     }
438e5c89e4eSSatish Balay 
439a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
440e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
441e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
442e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
443e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
444e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
44584e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
446e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
447e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
448e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
449e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
450e5c89e4eSSatish Balay   }
451e5c89e4eSSatish Balay 
452e5c89e4eSSatish Balay   /*
453e5c89e4eSSatish Balay        Print "higher-level" package help message
454e5c89e4eSSatish Balay   */
455e5c89e4eSSatish Balay   if (flg3){
456e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
457e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
458e5c89e4eSSatish Balay     }
459e5c89e4eSSatish Balay   }
460e5c89e4eSSatish Balay 
461e5c89e4eSSatish Balay   /*
462e5c89e4eSSatish Balay       Setup the error handling
463e5c89e4eSSatish Balay   */
46490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
465acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
466cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
46790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
468acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
469cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
47090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
471acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
472e5c89e4eSSatish Balay   if (flg1) {
473e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
474e5c89e4eSSatish Balay   }
47590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
476acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
477cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
47896cc47afSJed Brown   flg1 = PETSC_FALSE;
479acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
48096cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
481e5c89e4eSSatish Balay 
482e5c89e4eSSatish Balay   /*
483e5c89e4eSSatish Balay       Setup debugger information
484e5c89e4eSSatish Balay   */
485e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
486e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
487e5c89e4eSSatish Balay   if (flg1) {
488e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
489e5c89e4eSSatish Balay 
490e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
491e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
492e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
493e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
494e5c89e4eSSatish Balay   }
4955e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
4965e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
497e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
498e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
499e5c89e4eSSatish Balay   if (flg1 || flg2) {
500e5c89e4eSSatish Balay     PetscMPIInt    size;
501e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
502e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
503e5c89e4eSSatish Balay     /*
504e5c89e4eSSatish Balay        we have to make sure that all processors have opened
505e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
506e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
507e5c89e4eSSatish Balay        and kill the program.
508e5c89e4eSSatish Balay     */
509e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
510e5c89e4eSSatish Balay     if (size > 2) {
511533163c2SBarry Smith       PetscMPIInt dummy = 0;
512e5c89e4eSSatish Balay       MPI_Status  status;
513e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
514e5c89e4eSSatish Balay         if (rank != i) {
515e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
516e5c89e4eSSatish Balay         }
517e5c89e4eSSatish Balay       }
518e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
519e5c89e4eSSatish Balay         if (rank != i) {
520e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
521e5c89e4eSSatish Balay         }
522e5c89e4eSSatish Balay       }
523e5c89e4eSSatish Balay     }
524e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
525e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
526e5c89e4eSSatish Balay     lsize = size;
527e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
528e5c89e4eSSatish Balay     if (flag) {
529e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
530e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
531e5c89e4eSSatish Balay       }
532e5c89e4eSSatish Balay     }
533e5c89e4eSSatish Balay     if (!flag) {
534e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
535e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
536e5c89e4eSSatish Balay       if (flg1) {
537e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
538e5c89e4eSSatish Balay       } else {
539e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
540e5c89e4eSSatish Balay       }
541e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
542e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
543e5c89e4eSSatish Balay     }
544e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
545e5c89e4eSSatish Balay   }
546e5c89e4eSSatish Balay 
547e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
548cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
549e5c89e4eSSatish Balay 
55093ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
55122b84c2fSbcordonn   /*
55222b84c2fSbcordonn     Activates new sockets for zope if needed
55322b84c2fSbcordonn   */
55484ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
555d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5566dc8fec2Sbcordonn   if (flgz){
55722b84c2fSbcordonn     int  sockfd;
558f1384234SBarry Smith     char hostname[256];
55922b84c2fSbcordonn     char username[256];
5606dc8fec2Sbcordonn     int  remoteport = 9999;
5619c4c166aSBarry Smith 
56284ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
56384ab5442Sbcordonn     if (!hostname[0]){
5649c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5659c4c166aSBarry Smith     }
56622b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5679c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
56822b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
56922b84c2fSbcordonn     if (flgzout){
57022b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
571606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5726dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5739c4c166aSBarry Smith     } else {
574d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
575d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5769c4c166aSBarry Smith     }
5779c4c166aSBarry Smith   }
57893ba235fSBarry Smith #endif
579ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
580ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
581ffc871a5SBarry Smith   if (flgz){
582ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
583ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
584ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
585ffc871a5SBarry Smith   }
586ffc871a5SBarry Smith #endif
5876dc8fec2Sbcordonn 
588e5c89e4eSSatish Balay   /*
589e5c89e4eSSatish Balay         Setup profiling and logging
590e5c89e4eSSatish Balay   */
5916cf91177SBarry Smith #if defined (PETSC_USE_INFO)
5928bb29257SSatish Balay   {
593e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
5946cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
5958bb29257SSatish Balay     if (flg1 && logname[0]) {
596fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
5978bb29257SSatish Balay     } else if (flg1) {
598fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
599e5c89e4eSSatish Balay     }
600e5c89e4eSSatish Balay   }
601865f6aa8SSatish Balay #endif
602865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
603865f6aa8SSatish Balay   mname[0] = 0;
604f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
605865f6aa8SSatish Balay   if (flg1) {
606865f6aa8SSatish Balay     if (mname[0]) {
607f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
608865f6aa8SSatish Balay     } else {
609f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
610865f6aa8SSatish Balay     }
611865f6aa8SSatish Balay   }
612e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
61390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
614fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
615e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
616e5c89e4eSSatish Balay #endif
61790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
61890d69ab7SBarry Smith   flg2 = PETSC_FALSE;
61990d69ab7SBarry Smith   flg3 = PETSC_FALSE;
620acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
621acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
622d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6239f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
624e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6259f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
626e5c89e4eSSatish Balay 
627e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
628e5c89e4eSSatish Balay   if (flg1) {
629e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
630e5c89e4eSSatish Balay     FILE *file;
631e5c89e4eSSatish Balay     if (mname[0]) {
632e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
633e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
634e5c89e4eSSatish Balay       file = fopen(fname,"w");
635f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
636e5c89e4eSSatish Balay     } else {
637da9f1d6bSBarry Smith       file = PETSC_STDOUT;
638e5c89e4eSSatish Balay     }
639e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
640e5c89e4eSSatish Balay   }
641e5c89e4eSSatish Balay #endif
642e5c89e4eSSatish Balay 
643e5c89e4eSSatish Balay   /*
644e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
645e5c89e4eSSatish Balay   */
64663d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
647e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
648e5c89e4eSSatish Balay #endif
649e5c89e4eSSatish Balay 
650acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
651e5c89e4eSSatish Balay 
652ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
653af359df3SBarry Smith   /*
654af359df3SBarry Smith       Determine whether user specified maximum number of threads
655af359df3SBarry Smith    */
656af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
657af359df3SBarry Smith 
658b154f58aSKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-main",&flg1);CHKERRQ(ierr);
659b154f58aSKerry Stevens   if(flg1) {
660b154f58aSKerry Stevens     cpu_set_t mset;
661b154f58aSKerry Stevens     int icorr,ncorr = get_nprocs();
662b154f58aSKerry Stevens     ierr = PetscOptionsGetInt(PETSC_NULL,"-main",&icorr,PETSC_NULL);CHKERRQ(ierr);
663b154f58aSKerry Stevens     CPU_ZERO(&mset);
664b154f58aSKerry Stevens     CPU_SET(icorr%ncorr,&mset);
665b154f58aSKerry Stevens     sched_setaffinity(0,sizeof(cpu_set_t),&mset);
666b154f58aSKerry Stevens   }
667b154f58aSKerry Stevens 
668af359df3SBarry Smith   PetscInt N_CORES = get_nprocs();
669af359df3SBarry Smith   ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
670af359df3SBarry Smith   char tstr[9];
671af359df3SBarry Smith   char tbuf[2];
672af359df3SBarry Smith   strcpy(tstr,"-thread");
673af359df3SBarry Smith   for(i=0;i<PetscMaxThreads;i++) {
674af359df3SBarry Smith     ThreadCoreAffinity[i] = i;
675af359df3SBarry Smith     sprintf(tbuf,"%d",i);
676af359df3SBarry Smith     strcat(tstr,tbuf);
677af359df3SBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
678af359df3SBarry Smith     if(flg1) {
679af359df3SBarry Smith       ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
680af359df3SBarry Smith       ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
681af359df3SBarry Smith     }
682af359df3SBarry Smith     tstr[7] = '\0';
683af359df3SBarry Smith   }
684cfcfc605SKerry Stevens 
685cfcfc605SKerry Stevens   /*
686cfcfc605SKerry Stevens       Determine whether to use thread pool
687cfcfc605SKerry Stevens    */
688cfcfc605SKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
689cfcfc605SKerry Stevens   if (flg1) {
690cfcfc605SKerry Stevens     PetscUseThreadPool = PETSC_TRUE;
691af359df3SBarry Smith     /* get the thread pool type */
692af359df3SBarry Smith     PetscInt ipool = 0;
693af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
694af359df3SBarry Smith 
695af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
696af359df3SBarry Smith     switch(ipool) {
697af359df3SBarry Smith     case 1:
698af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
699af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
700af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
701af359df3SBarry Smith       MainWait              = &MainWait_Tree;
702af359df3SBarry Smith       MainJob               = &MainJob_Tree;
703af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
704af359df3SBarry Smith       break;
705af359df3SBarry Smith     case 2:
706af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
707af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
708af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
709af359df3SBarry Smith       MainWait              = &MainWait_Main;
710af359df3SBarry Smith       MainJob               = &MainJob_Main;
711af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
712af359df3SBarry Smith       break;
713af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
714af359df3SBarry Smith     case 3:
715af359df3SBarry Smith #else
716af359df3SBarry Smith     default:
717af359df3SBarry Smith #endif
718af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
719af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
720af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
721af359df3SBarry Smith       MainWait              = &MainWait_Chain;
722af359df3SBarry Smith       MainJob               = &MainJob_Chain;
723af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
724af359df3SBarry Smith       break;
725af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
726af359df3SBarry Smith     default:
727af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
728af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
729af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
730af359df3SBarry Smith       MainWait              = &MainWait_True;
731af359df3SBarry Smith       MainJob               = &MainJob_True;
732af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
733af359df3SBarry Smith       break;
734af359df3SBarry Smith #endif
735af359df3SBarry Smith     }
736af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
737683509dcSKerry Stevens   } else {
738683509dcSKerry Stevens     //need to define these in the case on 'no threads' or 'thread create/destroy'
739683509dcSKerry Stevens     //could take any of the above versions
740683509dcSKerry Stevens     MainJob               = &MainJob_Spawn;
741af359df3SBarry Smith   }
742af359df3SBarry Smith #endif
743e5c89e4eSSatish Balay   /*
744e5c89e4eSSatish Balay        Print basic help message
745e5c89e4eSSatish Balay   */
746e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
747e5c89e4eSSatish Balay   if (flg1) {
748e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
749301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
750301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
751301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
752e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
753e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
754e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
755e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
756e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
757e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
758e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
759e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
760e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
761e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
762e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
763e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
764e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
765e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
766e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
767e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
768e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
769e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
770e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7714161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7724161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
773e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
774e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
775e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
776e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
777e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
778e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
779a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
780e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
78140ab9619SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -server <port>: Run PETSc webserver (default port is 8080) see PetscWebServe()\n");CHKERRQ(ierr);
782e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
783e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
7847ef452c0SMatthew G Knepley     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary _summary_python]: logging objects and events\n");CHKERRQ(ierr);
785e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
786e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
787e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
788e5c89e4eSSatish Balay #endif
7896cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
790e5c89e4eSSatish Balay #endif
791e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
792e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
793e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
794e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
795e5c89e4eSSatish Balay   }
796e5c89e4eSSatish Balay 
797a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
798e5c89e4eSSatish Balay   if (flg1) {
799e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
800e5c89e4eSSatish Balay   }
801e5c89e4eSSatish Balay 
8026cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
803e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
804e5c89e4eSSatish Balay   if (f) {
8056cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
806e5c89e4eSSatish Balay   }
807827f890bSBarry Smith 
8088154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
809c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
81073113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
81173113deaSBarry Smith   else flg1 = PETSC_FALSE;
8128154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
8138154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
814bab1f7e6SVictor Minden #endif
815192daf7cSBarry Smith 
816e5c89e4eSSatish Balay   PetscFunctionReturn(0);
817e5c89e4eSSatish Balay }
818df413903SBarry Smith 
819ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
820ba61063dSBarry Smith 
82151d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
82251d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
82351d315f7SKerry Stevens   PetscErrorCode iterr;
82451d315f7SKerry Stevens   int icorr,ierr;
82551d315f7SKerry Stevens   int* pId = (int*)arg;
82651d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
82751d315f7SKerry Stevens   PetscBool PeeOn;
82851d315f7SKerry Stevens   cpu_set_t mset;
8299e800a48SKerry Stevens   //printf("Thread %d In Tree Thread Function\n",ThreadId);
83051d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
83151d315f7SKerry Stevens   CPU_ZERO(&mset);
83251d315f7SKerry Stevens   CPU_SET(icorr,&mset);
83351d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
83451d315f7SKerry Stevens 
83551d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
83651d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
83751d315f7SKerry Stevens   }
83851d315f7SKerry Stevens   else {
83951d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
84051d315f7SKerry Stevens   }
84151d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
842ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
84351d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
84451d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
84551d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
84651d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
84751d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
848ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
849ba61063dSBarry Smith            upon return, has the lock */
85051d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
85151d315f7SKerry Stevens         }
85251d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
85351d315f7SKerry Stevens       }
85451d315f7SKerry Stevens     }
855ba61063dSBarry Smith     /* your subordinates are now ready */
85651d315f7SKerry Stevens   }
85751d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
858ba61063dSBarry Smith   /* update your ready status */
85951d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
86051d315f7SKerry Stevens   if(ThreadId==0) {
86151d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
862ba61063dSBarry Smith     /* ignal main */
86351d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
86451d315f7SKerry Stevens   }
86551d315f7SKerry Stevens   else {
866ba61063dSBarry Smith     /* tell your boss that you're ready to work */
86751d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
86851d315f7SKerry Stevens   }
869ba61063dSBarry Smith   /* the while loop needs to have an exit
870ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
871ba61063dSBarry Smith    and calling FuncFinish */
87251d315f7SKerry Stevens   while(PetscThreadGo) {
873ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
874ba61063dSBarry Smith       waiting when you don't have to causes problems
875ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
87651d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
877ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
878ba61063dSBarry Smith        upon return, has the lock */
87951d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
88051d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
88151d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
88251d315f7SKerry Stevens     }
88351d315f7SKerry Stevens     if(ThreadId==0) {
88451d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
88551d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
88651d315f7SKerry Stevens     }
88751d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
88851d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
889ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
89051d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
89151d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
89251d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
89351d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
89451d315f7SKerry Stevens         }
89551d315f7SKerry Stevens       }
89651d315f7SKerry Stevens     }
897ba61063dSBarry Smith     /* do your job */
89851d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
89951d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
90051d315f7SKerry Stevens     }
90151d315f7SKerry Stevens     else {
90251d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
90351d315f7SKerry Stevens     }
90451d315f7SKerry Stevens     if(iterr!=0) {
90551d315f7SKerry Stevens       ithreaderr = 1;
90651d315f7SKerry Stevens     }
90751d315f7SKerry Stevens     if(PetscThreadGo) {
908ba61063dSBarry Smith       /* reset job, get ready for more */
90951d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
910ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
911ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
91251d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
91351d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
91451d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
91551d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
91651d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
917ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
918ba61063dSBarry Smith                upon return, has the lock */
91951d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
92051d315f7SKerry Stevens             }
92151d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
92251d315f7SKerry Stevens           }
92351d315f7SKerry Stevens 	}
924ba61063dSBarry Smith         /* your subordinates are now ready */
92551d315f7SKerry Stevens       }
92651d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
92751d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
92851d315f7SKerry Stevens       if(ThreadId==0) {
929ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
930ba61063dSBarry Smith         /* root thread signals 'main' */
93151d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
93251d315f7SKerry Stevens       }
93351d315f7SKerry Stevens       else {
934ba61063dSBarry Smith         /* signal your boss before you go to sleep */
93551d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
93651d315f7SKerry Stevens       }
93751d315f7SKerry Stevens     }
93851d315f7SKerry Stevens   }
93951d315f7SKerry Stevens   return NULL;
94051d315f7SKerry Stevens }
94151d315f7SKerry Stevens 
94251d315f7SKerry Stevens #undef __FUNCT__
94351d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
94451d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
94551d315f7SKerry Stevens   PetscInt i,ierr;
94651d315f7SKerry Stevens   int status;
94751d315f7SKerry Stevens 
94851d315f7SKerry Stevens   if(PetscUseThreadPool) {
94951d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
95051d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
95151d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
95251d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
95351d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
95451d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
95551d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
95651d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
95751d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95851d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95951d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
96051d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
961ba61063dSBarry Smith     /* initialize job structure */
96251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
96351d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
96451d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
96551d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
96651d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
96751d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
96851d315f7SKerry Stevens     }
96951d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
97051d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
97151d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
97251d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
97351d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
97451d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
97551d315f7SKerry Stevens     }
97651d315f7SKerry Stevens     job_tree.pfunc = NULL;
97751d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
97851d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
97951d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
98051d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
981ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
98251d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
983ba61063dSBarry Smith     /* create threads */
98451d315f7SKerry Stevens     for(i=0; i<N; i++) {
98551d315f7SKerry Stevens       pVal[i] = i;
98651d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
987ba61063dSBarry Smith       /* should check status */
98851d315f7SKerry Stevens     }
98951d315f7SKerry Stevens   }
99051d315f7SKerry Stevens   return NULL;
99151d315f7SKerry Stevens }
99251d315f7SKerry Stevens 
99351d315f7SKerry Stevens #undef __FUNCT__
99451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
99551d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
99651d315f7SKerry Stevens   int i,ierr;
99751d315f7SKerry Stevens   void* jstatus;
99851d315f7SKerry Stevens 
99951d315f7SKerry Stevens   PetscFunctionBegin;
100051d315f7SKerry Stevens 
1001ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1002ba61063dSBarry Smith   /* join the threads */
100351d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
100451d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1005ba61063dSBarry Smith     /* do error checking*/
100651d315f7SKerry Stevens   }
100751d315f7SKerry Stevens   free(PetscThreadPoint);
100851d315f7SKerry Stevens   free(arrmutex);
100951d315f7SKerry Stevens   free(arrcond1);
101051d315f7SKerry Stevens   free(arrcond2);
101151d315f7SKerry Stevens   free(arrstart);
101251d315f7SKerry Stevens   free(arrready);
101351d315f7SKerry Stevens   free(job_tree.pdata);
101451d315f7SKerry Stevens   free(pVal);
1015cfcfc605SKerry Stevens 
101651d315f7SKerry Stevens   PetscFunctionReturn(0);
101751d315f7SKerry Stevens }
101851d315f7SKerry Stevens 
101951d315f7SKerry Stevens #undef __FUNCT__
102051d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
102151d315f7SKerry Stevens void MainWait_Tree() {
102251d315f7SKerry Stevens   int ierr;
102351d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
102451d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
102551d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
102651d315f7SKerry Stevens   }
102751d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
102851d315f7SKerry Stevens }
102951d315f7SKerry Stevens 
103051d315f7SKerry Stevens #undef __FUNCT__
103151d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
103251d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
103351d315f7SKerry Stevens   int i,ierr;
103451d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
103536d20dc5SKerry Stevens 
103651d315f7SKerry Stevens   MainWait();
103751d315f7SKerry Stevens   job_tree.pfunc = pFunc;
103851d315f7SKerry Stevens   job_tree.pdata = data;
103951d315f7SKerry Stevens   job_tree.startJob = PETSC_TRUE;
104051d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
104151d315f7SKerry Stevens     *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
104251d315f7SKerry Stevens   }
104351d315f7SKerry Stevens   job_tree.eJobStat = JobInitiated;
104451d315f7SKerry Stevens   ierr = pthread_cond_signal(job_tree.cond2array[0]);
104551d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1046ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
104751d315f7SKerry Stevens   }
104836d20dc5SKerry Stevens 
104951d315f7SKerry Stevens   if(ithreaderr) {
105051d315f7SKerry Stevens     ijoberr = ithreaderr;
105151d315f7SKerry Stevens   }
105251d315f7SKerry Stevens   return ijoberr;
105351d315f7SKerry Stevens }
105451d315f7SKerry Stevens /****  ****/
105551d315f7SKerry Stevens 
105651d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
105751d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
105851d315f7SKerry Stevens   PetscErrorCode iterr;
105951d315f7SKerry Stevens   int icorr,ierr;
106051d315f7SKerry Stevens   int* pId = (int*)arg;
106151d315f7SKerry Stevens   int ThreadId = *pId;
106251d315f7SKerry Stevens   cpu_set_t mset;
10639e800a48SKerry Stevens   //printf("Thread %d In Main Thread Function\n",ThreadId);
106451d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
106551d315f7SKerry Stevens   CPU_ZERO(&mset);
106651d315f7SKerry Stevens   CPU_SET(icorr,&mset);
106751d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
106851d315f7SKerry Stevens 
106951d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1070ba61063dSBarry Smith   /* update your ready status */
107151d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1072ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
107351d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
107451d315f7SKerry Stevens 
1075ba61063dSBarry Smith   /* the while loop needs to have an exit
1076ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1077ba61063dSBarry Smith      and calling FuncFinish */
107851d315f7SKerry Stevens   while(PetscThreadGo) {
1079ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1080ba61063dSBarry Smith        waiting when you don't have to causes problems
1081ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
108251d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1083ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1084ba61063dSBarry Smith        upon return, has the lock */
108551d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1086ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
108751d315f7SKerry Stevens     }
108851d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
108951d315f7SKerry Stevens     if(job_main.pdata==NULL) {
109051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
109151d315f7SKerry Stevens     }
109251d315f7SKerry Stevens     else {
109351d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
109451d315f7SKerry Stevens     }
109551d315f7SKerry Stevens     if(iterr!=0) {
109651d315f7SKerry Stevens       ithreaderr = 1;
109751d315f7SKerry Stevens     }
109851d315f7SKerry Stevens     if(PetscThreadGo) {
1099ba61063dSBarry Smith       /* reset job, get ready for more */
110051d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
110151d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1102ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
110351d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
110451d315f7SKerry Stevens     }
110551d315f7SKerry Stevens   }
110651d315f7SKerry Stevens   return NULL;
110751d315f7SKerry Stevens }
110851d315f7SKerry Stevens 
110951d315f7SKerry Stevens #undef __FUNCT__
111051d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
111151d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
111251d315f7SKerry Stevens   PetscInt i,ierr;
111351d315f7SKerry Stevens   int status;
111451d315f7SKerry Stevens 
111551d315f7SKerry Stevens   if(PetscUseThreadPool) {
111651d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
111751d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
111851d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
111951d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
112051d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
112151d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
112251d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
112351d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
112451d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112551d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112651d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1127ba61063dSBarry Smith     /* initialize job structure */
112851d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
112951d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
113051d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
113151d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
113251d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
113351d315f7SKerry Stevens     }
113451d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
113551d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
113651d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
113751d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
113851d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
113951d315f7SKerry Stevens     }
114051d315f7SKerry Stevens     job_main.pfunc = NULL;
114151d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
114251d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1143ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
114451d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1145ba61063dSBarry Smith     /* create threads */
114651d315f7SKerry Stevens     for(i=0; i<N; i++) {
114751d315f7SKerry Stevens       pVal[i] = i;
114851d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1149ba61063dSBarry Smith       /* error check */
115051d315f7SKerry Stevens     }
115151d315f7SKerry Stevens   }
115251d315f7SKerry Stevens   else {
115351d315f7SKerry Stevens   }
115451d315f7SKerry Stevens   return NULL;
115551d315f7SKerry Stevens }
115651d315f7SKerry Stevens 
115751d315f7SKerry Stevens #undef __FUNCT__
115851d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
115951d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
116051d315f7SKerry Stevens   int i,ierr;
116151d315f7SKerry Stevens   void* jstatus;
116251d315f7SKerry Stevens 
116351d315f7SKerry Stevens   PetscFunctionBegin;
116451d315f7SKerry Stevens 
1165ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1166ba61063dSBarry Smith   /* join the threads */
116751d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1168ba61063dSBarry Smith     ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
116951d315f7SKerry Stevens   }
117051d315f7SKerry Stevens   free(PetscThreadPoint);
117151d315f7SKerry Stevens   free(arrmutex);
117251d315f7SKerry Stevens   free(arrcond1);
117351d315f7SKerry Stevens   free(arrcond2);
117451d315f7SKerry Stevens   free(arrstart);
117551d315f7SKerry Stevens   free(arrready);
117651d315f7SKerry Stevens   free(job_main.pdata);
117751d315f7SKerry Stevens   free(pVal);
1178cfcfc605SKerry Stevens 
117951d315f7SKerry Stevens   PetscFunctionReturn(0);
118051d315f7SKerry Stevens }
118151d315f7SKerry Stevens 
118251d315f7SKerry Stevens #undef __FUNCT__
118351d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
118451d315f7SKerry Stevens void MainWait_Main() {
118551d315f7SKerry Stevens   int i,ierr;
118651d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
118751d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
118851d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
118951d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
119051d315f7SKerry Stevens     }
119151d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
119251d315f7SKerry Stevens   }
119351d315f7SKerry Stevens }
119451d315f7SKerry Stevens 
119551d315f7SKerry Stevens #undef __FUNCT__
119651d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
119751d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
119851d315f7SKerry Stevens   int i,ierr;
119951d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
120036d20dc5SKerry Stevens 
1201ba61063dSBarry Smith   MainWait(); /* you know everyone is waiting to be signalled! */
120251d315f7SKerry Stevens   job_main.pfunc = pFunc;
120351d315f7SKerry Stevens   job_main.pdata = data;
120451d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1205ba61063dSBarry Smith     *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
120651d315f7SKerry Stevens   }
1207ba61063dSBarry Smith   /* tell the threads to go to work */
120851d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
120951d315f7SKerry Stevens     ierr = pthread_cond_signal(job_main.cond2array[i]);
121051d315f7SKerry Stevens   }
121151d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1212ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
121351d315f7SKerry Stevens   }
121436d20dc5SKerry Stevens 
121551d315f7SKerry Stevens   if(ithreaderr) {
121651d315f7SKerry Stevens     ijoberr = ithreaderr;
121751d315f7SKerry Stevens   }
121851d315f7SKerry Stevens   return ijoberr;
121951d315f7SKerry Stevens }
122051d315f7SKerry Stevens /****  ****/
122151d315f7SKerry Stevens 
122251d315f7SKerry Stevens /**** Chain Thread Functions ****/
122351d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
122451d315f7SKerry Stevens   PetscErrorCode iterr;
122551d315f7SKerry Stevens   int icorr,ierr;
122651d315f7SKerry Stevens   int* pId = (int*)arg;
122751d315f7SKerry Stevens   int ThreadId = *pId;
122851d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
122951d315f7SKerry Stevens   PetscBool PeeOn;
123051d315f7SKerry Stevens   cpu_set_t mset;
12319e800a48SKerry Stevens   //printf("Thread %d In Chain Thread Function\n",ThreadId);
123251d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
123351d315f7SKerry Stevens   CPU_ZERO(&mset);
123451d315f7SKerry Stevens   CPU_SET(icorr,&mset);
123551d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
123651d315f7SKerry Stevens 
123751d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
123851d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
123951d315f7SKerry Stevens   }
124051d315f7SKerry Stevens   else {
124151d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
124251d315f7SKerry Stevens   }
124351d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1244ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
124551d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
124651d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1247ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1248ba61063dSBarry Smith        upon return, has the lock */
124951d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
125051d315f7SKerry Stevens     }
125151d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1252ba61063dSBarry Smith     /* your subordinate is now ready*/
125351d315f7SKerry Stevens   }
125451d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1255ba61063dSBarry Smith   /* update your ready status */
125651d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
125751d315f7SKerry Stevens   if(ThreadId==0) {
125851d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1259ba61063dSBarry Smith     /* signal main */
126051d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
126151d315f7SKerry Stevens   }
126251d315f7SKerry Stevens   else {
1263ba61063dSBarry Smith     /* tell your boss that you're ready to work */
126451d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
126551d315f7SKerry Stevens   }
1266ba61063dSBarry Smith   /*  the while loop needs to have an exit
1267ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1268ba61063dSBarry Smith    and calling FuncFinish */
126951d315f7SKerry Stevens   while(PetscThreadGo) {
1270ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1271ba61063dSBarry Smith        waiting when you don't have to causes problems
1272ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
127351d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1274ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1275ba61063dSBarry Smith        upon return, has the lock */
127651d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
127751d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
127851d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
127951d315f7SKerry Stevens     }
128051d315f7SKerry Stevens     if(ThreadId==0) {
128151d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
128251d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
128351d315f7SKerry Stevens     }
128451d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
128551d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1286ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
128751d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
128851d315f7SKerry Stevens     }
1289ba61063dSBarry Smith     /* do your job */
129051d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
129151d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
129251d315f7SKerry Stevens     }
129351d315f7SKerry Stevens     else {
129451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
129551d315f7SKerry Stevens     }
129651d315f7SKerry Stevens     if(iterr!=0) {
129751d315f7SKerry Stevens       ithreaderr = 1;
129851d315f7SKerry Stevens     }
129951d315f7SKerry Stevens     if(PetscThreadGo) {
1300ba61063dSBarry Smith       /* reset job, get ready for more */
130151d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1302ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1303ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
130451d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
130551d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1306ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1307ba61063dSBarry Smith            upon return, has the lock */
130851d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
130951d315f7SKerry Stevens         }
131051d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1311ba61063dSBarry Smith         /* your subordinate is now ready */
131251d315f7SKerry Stevens       }
131351d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
131451d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
131551d315f7SKerry Stevens       if(ThreadId==0) {
1316ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1317ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
131851d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
131951d315f7SKerry Stevens       }
132051d315f7SKerry Stevens       else {
1321ba61063dSBarry Smith         /* signal your boss before you go to sleep */
132251d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
132351d315f7SKerry Stevens       }
132451d315f7SKerry Stevens     }
132551d315f7SKerry Stevens   }
132651d315f7SKerry Stevens   return NULL;
132751d315f7SKerry Stevens }
132851d315f7SKerry Stevens 
132951d315f7SKerry Stevens #undef __FUNCT__
133051d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
133151d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
133251d315f7SKerry Stevens   PetscInt i,ierr;
133351d315f7SKerry Stevens   int status;
133451d315f7SKerry Stevens 
133551d315f7SKerry Stevens   if(PetscUseThreadPool) {
133651d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
133751d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
133851d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
133951d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
134051d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
134151d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
134251d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
134351d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
134451d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134551d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134651d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
134751d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1348ba61063dSBarry Smith     /* initialize job structure */
134951d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
135051d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
135151d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
135251d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
135351d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
135451d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
135551d315f7SKerry Stevens     }
135651d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
135751d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
135851d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
135951d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
136051d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
136151d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
136251d315f7SKerry Stevens     }
136351d315f7SKerry Stevens     job_chain.pfunc = NULL;
136451d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
136551d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
136651d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
136751d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1368ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
136951d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1370ba61063dSBarry Smith     /* create threads */
137151d315f7SKerry Stevens     for(i=0; i<N; i++) {
137251d315f7SKerry Stevens       pVal[i] = i;
137351d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1374ba61063dSBarry Smith       /* should check error */
137551d315f7SKerry Stevens     }
137651d315f7SKerry Stevens   }
137751d315f7SKerry Stevens   else {
137851d315f7SKerry Stevens   }
137951d315f7SKerry Stevens   return NULL;
138051d315f7SKerry Stevens }
138151d315f7SKerry Stevens 
138251d315f7SKerry Stevens 
138351d315f7SKerry Stevens #undef __FUNCT__
138451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
138551d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
138651d315f7SKerry Stevens   int i,ierr;
138751d315f7SKerry Stevens   void* jstatus;
138851d315f7SKerry Stevens 
138951d315f7SKerry Stevens   PetscFunctionBegin;
139051d315f7SKerry Stevens 
1391ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1392ba61063dSBarry Smith   /* join the threads */
139351d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
139451d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1395ba61063dSBarry Smith     /* should check error */
139651d315f7SKerry Stevens   }
139751d315f7SKerry Stevens   free(PetscThreadPoint);
139851d315f7SKerry Stevens   free(arrmutex);
139951d315f7SKerry Stevens   free(arrcond1);
140051d315f7SKerry Stevens   free(arrcond2);
140151d315f7SKerry Stevens   free(arrstart);
140251d315f7SKerry Stevens   free(arrready);
140351d315f7SKerry Stevens   free(job_chain.pdata);
140451d315f7SKerry Stevens   free(pVal);
1405cfcfc605SKerry Stevens 
140651d315f7SKerry Stevens   PetscFunctionReturn(0);
140751d315f7SKerry Stevens }
140851d315f7SKerry Stevens 
140951d315f7SKerry Stevens #undef __FUNCT__
141051d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
141151d315f7SKerry Stevens void MainWait_Chain() {
141251d315f7SKerry Stevens   int ierr;
141351d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
141451d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
141551d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
141651d315f7SKerry Stevens   }
141751d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
141851d315f7SKerry Stevens }
141951d315f7SKerry Stevens 
142051d315f7SKerry Stevens #undef __FUNCT__
142151d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
142251d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
142351d315f7SKerry Stevens   int i,ierr;
142451d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
142536d20dc5SKerry Stevens 
142651d315f7SKerry Stevens   MainWait();
142751d315f7SKerry Stevens   job_chain.pfunc = pFunc;
142851d315f7SKerry Stevens   job_chain.pdata = data;
142951d315f7SKerry Stevens   job_chain.startJob = PETSC_TRUE;
143051d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
143151d315f7SKerry Stevens     *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
143251d315f7SKerry Stevens   }
143351d315f7SKerry Stevens   job_chain.eJobStat = JobInitiated;
143451d315f7SKerry Stevens   ierr = pthread_cond_signal(job_chain.cond2array[0]);
143551d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1436ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
143751d315f7SKerry Stevens   }
143836d20dc5SKerry 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;
14549e800a48SKerry Stevens   //printf("Thread %d In True Pool Thread Function\n",ThreadId);
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 */
14756c9b208dSKerry Stevens       //printf("Thread %d Going to Sleep!\n",ThreadId);
147651d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
147751dcc849SKerry Stevens     }
147851d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
147951d315f7SKerry Stevens     job_true.iNumJobThreads--;
148051d315f7SKerry Stevens     job_true.iNumReadyThreads--;
148151d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
148251d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
148351d315f7SKerry Stevens     if(job_true.pdata==NULL) {
148451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
148551dcc849SKerry Stevens     }
148651dcc849SKerry Stevens     else {
148751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
148851dcc849SKerry Stevens     }
14890ca81413SKerry Stevens     if(iterr!=0) {
14900ca81413SKerry Stevens       ithreaderr = 1;
14910ca81413SKerry Stevens     }
14926c9b208dSKerry Stevens     //printf("Thread %d Finished Job\n",ThreadId);
1493ba61063dSBarry Smith     /* the barrier is necessary BECAUSE: look at job_true.iNumReadyThreads
1494ba61063dSBarry Smith       what happens if a thread finishes before they all start? BAD!
1495ba61063dSBarry Smith      what happens if a thread finishes before any else start? BAD! */
1496ba61063dSBarry Smith     pthread_barrier_wait(job_true.pbarr); /* ensures all threads are finished */
1497ba61063dSBarry Smith     /* reset job */
149851dcc849SKerry Stevens     if(PetscThreadGo) {
149951d315f7SKerry Stevens       pthread_mutex_lock(&job_true.mutex);
150051d315f7SKerry Stevens       job_true.iNumReadyThreads++;
150151d315f7SKerry Stevens       if(job_true.iNumReadyThreads==PetscMaxThreads) {
1502ba61063dSBarry Smith 	/* signal the 'main' thread that the job is done! (only done once) */
150351dcc849SKerry Stevens 	ierr = pthread_cond_signal(&main_cond);
150451dcc849SKerry Stevens       }
150551dcc849SKerry Stevens     }
150651dcc849SKerry Stevens   }
150751dcc849SKerry Stevens   return NULL;
150851dcc849SKerry Stevens }
150951dcc849SKerry Stevens 
1510f09cb4aaSKerry Stevens #undef __FUNCT__
151151d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_True"
151251d315f7SKerry Stevens void* PetscThreadInitialize_True(PetscInt N) {
151351dcc849SKerry Stevens   PetscInt i;
151451dcc849SKerry Stevens   int status;
15150ca81413SKerry Stevens 
1516f09cb4aaSKerry Stevens   pVal = (int*)malloc(N*sizeof(int));
1517ba61063dSBarry Smith   /* allocate memory in the heap for the thread structure */
151851dcc849SKerry Stevens   PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1519ba61063dSBarry Smith   BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
152051d315f7SKerry Stevens   job_true.pdata = (void**)malloc(N*sizeof(void*));
152151dcc849SKerry Stevens   for(i=0; i<N; i++) {
1522f09cb4aaSKerry Stevens     pVal[i] = i;
1523f09cb4aaSKerry Stevens     status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1524ba61063dSBarry Smith     /* error check to ensure proper thread creation */
152551dcc849SKerry Stevens     status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1526ba61063dSBarry Smith     /* should check error */
152751dcc849SKerry Stevens   }
15286c9b208dSKerry Stevens   //printf("Finished True Thread Pool Initialization\n");
152951dcc849SKerry Stevens   return NULL;
153051dcc849SKerry Stevens }
153151dcc849SKerry Stevens 
1532f09cb4aaSKerry Stevens 
1533f09cb4aaSKerry Stevens #undef __FUNCT__
153451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
153551d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
153651dcc849SKerry Stevens   int i,ierr;
153751dcc849SKerry Stevens   void* jstatus;
153851dcc849SKerry Stevens 
153951dcc849SKerry Stevens   PetscFunctionBegin;
1540cfcfc605SKerry Stevens 
1541ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1542ba61063dSBarry Smith   /* join the threads */
154351dcc849SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
154451dcc849SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
154551dcc849SKerry Stevens   }
154651dcc849SKerry Stevens   free(BarrPoint);
154751dcc849SKerry Stevens   free(PetscThreadPoint);
1548cfcfc605SKerry Stevens 
154951dcc849SKerry Stevens   PetscFunctionReturn(0);
155051dcc849SKerry Stevens }
155151dcc849SKerry Stevens 
1552f09cb4aaSKerry Stevens #undef __FUNCT__
155351d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
155451d315f7SKerry Stevens void MainWait_True() {
155551dcc849SKerry Stevens   int ierr;
15566c9b208dSKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
155751d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
155851d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
155951dcc849SKerry Stevens   }
156051d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
156151dcc849SKerry Stevens }
156251dcc849SKerry Stevens 
1563f09cb4aaSKerry Stevens #undef __FUNCT__
156451d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
156551d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
156651dcc849SKerry Stevens   int ierr;
15670ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
156836d20dc5SKerry Stevens 
15690ca81413SKerry Stevens   MainWait();
157051d315f7SKerry Stevens   job_true.pfunc = pFunc;
157151d315f7SKerry Stevens   job_true.pdata = data;
157251d315f7SKerry Stevens   job_true.pbarr = &BarrPoint[n];
157351d315f7SKerry Stevens   job_true.iNumJobThreads = n;
157451d315f7SKerry Stevens   job_true.startJob = PETSC_TRUE;
157551d315f7SKerry Stevens   ierr = pthread_cond_broadcast(&job_true.cond);
15760ca81413SKerry Stevens   if(pFunc!=FuncFinish) {
1577ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done */
15780ca81413SKerry Stevens   }
157936d20dc5SKerry Stevens 
15800ca81413SKerry Stevens   if(ithreaderr) {
15810ca81413SKerry Stevens     ijoberr = ithreaderr;
15820ca81413SKerry Stevens   }
15830ca81413SKerry Stevens   return ijoberr;
158451dcc849SKerry Stevens }
1585683509dcSKerry Stevens /**** NO THREAD POOL FUNCTION ****/
1586683509dcSKerry Stevens #undef __FUNCT__
1587683509dcSKerry Stevens #define __FUNCT__ "MainJob_Spawn"
1588683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void** data,PetscInt n) {
1589683509dcSKerry Stevens   PetscErrorCode ijoberr = 0;
1590683509dcSKerry Stevens 
1591683509dcSKerry Stevens   pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
1592cfcfc605SKerry Stevens   PetscThreadPoint = apThread; /* point to same place */
1593683509dcSKerry Stevens   PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1594683509dcSKerry Stevens   PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
1595683509dcSKerry Stevens   free(apThread);
1596683509dcSKerry Stevens 
1597683509dcSKerry Stevens   return ijoberr;
1598683509dcSKerry Stevens }
159951d315f7SKerry Stevens /****  ****/
1600ba61063dSBarry Smith #endif
160151dcc849SKerry Stevens 
160251dcc849SKerry Stevens void* FuncFinish(void* arg) {
160351dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16040ca81413SKerry Stevens   return(0);
160551dcc849SKerry Stevens }
1606ba61063dSBarry Smith 
1607ba61063dSBarry Smith #endif
1608