xref: /petsc/src/sys/objects/init.c (revision 24edbb7e7869e61eb26a18dba63dc481248fd073)
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)
15*24edbb7eSSatish Balay #ifndef __USE_GNU
16*24edbb7eSSatish Balay #define __USE_GNU
17*24edbb7eSSatish Balay #endif
188f54b378SBarry Smith #include <sched.h>
19cd723cd1SBarry Smith #endif
20bda8bf91SBarry Smith #if defined(PETSC_HAVE_PTHREAD_H)
2151dcc849SKerry Stevens #include <pthread.h>
22ba61063dSBarry Smith #endif
23ef386f4bSSatish Balay 
24ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
2551d315f7SKerry Stevens #include <sys/sysinfo.h>
26ba61063dSBarry Smith #endif
27121deb67SSatish Balay #if defined(PETSC_HAVE_UNISTD_H)
2851d315f7SKerry Stevens #include <unistd.h>
29121deb67SSatish Balay #endif
30e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
31e5c89e4eSSatish Balay #include <stdlib.h>
32e5c89e4eSSatish Balay #endif
33e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
34e5c89e4eSSatish Balay #include <malloc.h>
35e5c89e4eSSatish Balay #endif
36555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
37555d055bSBarry Smith #include <valgrind/valgrind.h>
38555d055bSBarry Smith #endif
39555d055bSBarry Smith 
40e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
41e5c89e4eSSatish Balay /*
42e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
43e5c89e4eSSatish Balay    already started before PETSc was initialized.
44e5c89e4eSSatish Balay */
457087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
467087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
477087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
4851dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
4951dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
507087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
517087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
52ba61063dSBarry Smith 
53ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
5451dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
5551dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
56af359df3SBarry Smith #define PETSC_HAVE_PTHREAD_BARRIER
57ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
58ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
59ba61063dSBarry Smith #endif
6051d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
61f09cb4aaSKerry Stevens int*         pVal;
6251dcc849SKerry Stevens 
63ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
6451d315f7SKerry Stevens int* ThreadCoreAffinity;
6551d315f7SKerry Stevens 
66ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
6751d315f7SKerry Stevens 
6851d315f7SKerry Stevens typedef struct {
6951d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
7051d315f7SKerry Stevens   pthread_cond_t**  cond1array;
7151d315f7SKerry Stevens   pthread_cond_t** cond2array;
7251d315f7SKerry Stevens   void* (*pfunc)(void*);
7351d315f7SKerry Stevens   void** pdata;
7451d315f7SKerry Stevens   PetscBool startJob;
7551d315f7SKerry Stevens   estat eJobStat;
7651d315f7SKerry Stevens   PetscBool** arrThreadStarted;
7751d315f7SKerry Stevens   PetscBool** arrThreadReady;
7851d315f7SKerry Stevens } sjob_tree;
7951d315f7SKerry Stevens sjob_tree job_tree;
8051d315f7SKerry Stevens typedef struct {
8151d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8251d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8351d315f7SKerry Stevens   pthread_cond_t** cond2array;
8451d315f7SKerry Stevens   void* (*pfunc)(void*);
8551d315f7SKerry Stevens   void** pdata;
8651d315f7SKerry Stevens   PetscBool** arrThreadReady;
8751d315f7SKerry Stevens } sjob_main;
8851d315f7SKerry Stevens sjob_main job_main;
8951d315f7SKerry Stevens typedef struct {
9051d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
9151d315f7SKerry Stevens   pthread_cond_t**  cond1array;
9251d315f7SKerry Stevens   pthread_cond_t** cond2array;
9351d315f7SKerry Stevens   void* (*pfunc)(void*);
9451d315f7SKerry Stevens   void** pdata;
9551d315f7SKerry Stevens   PetscBool startJob;
9651d315f7SKerry Stevens   estat eJobStat;
9751d315f7SKerry Stevens   PetscBool** arrThreadStarted;
9851d315f7SKerry Stevens   PetscBool** arrThreadReady;
9951d315f7SKerry Stevens } sjob_chain;
10051d315f7SKerry Stevens sjob_chain job_chain;
101ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
10251dcc849SKerry Stevens typedef struct {
10351dcc849SKerry Stevens   pthread_mutex_t mutex;
10451dcc849SKerry Stevens   pthread_cond_t cond;
10551dcc849SKerry Stevens   void* (*pfunc)(void*);
10651dcc849SKerry Stevens   void** pdata;
10751dcc849SKerry Stevens   pthread_barrier_t* pbarr;
10851dcc849SKerry Stevens   int iNumJobThreads;
10951dcc849SKerry Stevens   int iNumReadyThreads;
11051dcc849SKerry Stevens   PetscBool startJob;
11151d315f7SKerry Stevens } sjob_true;
11251d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
113ba61063dSBarry Smith #endif
11451dcc849SKerry Stevens 
115ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
116ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
117ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
118ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
119ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
120ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
12151dcc849SKerry Stevens 
12251d315f7SKerry Stevens /* Function Pointers */
12351d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
12451d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
12551d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
12651d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
12751d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
12836d20dc5SKerry Stevens /**** Tree Thread Pool Functions ****/
12951d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
13051d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
13151d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
13251d315f7SKerry Stevens void           MainWait_Tree(void);
13351d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
13436d20dc5SKerry Stevens /**** Main Thread Pool Functions ****/
13551d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
13651d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
13751d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
13851d315f7SKerry Stevens void           MainWait_Main(void);
13951d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
14036d20dc5SKerry Stevens /**** Chain Thread Pool Functions ****/
14151d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
14251d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
14351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
14451d315f7SKerry Stevens void           MainWait_Chain(void);
14551d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
14636d20dc5SKerry Stevens /**** True Thread Pool Functions ****/
14751d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
14851d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
14951d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
15051d315f7SKerry Stevens void           MainWait_True(void);
15151d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
15236d20dc5SKerry Stevens /**** NO Thread Pool Function  ****/
153683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void**,PetscInt);
15436d20dc5SKerry Stevens /****  ****/
15551dcc849SKerry Stevens void* FuncFinish(void*);
1560ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1570ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
158ba61063dSBarry Smith #endif
159e5c89e4eSSatish Balay 
160e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
161e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
162e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
163e5c89e4eSSatish Balay #endif
1642c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
165500d8756SSatish Balay MPI_Datatype   MPIU_C_DOUBLE_COMPLEX;
166500d8756SSatish Balay MPI_Datatype   MPIU_C_COMPLEX;
1672c876bd9SBarry Smith #endif
1687087cfbeSBarry Smith PetscScalar    PETSC_i;
169e5c89e4eSSatish Balay #else
1707087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
171e5c89e4eSSatish Balay #endif
172ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
173c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
174c90a1750SBarry Smith #endif
1757087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1767087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
17775567043SBarry Smith 
178e5c89e4eSSatish Balay /*
179e5c89e4eSSatish Balay      These are needed by petscbt.h
180e5c89e4eSSatish Balay */
181c6db04a5SJed Brown #include <petscbt.h>
1827087cfbeSBarry Smith char      _BT_mask = ' ';
1837087cfbeSBarry Smith char      _BT_c = ' ';
1847087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
185e5c89e4eSSatish Balay 
186e5c89e4eSSatish Balay /*
187e5c89e4eSSatish Balay        Function that is called to display all error messages
188e5c89e4eSSatish Balay */
1897087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1907087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
191238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1927087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
193238ccf28SShri Abhyankar #else
1947087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
195238ccf28SShri Abhyankar #endif
196bab1f7e6SVictor Minden /*
1978154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1988154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
199bab1f7e6SVictor Minden 
200e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
201e5c89e4eSSatish Balay /*
202e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
203e5c89e4eSSatish Balay */
204e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
205e5c89e4eSSatish Balay 
206e5c89e4eSSatish Balay #undef __FUNCT__
207f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
2087087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
209e5c89e4eSSatish Balay {
210e5c89e4eSSatish Balay   PetscErrorCode ierr;
211e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
212e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
213e5c89e4eSSatish Balay   char           version[256];
214e5c89e4eSSatish Balay 
215e5c89e4eSSatish Balay   PetscFunctionBegin;
216e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
217e5c89e4eSSatish Balay   if (!rank) {
218e5c89e4eSSatish Balay     char        arch[10];
219f56c2debSBarry Smith     int         err;
22088c29154SBarry Smith     PetscViewer viewer;
221f56c2debSBarry Smith 
222e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
223e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
224a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
225e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
226e5c89e4eSSatish Balay     if (filename) {
227e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
228e5c89e4eSSatish Balay     } else {
229e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
230e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
231e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
232e5c89e4eSSatish Balay     }
233e5c89e4eSSatish Balay 
234e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
235e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
236e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
237e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
238e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
23988c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
24088c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2416bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
242e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
243f56c2debSBarry Smith     err = fflush(*fd);
244e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
245e5c89e4eSSatish Balay   }
246e5c89e4eSSatish Balay   PetscFunctionReturn(0);
247e5c89e4eSSatish Balay }
248e5c89e4eSSatish Balay 
249e5c89e4eSSatish Balay #undef __FUNCT__
250f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2517087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
252e5c89e4eSSatish Balay {
253e5c89e4eSSatish Balay   PetscErrorCode ierr;
254e5c89e4eSSatish Balay   PetscMPIInt    rank;
255e5c89e4eSSatish Balay   char           date[64];
256f56c2debSBarry Smith   int            err;
257e5c89e4eSSatish Balay 
258e5c89e4eSSatish Balay   PetscFunctionBegin;
259e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
260e5c89e4eSSatish Balay   if (!rank) {
261e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
262e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
263e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
264e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
265f56c2debSBarry Smith     err = fflush(*fd);
266e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
267f56c2debSBarry Smith     err = fclose(*fd);
268e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
269e5c89e4eSSatish Balay   }
270e5c89e4eSSatish Balay   PetscFunctionReturn(0);
271e5c89e4eSSatish Balay }
272e5c89e4eSSatish Balay 
273e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
274e5c89e4eSSatish Balay 
275e5c89e4eSSatish Balay /*
276e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
277e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
278e5c89e4eSSatish Balay 
279e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2803c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2813c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
282e5c89e4eSSatish Balay */
283e5c89e4eSSatish Balay 
284e5c89e4eSSatish Balay #undef __FUNCT__
285e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
286e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
287e5c89e4eSSatish Balay {
288e5c89e4eSSatish Balay   PetscFunctionBegin;
2893c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
290e5c89e4eSSatish Balay   abort();
291e5c89e4eSSatish Balay }
292e5c89e4eSSatish Balay 
293e5c89e4eSSatish Balay #undef __FUNCT__
294e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
295e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
296e5c89e4eSSatish Balay {
297e5c89e4eSSatish Balay   PetscErrorCode ierr;
298e5c89e4eSSatish Balay 
299e5c89e4eSSatish Balay   PetscFunctionBegin;
3003c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
301e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
302e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
3033c311c98SBarry Smith     MPI_Abort(*comm,*flag);
304e5c89e4eSSatish Balay   }
305e5c89e4eSSatish Balay }
306e5c89e4eSSatish Balay 
307e5c89e4eSSatish Balay #undef __FUNCT__
308e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
309e5c89e4eSSatish Balay /*@C
310e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
311e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
312e5c89e4eSSatish Balay 
313e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
314e5c89e4eSSatish Balay 
315e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
316e5c89e4eSSatish Balay 
317e5c89e4eSSatish Balay    Level: advanced
318e5c89e4eSSatish Balay 
319e5c89e4eSSatish Balay    Note:
320e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
321e5c89e4eSSatish Balay 
32288c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
323e5c89e4eSSatish Balay @*/
3247087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
325e5c89e4eSSatish Balay {
326e5c89e4eSSatish Balay   PetscFunctionBegin;
327e5c89e4eSSatish Balay   PetscFinalize();
328e5c89e4eSSatish Balay   exit(0);
329e5c89e4eSSatish Balay   return 0;
330e5c89e4eSSatish Balay }
331e5c89e4eSSatish Balay 
332ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
33309573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
334ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
335e5c89e4eSSatish Balay static char       emacsmachinename[256];
336e5c89e4eSSatish Balay 
337e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
338e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
339e5c89e4eSSatish Balay 
340e5c89e4eSSatish Balay #undef __FUNCT__
341e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
342e5c89e4eSSatish Balay /*@C
343e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
344e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
345e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
346e5c89e4eSSatish Balay 
347e5c89e4eSSatish Balay    Input Parameter:
348e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
349da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
350e5c89e4eSSatish Balay 
351e5c89e4eSSatish Balay    Level: developer
352e5c89e4eSSatish Balay 
353e5c89e4eSSatish Balay    Concepts: package help message
354e5c89e4eSSatish Balay 
355e5c89e4eSSatish Balay @*/
3567087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
357e5c89e4eSSatish Balay {
358e5c89e4eSSatish Balay   PetscFunctionBegin;
359e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
360e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
361e5c89e4eSSatish Balay   PetscFunctionReturn(0);
362e5c89e4eSSatish Balay }
363e5c89e4eSSatish Balay 
364e5c89e4eSSatish Balay #undef __FUNCT__
365e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3667087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
367e5c89e4eSSatish Balay {
368e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
369e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
370ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
371e5c89e4eSSatish Balay   PetscErrorCode ierr;
372a6d0e24fSJed Brown   PetscReal      si;
373e5c89e4eSSatish Balay   int            i;
374e5c89e4eSSatish Balay   PetscMPIInt    rank;
375e5c89e4eSSatish Balay   char           version[256];
376e5c89e4eSSatish Balay 
377e5c89e4eSSatish Balay   PetscFunctionBegin;
378e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
379e5c89e4eSSatish Balay 
380e5c89e4eSSatish Balay   /*
381e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
382e5c89e4eSSatish Balay   */
3838bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
38481b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
385acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
386e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
387555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
388555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
389555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
390555d055bSBarry Smith #endif
391e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
392555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
393555d055bSBarry Smith     }
394555d055bSBarry Smith #endif
395e5c89e4eSSatish Balay   }
396e5c89e4eSSatish Balay #else
397acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
398acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
399e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
400e5c89e4eSSatish Balay #endif
401e5c89e4eSSatish Balay   if (flg3) {
402e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
403e5c89e4eSSatish Balay   }
40490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
405acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
406e5c89e4eSSatish Balay   if (flg1) {
407e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
408e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
409e5c89e4eSSatish Balay   }
410e5c89e4eSSatish Balay 
41190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
412acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4137783f70dSSatish Balay   if (!flg1) {
41490d69ab7SBarry Smith     flg1 = PETSC_FALSE;
415acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4167783f70dSSatish Balay   }
417e5c89e4eSSatish Balay   if (flg1) {
418e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
419e5c89e4eSSatish Balay   }
420e5c89e4eSSatish Balay 
421e5c89e4eSSatish Balay   /*
422e5c89e4eSSatish Balay       Set the display variable for graphics
423e5c89e4eSSatish Balay   */
424e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
425e5c89e4eSSatish Balay 
42651dcc849SKerry Stevens 
42751dcc849SKerry Stevens   /*
428e5c89e4eSSatish Balay       Print the PETSc version information
429e5c89e4eSSatish Balay   */
430e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
431e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
432e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
433e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
434e5c89e4eSSatish Balay 
435e5c89e4eSSatish Balay     /*
436e5c89e4eSSatish Balay        Print "higher-level" package version message
437e5c89e4eSSatish Balay     */
438e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
439e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
440e5c89e4eSSatish Balay     }
441e5c89e4eSSatish Balay 
442a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
443e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
444e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
445e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
446e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
447e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
44884e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
449e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
450e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
451e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
452e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
453e5c89e4eSSatish Balay   }
454e5c89e4eSSatish Balay 
455e5c89e4eSSatish Balay   /*
456e5c89e4eSSatish Balay        Print "higher-level" package help message
457e5c89e4eSSatish Balay   */
458e5c89e4eSSatish Balay   if (flg3){
459e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
460e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
461e5c89e4eSSatish Balay     }
462e5c89e4eSSatish Balay   }
463e5c89e4eSSatish Balay 
464e5c89e4eSSatish Balay   /*
465e5c89e4eSSatish Balay       Setup the error handling
466e5c89e4eSSatish Balay   */
46790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
468acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
469cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
47090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
471acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
472cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
47390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
474acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
475e5c89e4eSSatish Balay   if (flg1) {
476e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
477e5c89e4eSSatish Balay   }
47890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
479acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
480cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
48196cc47afSJed Brown   flg1 = PETSC_FALSE;
482acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
48396cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
484e5c89e4eSSatish Balay 
485e5c89e4eSSatish Balay   /*
486e5c89e4eSSatish Balay       Setup debugger information
487e5c89e4eSSatish Balay   */
488e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
489e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
490e5c89e4eSSatish Balay   if (flg1) {
491e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
492e5c89e4eSSatish Balay 
493e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
494e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
495e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
496e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
497e5c89e4eSSatish Balay   }
4985e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
4995e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
500e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
501e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
502e5c89e4eSSatish Balay   if (flg1 || flg2) {
503e5c89e4eSSatish Balay     PetscMPIInt    size;
504e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
505e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
506e5c89e4eSSatish Balay     /*
507e5c89e4eSSatish Balay        we have to make sure that all processors have opened
508e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
509e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
510e5c89e4eSSatish Balay        and kill the program.
511e5c89e4eSSatish Balay     */
512e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
513e5c89e4eSSatish Balay     if (size > 2) {
514533163c2SBarry Smith       PetscMPIInt dummy = 0;
515e5c89e4eSSatish Balay       MPI_Status  status;
516e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
517e5c89e4eSSatish Balay         if (rank != i) {
518e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
519e5c89e4eSSatish Balay         }
520e5c89e4eSSatish Balay       }
521e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
522e5c89e4eSSatish Balay         if (rank != i) {
523e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
524e5c89e4eSSatish Balay         }
525e5c89e4eSSatish Balay       }
526e5c89e4eSSatish Balay     }
527e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
528e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
529e5c89e4eSSatish Balay     lsize = size;
530e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
531e5c89e4eSSatish Balay     if (flag) {
532e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
533e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
534e5c89e4eSSatish Balay       }
535e5c89e4eSSatish Balay     }
536e5c89e4eSSatish Balay     if (!flag) {
537e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
538e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
539e5c89e4eSSatish Balay       if (flg1) {
540e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
541e5c89e4eSSatish Balay       } else {
542e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
543e5c89e4eSSatish Balay       }
544e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
545e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
546e5c89e4eSSatish Balay     }
547e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
548e5c89e4eSSatish Balay   }
549e5c89e4eSSatish Balay 
550e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
551cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
552e5c89e4eSSatish Balay 
55393ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
55422b84c2fSbcordonn   /*
55522b84c2fSbcordonn     Activates new sockets for zope if needed
55622b84c2fSbcordonn   */
55784ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
558d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5596dc8fec2Sbcordonn   if (flgz){
56022b84c2fSbcordonn     int  sockfd;
561f1384234SBarry Smith     char hostname[256];
56222b84c2fSbcordonn     char username[256];
5636dc8fec2Sbcordonn     int  remoteport = 9999;
5649c4c166aSBarry Smith 
56584ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
56684ab5442Sbcordonn     if (!hostname[0]){
5679c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5689c4c166aSBarry Smith     }
56922b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5709c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
57122b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
57222b84c2fSbcordonn     if (flgzout){
57322b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
574606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5756dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5769c4c166aSBarry Smith     } else {
577d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
578d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5799c4c166aSBarry Smith     }
5809c4c166aSBarry Smith   }
58193ba235fSBarry Smith #endif
582ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
583ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
584ffc871a5SBarry Smith   if (flgz){
585ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
586ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
587ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
588ffc871a5SBarry Smith   }
589ffc871a5SBarry Smith #endif
5906dc8fec2Sbcordonn 
591e5c89e4eSSatish Balay   /*
592e5c89e4eSSatish Balay         Setup profiling and logging
593e5c89e4eSSatish Balay   */
5946cf91177SBarry Smith #if defined (PETSC_USE_INFO)
5958bb29257SSatish Balay   {
596e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
5976cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
5988bb29257SSatish Balay     if (flg1 && logname[0]) {
599fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
6008bb29257SSatish Balay     } else if (flg1) {
601fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
602e5c89e4eSSatish Balay     }
603e5c89e4eSSatish Balay   }
604865f6aa8SSatish Balay #endif
605865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
606865f6aa8SSatish Balay   mname[0] = 0;
607f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
608865f6aa8SSatish Balay   if (flg1) {
609865f6aa8SSatish Balay     if (mname[0]) {
610f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
611865f6aa8SSatish Balay     } else {
612f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
613865f6aa8SSatish Balay     }
614865f6aa8SSatish Balay   }
615e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
61690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
617fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
618e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
619e5c89e4eSSatish Balay #endif
62090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
62190d69ab7SBarry Smith   flg2 = PETSC_FALSE;
62290d69ab7SBarry Smith   flg3 = PETSC_FALSE;
623acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
624acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
625d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6269f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
627e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6289f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
629e5c89e4eSSatish Balay 
630e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
631e5c89e4eSSatish Balay   if (flg1) {
632e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
633e5c89e4eSSatish Balay     FILE *file;
634e5c89e4eSSatish Balay     if (mname[0]) {
635e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
636e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
637e5c89e4eSSatish Balay       file = fopen(fname,"w");
638f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
639e5c89e4eSSatish Balay     } else {
640da9f1d6bSBarry Smith       file = PETSC_STDOUT;
641e5c89e4eSSatish Balay     }
642e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
643e5c89e4eSSatish Balay   }
644e5c89e4eSSatish Balay #endif
645e5c89e4eSSatish Balay 
646e5c89e4eSSatish Balay   /*
647e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
648e5c89e4eSSatish Balay   */
64963d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
650e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
651e5c89e4eSSatish Balay #endif
652e5c89e4eSSatish Balay 
653acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
654e5c89e4eSSatish Balay 
655ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
656af359df3SBarry Smith   /*
657af359df3SBarry Smith       Determine whether user specified maximum number of threads
658af359df3SBarry Smith    */
659af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
660af359df3SBarry Smith 
661b154f58aSKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-main",&flg1);CHKERRQ(ierr);
662b154f58aSKerry Stevens   if(flg1) {
663b154f58aSKerry Stevens     cpu_set_t mset;
664b154f58aSKerry Stevens     int icorr,ncorr = get_nprocs();
665b154f58aSKerry Stevens     ierr = PetscOptionsGetInt(PETSC_NULL,"-main",&icorr,PETSC_NULL);CHKERRQ(ierr);
666b154f58aSKerry Stevens     CPU_ZERO(&mset);
667b154f58aSKerry Stevens     CPU_SET(icorr%ncorr,&mset);
668b154f58aSKerry Stevens     sched_setaffinity(0,sizeof(cpu_set_t),&mset);
669b154f58aSKerry Stevens   }
670b154f58aSKerry Stevens 
671af359df3SBarry Smith   PetscInt N_CORES = get_nprocs();
672af359df3SBarry Smith   ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
673af359df3SBarry Smith   char tstr[9];
674af359df3SBarry Smith   char tbuf[2];
675af359df3SBarry Smith   strcpy(tstr,"-thread");
676af359df3SBarry Smith   for(i=0;i<PetscMaxThreads;i++) {
677af359df3SBarry Smith     ThreadCoreAffinity[i] = i;
678af359df3SBarry Smith     sprintf(tbuf,"%d",i);
679af359df3SBarry Smith     strcat(tstr,tbuf);
680af359df3SBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
681af359df3SBarry Smith     if(flg1) {
682af359df3SBarry Smith       ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
683af359df3SBarry Smith       ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
684af359df3SBarry Smith     }
685af359df3SBarry Smith     tstr[7] = '\0';
686af359df3SBarry Smith   }
687cfcfc605SKerry Stevens 
688cfcfc605SKerry Stevens   /*
689cfcfc605SKerry Stevens       Determine whether to use thread pool
690cfcfc605SKerry Stevens    */
691cfcfc605SKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
692cfcfc605SKerry Stevens   if (flg1) {
693cfcfc605SKerry Stevens     PetscUseThreadPool = PETSC_TRUE;
694af359df3SBarry Smith     /* get the thread pool type */
695af359df3SBarry Smith     PetscInt ipool = 0;
696af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
697af359df3SBarry Smith 
698af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
699af359df3SBarry Smith     switch(ipool) {
700af359df3SBarry Smith     case 1:
701af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
702af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
703af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
704af359df3SBarry Smith       MainWait              = &MainWait_Tree;
705af359df3SBarry Smith       MainJob               = &MainJob_Tree;
706af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
707af359df3SBarry Smith       break;
708af359df3SBarry Smith     case 2:
709af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
710af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
711af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
712af359df3SBarry Smith       MainWait              = &MainWait_Main;
713af359df3SBarry Smith       MainJob               = &MainJob_Main;
714af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
715af359df3SBarry Smith       break;
716af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
717af359df3SBarry Smith     case 3:
718af359df3SBarry Smith #else
719af359df3SBarry Smith     default:
720af359df3SBarry Smith #endif
721af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
722af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
723af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
724af359df3SBarry Smith       MainWait              = &MainWait_Chain;
725af359df3SBarry Smith       MainJob               = &MainJob_Chain;
726af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
727af359df3SBarry Smith       break;
728af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
729af359df3SBarry Smith     default:
730af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
731af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
732af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
733af359df3SBarry Smith       MainWait              = &MainWait_True;
734af359df3SBarry Smith       MainJob               = &MainJob_True;
735af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
736af359df3SBarry Smith       break;
737af359df3SBarry Smith #endif
738af359df3SBarry Smith     }
739af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
740683509dcSKerry Stevens   } else {
741683509dcSKerry Stevens     //need to define these in the case on 'no threads' or 'thread create/destroy'
742683509dcSKerry Stevens     //could take any of the above versions
743683509dcSKerry Stevens     MainJob               = &MainJob_Spawn;
744af359df3SBarry Smith   }
745af359df3SBarry Smith #endif
746e5c89e4eSSatish Balay   /*
747e5c89e4eSSatish Balay        Print basic help message
748e5c89e4eSSatish Balay   */
749e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
750e5c89e4eSSatish Balay   if (flg1) {
751e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
752301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
753301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
754301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
755e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
756e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
757e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
758e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
759e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
760e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
761e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
762e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
763e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
764e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
765e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
766e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
767e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
768e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
769e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
770e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
771e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
772e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
773e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7744161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7754161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
776e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
777e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
778e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
779e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
780e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
781e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
782a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
783e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
78440ab9619SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -server <port>: Run PETSc webserver (default port is 8080) see PetscWebServe()\n");CHKERRQ(ierr);
785e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
786e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
7877ef452c0SMatthew G Knepley     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary _summary_python]: logging objects and events\n");CHKERRQ(ierr);
788e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
789e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
790e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
791e5c89e4eSSatish Balay #endif
7926cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
793e5c89e4eSSatish Balay #endif
794e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
795e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
796e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
797e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
798e5c89e4eSSatish Balay   }
799e5c89e4eSSatish Balay 
800a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
801e5c89e4eSSatish Balay   if (flg1) {
802e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
803e5c89e4eSSatish Balay   }
804e5c89e4eSSatish Balay 
8056cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
806e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
807e5c89e4eSSatish Balay   if (f) {
8086cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
809e5c89e4eSSatish Balay   }
810827f890bSBarry Smith 
8118154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
812c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
81373113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
81473113deaSBarry Smith   else flg1 = PETSC_FALSE;
8158154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
8168154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
817bab1f7e6SVictor Minden #endif
818192daf7cSBarry Smith 
819e5c89e4eSSatish Balay   PetscFunctionReturn(0);
820e5c89e4eSSatish Balay }
821df413903SBarry Smith 
822ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
823ba61063dSBarry Smith 
82451d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
82551d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
82651d315f7SKerry Stevens   PetscErrorCode iterr;
82751d315f7SKerry Stevens   int icorr,ierr;
82851d315f7SKerry Stevens   int* pId = (int*)arg;
82951d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
83051d315f7SKerry Stevens   PetscBool PeeOn;
83151d315f7SKerry Stevens   cpu_set_t mset;
8329e800a48SKerry Stevens   //printf("Thread %d In Tree Thread Function\n",ThreadId);
83351d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
83451d315f7SKerry Stevens   CPU_ZERO(&mset);
83551d315f7SKerry Stevens   CPU_SET(icorr,&mset);
83651d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
83751d315f7SKerry Stevens 
83851d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
83951d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
84051d315f7SKerry Stevens   }
84151d315f7SKerry Stevens   else {
84251d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
84351d315f7SKerry Stevens   }
84451d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
845ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
84651d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
84751d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
84851d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
84951d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
85051d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
851ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
852ba61063dSBarry Smith            upon return, has the lock */
85351d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
85451d315f7SKerry Stevens         }
85551d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
85651d315f7SKerry Stevens       }
85751d315f7SKerry Stevens     }
858ba61063dSBarry Smith     /* your subordinates are now ready */
85951d315f7SKerry Stevens   }
86051d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
861ba61063dSBarry Smith   /* update your ready status */
86251d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
86351d315f7SKerry Stevens   if(ThreadId==0) {
86451d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
865ba61063dSBarry Smith     /* ignal main */
86651d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
86751d315f7SKerry Stevens   }
86851d315f7SKerry Stevens   else {
869ba61063dSBarry Smith     /* tell your boss that you're ready to work */
87051d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
87151d315f7SKerry Stevens   }
872ba61063dSBarry Smith   /* the while loop needs to have an exit
873ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
874ba61063dSBarry Smith    and calling FuncFinish */
87551d315f7SKerry Stevens   while(PetscThreadGo) {
876ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
877ba61063dSBarry Smith       waiting when you don't have to causes problems
878ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
87951d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
880ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
881ba61063dSBarry Smith        upon return, has the lock */
88251d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
88351d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
88451d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
88551d315f7SKerry Stevens     }
88651d315f7SKerry Stevens     if(ThreadId==0) {
88751d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
88851d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
88951d315f7SKerry Stevens     }
89051d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
89151d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
892ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
89351d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
89451d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
89551d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
89651d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
89751d315f7SKerry Stevens         }
89851d315f7SKerry Stevens       }
89951d315f7SKerry Stevens     }
900ba61063dSBarry Smith     /* do your job */
90151d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
90251d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
90351d315f7SKerry Stevens     }
90451d315f7SKerry Stevens     else {
90551d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
90651d315f7SKerry Stevens     }
90751d315f7SKerry Stevens     if(iterr!=0) {
90851d315f7SKerry Stevens       ithreaderr = 1;
90951d315f7SKerry Stevens     }
91051d315f7SKerry Stevens     if(PetscThreadGo) {
911ba61063dSBarry Smith       /* reset job, get ready for more */
91251d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
913ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
914ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
91551d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
91651d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
91751d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
91851d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
91951d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
920ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
921ba61063dSBarry Smith                upon return, has the lock */
92251d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
92351d315f7SKerry Stevens             }
92451d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
92551d315f7SKerry Stevens           }
92651d315f7SKerry Stevens 	}
927ba61063dSBarry Smith         /* your subordinates are now ready */
92851d315f7SKerry Stevens       }
92951d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
93051d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
93151d315f7SKerry Stevens       if(ThreadId==0) {
932ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
933ba61063dSBarry Smith         /* root thread signals 'main' */
93451d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
93551d315f7SKerry Stevens       }
93651d315f7SKerry Stevens       else {
937ba61063dSBarry Smith         /* signal your boss before you go to sleep */
93851d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
93951d315f7SKerry Stevens       }
94051d315f7SKerry Stevens     }
94151d315f7SKerry Stevens   }
94251d315f7SKerry Stevens   return NULL;
94351d315f7SKerry Stevens }
94451d315f7SKerry Stevens 
94551d315f7SKerry Stevens #undef __FUNCT__
94651d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
94751d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
94851d315f7SKerry Stevens   PetscInt i,ierr;
94951d315f7SKerry Stevens   int status;
95051d315f7SKerry Stevens 
95151d315f7SKerry Stevens   if(PetscUseThreadPool) {
95251d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
95351d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
95451d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
95551d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
95651d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
95751d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
95851d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
95951d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
96051d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
96151d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
96251d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
96351d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
964ba61063dSBarry Smith     /* initialize job structure */
96551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
96651d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
96751d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
96851d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
96951d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
97051d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
97151d315f7SKerry Stevens     }
97251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
97351d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
97451d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
97551d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
97651d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
97751d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
97851d315f7SKerry Stevens     }
97951d315f7SKerry Stevens     job_tree.pfunc = NULL;
98051d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
98151d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
98251d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
98351d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
984ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
98551d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
986ba61063dSBarry Smith     /* create threads */
98751d315f7SKerry Stevens     for(i=0; i<N; i++) {
98851d315f7SKerry Stevens       pVal[i] = i;
98951d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
990ba61063dSBarry Smith       /* should check status */
99151d315f7SKerry Stevens     }
99251d315f7SKerry Stevens   }
99351d315f7SKerry Stevens   return NULL;
99451d315f7SKerry Stevens }
99551d315f7SKerry Stevens 
99651d315f7SKerry Stevens #undef __FUNCT__
99751d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
99851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
99951d315f7SKerry Stevens   int i,ierr;
100051d315f7SKerry Stevens   void* jstatus;
100151d315f7SKerry Stevens 
100251d315f7SKerry Stevens   PetscFunctionBegin;
100351d315f7SKerry Stevens 
1004ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1005ba61063dSBarry Smith   /* join the threads */
100651d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
100751d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1008ba61063dSBarry Smith     /* do error checking*/
100951d315f7SKerry Stevens   }
101051d315f7SKerry Stevens   free(PetscThreadPoint);
101151d315f7SKerry Stevens   free(arrmutex);
101251d315f7SKerry Stevens   free(arrcond1);
101351d315f7SKerry Stevens   free(arrcond2);
101451d315f7SKerry Stevens   free(arrstart);
101551d315f7SKerry Stevens   free(arrready);
101651d315f7SKerry Stevens   free(job_tree.pdata);
101751d315f7SKerry Stevens   free(pVal);
1018cfcfc605SKerry Stevens 
101951d315f7SKerry Stevens   PetscFunctionReturn(0);
102051d315f7SKerry Stevens }
102151d315f7SKerry Stevens 
102251d315f7SKerry Stevens #undef __FUNCT__
102351d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
102451d315f7SKerry Stevens void MainWait_Tree() {
102551d315f7SKerry Stevens   int ierr;
102651d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
102751d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
102851d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
102951d315f7SKerry Stevens   }
103051d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
103151d315f7SKerry Stevens }
103251d315f7SKerry Stevens 
103351d315f7SKerry Stevens #undef __FUNCT__
103451d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
103551d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
103651d315f7SKerry Stevens   int i,ierr;
103751d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
103836d20dc5SKerry Stevens 
103951d315f7SKerry Stevens   MainWait();
104051d315f7SKerry Stevens   job_tree.pfunc = pFunc;
104151d315f7SKerry Stevens   job_tree.pdata = data;
104251d315f7SKerry Stevens   job_tree.startJob = PETSC_TRUE;
104351d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
104451d315f7SKerry Stevens     *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
104551d315f7SKerry Stevens   }
104651d315f7SKerry Stevens   job_tree.eJobStat = JobInitiated;
104751d315f7SKerry Stevens   ierr = pthread_cond_signal(job_tree.cond2array[0]);
104851d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1049ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
105051d315f7SKerry Stevens   }
105136d20dc5SKerry Stevens 
105251d315f7SKerry Stevens   if(ithreaderr) {
105351d315f7SKerry Stevens     ijoberr = ithreaderr;
105451d315f7SKerry Stevens   }
105551d315f7SKerry Stevens   return ijoberr;
105651d315f7SKerry Stevens }
105751d315f7SKerry Stevens /****  ****/
105851d315f7SKerry Stevens 
105951d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
106051d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
106151d315f7SKerry Stevens   PetscErrorCode iterr;
106251d315f7SKerry Stevens   int icorr,ierr;
106351d315f7SKerry Stevens   int* pId = (int*)arg;
106451d315f7SKerry Stevens   int ThreadId = *pId;
106551d315f7SKerry Stevens   cpu_set_t mset;
10669e800a48SKerry Stevens   //printf("Thread %d In Main Thread Function\n",ThreadId);
106751d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
106851d315f7SKerry Stevens   CPU_ZERO(&mset);
106951d315f7SKerry Stevens   CPU_SET(icorr,&mset);
107051d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
107151d315f7SKerry Stevens 
107251d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1073ba61063dSBarry Smith   /* update your ready status */
107451d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1075ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
107651d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
107751d315f7SKerry Stevens 
1078ba61063dSBarry Smith   /* the while loop needs to have an exit
1079ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1080ba61063dSBarry Smith      and calling FuncFinish */
108151d315f7SKerry Stevens   while(PetscThreadGo) {
1082ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1083ba61063dSBarry Smith        waiting when you don't have to causes problems
1084ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
108551d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1086ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1087ba61063dSBarry Smith        upon return, has the lock */
108851d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1089ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
109051d315f7SKerry Stevens     }
109151d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
109251d315f7SKerry Stevens     if(job_main.pdata==NULL) {
109351d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
109451d315f7SKerry Stevens     }
109551d315f7SKerry Stevens     else {
109651d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
109751d315f7SKerry Stevens     }
109851d315f7SKerry Stevens     if(iterr!=0) {
109951d315f7SKerry Stevens       ithreaderr = 1;
110051d315f7SKerry Stevens     }
110151d315f7SKerry Stevens     if(PetscThreadGo) {
1102ba61063dSBarry Smith       /* reset job, get ready for more */
110351d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
110451d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1105ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
110651d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
110751d315f7SKerry Stevens     }
110851d315f7SKerry Stevens   }
110951d315f7SKerry Stevens   return NULL;
111051d315f7SKerry Stevens }
111151d315f7SKerry Stevens 
111251d315f7SKerry Stevens #undef __FUNCT__
111351d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
111451d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
111551d315f7SKerry Stevens   PetscInt i,ierr;
111651d315f7SKerry Stevens   int status;
111751d315f7SKerry Stevens 
111851d315f7SKerry Stevens   if(PetscUseThreadPool) {
111951d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
112051d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
112151d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
112251d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
112351d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
112451d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
112551d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
112651d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
112751d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112851d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112951d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1130ba61063dSBarry Smith     /* initialize job structure */
113151d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
113251d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
113351d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
113451d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
113551d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
113651d315f7SKerry Stevens     }
113751d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
113851d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
113951d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
114051d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
114151d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
114251d315f7SKerry Stevens     }
114351d315f7SKerry Stevens     job_main.pfunc = NULL;
114451d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
114551d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1146ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
114751d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1148ba61063dSBarry Smith     /* create threads */
114951d315f7SKerry Stevens     for(i=0; i<N; i++) {
115051d315f7SKerry Stevens       pVal[i] = i;
115151d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1152ba61063dSBarry Smith       /* error check */
115351d315f7SKerry Stevens     }
115451d315f7SKerry Stevens   }
115551d315f7SKerry Stevens   else {
115651d315f7SKerry Stevens   }
115751d315f7SKerry Stevens   return NULL;
115851d315f7SKerry Stevens }
115951d315f7SKerry Stevens 
116051d315f7SKerry Stevens #undef __FUNCT__
116151d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
116251d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
116351d315f7SKerry Stevens   int i,ierr;
116451d315f7SKerry Stevens   void* jstatus;
116551d315f7SKerry Stevens 
116651d315f7SKerry Stevens   PetscFunctionBegin;
116751d315f7SKerry Stevens 
1168ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1169ba61063dSBarry Smith   /* join the threads */
117051d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1171ba61063dSBarry Smith     ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
117251d315f7SKerry Stevens   }
117351d315f7SKerry Stevens   free(PetscThreadPoint);
117451d315f7SKerry Stevens   free(arrmutex);
117551d315f7SKerry Stevens   free(arrcond1);
117651d315f7SKerry Stevens   free(arrcond2);
117751d315f7SKerry Stevens   free(arrstart);
117851d315f7SKerry Stevens   free(arrready);
117951d315f7SKerry Stevens   free(job_main.pdata);
118051d315f7SKerry Stevens   free(pVal);
1181cfcfc605SKerry Stevens 
118251d315f7SKerry Stevens   PetscFunctionReturn(0);
118351d315f7SKerry Stevens }
118451d315f7SKerry Stevens 
118551d315f7SKerry Stevens #undef __FUNCT__
118651d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
118751d315f7SKerry Stevens void MainWait_Main() {
118851d315f7SKerry Stevens   int i,ierr;
118951d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
119051d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
119151d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
119251d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
119351d315f7SKerry Stevens     }
119451d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
119551d315f7SKerry Stevens   }
119651d315f7SKerry Stevens }
119751d315f7SKerry Stevens 
119851d315f7SKerry Stevens #undef __FUNCT__
119951d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
120051d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
120151d315f7SKerry Stevens   int i,ierr;
120251d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
120336d20dc5SKerry Stevens 
1204ba61063dSBarry Smith   MainWait(); /* you know everyone is waiting to be signalled! */
120551d315f7SKerry Stevens   job_main.pfunc = pFunc;
120651d315f7SKerry Stevens   job_main.pdata = data;
120751d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1208ba61063dSBarry Smith     *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
120951d315f7SKerry Stevens   }
1210ba61063dSBarry Smith   /* tell the threads to go to work */
121151d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
121251d315f7SKerry Stevens     ierr = pthread_cond_signal(job_main.cond2array[i]);
121351d315f7SKerry Stevens   }
121451d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1215ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
121651d315f7SKerry Stevens   }
121736d20dc5SKerry Stevens 
121851d315f7SKerry Stevens   if(ithreaderr) {
121951d315f7SKerry Stevens     ijoberr = ithreaderr;
122051d315f7SKerry Stevens   }
122151d315f7SKerry Stevens   return ijoberr;
122251d315f7SKerry Stevens }
122351d315f7SKerry Stevens /****  ****/
122451d315f7SKerry Stevens 
122551d315f7SKerry Stevens /**** Chain Thread Functions ****/
122651d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
122751d315f7SKerry Stevens   PetscErrorCode iterr;
122851d315f7SKerry Stevens   int icorr,ierr;
122951d315f7SKerry Stevens   int* pId = (int*)arg;
123051d315f7SKerry Stevens   int ThreadId = *pId;
123151d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
123251d315f7SKerry Stevens   PetscBool PeeOn;
123351d315f7SKerry Stevens   cpu_set_t mset;
12349e800a48SKerry Stevens   //printf("Thread %d In Chain Thread Function\n",ThreadId);
123551d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
123651d315f7SKerry Stevens   CPU_ZERO(&mset);
123751d315f7SKerry Stevens   CPU_SET(icorr,&mset);
123851d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
123951d315f7SKerry Stevens 
124051d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
124151d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
124251d315f7SKerry Stevens   }
124351d315f7SKerry Stevens   else {
124451d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
124551d315f7SKerry Stevens   }
124651d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1247ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
124851d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
124951d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1250ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1251ba61063dSBarry Smith        upon return, has the lock */
125251d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
125351d315f7SKerry Stevens     }
125451d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1255ba61063dSBarry Smith     /* your subordinate is now ready*/
125651d315f7SKerry Stevens   }
125751d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1258ba61063dSBarry Smith   /* update your ready status */
125951d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
126051d315f7SKerry Stevens   if(ThreadId==0) {
126151d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1262ba61063dSBarry Smith     /* signal main */
126351d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
126451d315f7SKerry Stevens   }
126551d315f7SKerry Stevens   else {
1266ba61063dSBarry Smith     /* tell your boss that you're ready to work */
126751d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
126851d315f7SKerry Stevens   }
1269ba61063dSBarry Smith   /*  the while loop needs to have an exit
1270ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1271ba61063dSBarry Smith    and calling FuncFinish */
127251d315f7SKerry Stevens   while(PetscThreadGo) {
1273ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1274ba61063dSBarry Smith        waiting when you don't have to causes problems
1275ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
127651d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1277ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1278ba61063dSBarry Smith        upon return, has the lock */
127951d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
128051d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
128151d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
128251d315f7SKerry Stevens     }
128351d315f7SKerry Stevens     if(ThreadId==0) {
128451d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
128551d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
128651d315f7SKerry Stevens     }
128751d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
128851d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1289ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
129051d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
129151d315f7SKerry Stevens     }
1292ba61063dSBarry Smith     /* do your job */
129351d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
129451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
129551d315f7SKerry Stevens     }
129651d315f7SKerry Stevens     else {
129751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
129851d315f7SKerry Stevens     }
129951d315f7SKerry Stevens     if(iterr!=0) {
130051d315f7SKerry Stevens       ithreaderr = 1;
130151d315f7SKerry Stevens     }
130251d315f7SKerry Stevens     if(PetscThreadGo) {
1303ba61063dSBarry Smith       /* reset job, get ready for more */
130451d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1305ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1306ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
130751d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
130851d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1309ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1310ba61063dSBarry Smith            upon return, has the lock */
131151d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
131251d315f7SKerry Stevens         }
131351d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1314ba61063dSBarry Smith         /* your subordinate is now ready */
131551d315f7SKerry Stevens       }
131651d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
131751d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
131851d315f7SKerry Stevens       if(ThreadId==0) {
1319ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1320ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
132151d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
132251d315f7SKerry Stevens       }
132351d315f7SKerry Stevens       else {
1324ba61063dSBarry Smith         /* signal your boss before you go to sleep */
132551d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
132651d315f7SKerry Stevens       }
132751d315f7SKerry Stevens     }
132851d315f7SKerry Stevens   }
132951d315f7SKerry Stevens   return NULL;
133051d315f7SKerry Stevens }
133151d315f7SKerry Stevens 
133251d315f7SKerry Stevens #undef __FUNCT__
133351d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
133451d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
133551d315f7SKerry Stevens   PetscInt i,ierr;
133651d315f7SKerry Stevens   int status;
133751d315f7SKerry Stevens 
133851d315f7SKerry Stevens   if(PetscUseThreadPool) {
133951d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
134051d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
134151d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
134251d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
134351d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
134451d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
134551d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
134651d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
134751d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134851d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134951d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
135051d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1351ba61063dSBarry Smith     /* initialize job structure */
135251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
135351d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
135451d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
135551d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
135651d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
135751d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
135851d315f7SKerry Stevens     }
135951d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
136051d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
136151d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
136251d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
136351d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
136451d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
136551d315f7SKerry Stevens     }
136651d315f7SKerry Stevens     job_chain.pfunc = NULL;
136751d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
136851d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
136951d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
137051d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1371ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
137251d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1373ba61063dSBarry Smith     /* create threads */
137451d315f7SKerry Stevens     for(i=0; i<N; i++) {
137551d315f7SKerry Stevens       pVal[i] = i;
137651d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1377ba61063dSBarry Smith       /* should check error */
137851d315f7SKerry Stevens     }
137951d315f7SKerry Stevens   }
138051d315f7SKerry Stevens   else {
138151d315f7SKerry Stevens   }
138251d315f7SKerry Stevens   return NULL;
138351d315f7SKerry Stevens }
138451d315f7SKerry Stevens 
138551d315f7SKerry Stevens 
138651d315f7SKerry Stevens #undef __FUNCT__
138751d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
138851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
138951d315f7SKerry Stevens   int i,ierr;
139051d315f7SKerry Stevens   void* jstatus;
139151d315f7SKerry Stevens 
139251d315f7SKerry Stevens   PetscFunctionBegin;
139351d315f7SKerry Stevens 
1394ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1395ba61063dSBarry Smith   /* join the threads */
139651d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
139751d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1398ba61063dSBarry Smith     /* should check error */
139951d315f7SKerry Stevens   }
140051d315f7SKerry Stevens   free(PetscThreadPoint);
140151d315f7SKerry Stevens   free(arrmutex);
140251d315f7SKerry Stevens   free(arrcond1);
140351d315f7SKerry Stevens   free(arrcond2);
140451d315f7SKerry Stevens   free(arrstart);
140551d315f7SKerry Stevens   free(arrready);
140651d315f7SKerry Stevens   free(job_chain.pdata);
140751d315f7SKerry Stevens   free(pVal);
1408cfcfc605SKerry Stevens 
140951d315f7SKerry Stevens   PetscFunctionReturn(0);
141051d315f7SKerry Stevens }
141151d315f7SKerry Stevens 
141251d315f7SKerry Stevens #undef __FUNCT__
141351d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
141451d315f7SKerry Stevens void MainWait_Chain() {
141551d315f7SKerry Stevens   int ierr;
141651d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
141751d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
141851d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
141951d315f7SKerry Stevens   }
142051d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
142151d315f7SKerry Stevens }
142251d315f7SKerry Stevens 
142351d315f7SKerry Stevens #undef __FUNCT__
142451d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
142551d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
142651d315f7SKerry Stevens   int i,ierr;
142751d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
142836d20dc5SKerry Stevens 
142951d315f7SKerry Stevens   MainWait();
143051d315f7SKerry Stevens   job_chain.pfunc = pFunc;
143151d315f7SKerry Stevens   job_chain.pdata = data;
143251d315f7SKerry Stevens   job_chain.startJob = PETSC_TRUE;
143351d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
143451d315f7SKerry Stevens     *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
143551d315f7SKerry Stevens   }
143651d315f7SKerry Stevens   job_chain.eJobStat = JobInitiated;
143751d315f7SKerry Stevens   ierr = pthread_cond_signal(job_chain.cond2array[0]);
143851d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1439ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
144051d315f7SKerry Stevens   }
144136d20dc5SKerry Stevens 
144251d315f7SKerry Stevens   if(ithreaderr) {
144351d315f7SKerry Stevens     ijoberr = ithreaderr;
144451d315f7SKerry Stevens   }
144551d315f7SKerry Stevens   return ijoberr;
144651d315f7SKerry Stevens }
144751d315f7SKerry Stevens /****  ****/
144851d315f7SKerry Stevens 
1449ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
145051d315f7SKerry Stevens /**** True Thread Functions ****/
145151d315f7SKerry Stevens void* PetscThreadFunc_True(void* arg) {
145251d315f7SKerry Stevens   int icorr,ierr,iVal;
145351dcc849SKerry Stevens   int* pId = (int*)arg;
145451dcc849SKerry Stevens   int ThreadId = *pId;
14550ca81413SKerry Stevens   PetscErrorCode iterr;
145651d315f7SKerry Stevens   cpu_set_t mset;
14579e800a48SKerry Stevens   //printf("Thread %d In True Pool Thread Function\n",ThreadId);
145851d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
145951d315f7SKerry Stevens   CPU_ZERO(&mset);
146051d315f7SKerry Stevens   CPU_SET(icorr,&mset);
146151d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
146251d315f7SKerry Stevens 
146351d315f7SKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
146451d315f7SKerry Stevens   job_true.iNumReadyThreads++;
146551d315f7SKerry Stevens   if(job_true.iNumReadyThreads==PetscMaxThreads) {
146651dcc849SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
146751dcc849SKerry Stevens   }
1468ba61063dSBarry Smith   /*the while loop needs to have an exit
1469ba61063dSBarry Smith     the 'main' thread can terminate all the threads by performing a broadcast
1470ba61063dSBarry Smith    and calling FuncFinish */
147151dcc849SKerry Stevens   while(PetscThreadGo) {
1472ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
1473ba61063dSBarry Smith       waiting when you don't have to causes problems
1474ba61063dSBarry Smith      also need to wait if another thread sneaks in and messes with the predicate */
147551d315f7SKerry Stevens     while(job_true.startJob==PETSC_FALSE&&job_true.iNumJobThreads==0) {
1476ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1477ba61063dSBarry Smith        upon return, has the lock */
14786c9b208dSKerry Stevens       //printf("Thread %d Going to Sleep!\n",ThreadId);
147951d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
148051dcc849SKerry Stevens     }
148151d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
148251d315f7SKerry Stevens     job_true.iNumJobThreads--;
148351d315f7SKerry Stevens     job_true.iNumReadyThreads--;
148451d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
148551d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
148651d315f7SKerry Stevens     if(job_true.pdata==NULL) {
148751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
148851dcc849SKerry Stevens     }
148951dcc849SKerry Stevens     else {
149051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
149151dcc849SKerry Stevens     }
14920ca81413SKerry Stevens     if(iterr!=0) {
14930ca81413SKerry Stevens       ithreaderr = 1;
14940ca81413SKerry Stevens     }
14956c9b208dSKerry Stevens     //printf("Thread %d Finished Job\n",ThreadId);
1496ba61063dSBarry Smith     /* the barrier is necessary BECAUSE: look at job_true.iNumReadyThreads
1497ba61063dSBarry Smith       what happens if a thread finishes before they all start? BAD!
1498ba61063dSBarry Smith      what happens if a thread finishes before any else start? BAD! */
1499ba61063dSBarry Smith     pthread_barrier_wait(job_true.pbarr); /* ensures all threads are finished */
1500ba61063dSBarry Smith     /* reset job */
150151dcc849SKerry Stevens     if(PetscThreadGo) {
150251d315f7SKerry Stevens       pthread_mutex_lock(&job_true.mutex);
150351d315f7SKerry Stevens       job_true.iNumReadyThreads++;
150451d315f7SKerry Stevens       if(job_true.iNumReadyThreads==PetscMaxThreads) {
1505ba61063dSBarry Smith 	/* signal the 'main' thread that the job is done! (only done once) */
150651dcc849SKerry Stevens 	ierr = pthread_cond_signal(&main_cond);
150751dcc849SKerry Stevens       }
150851dcc849SKerry Stevens     }
150951dcc849SKerry Stevens   }
151051dcc849SKerry Stevens   return NULL;
151151dcc849SKerry Stevens }
151251dcc849SKerry Stevens 
1513f09cb4aaSKerry Stevens #undef __FUNCT__
151451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_True"
151551d315f7SKerry Stevens void* PetscThreadInitialize_True(PetscInt N) {
151651dcc849SKerry Stevens   PetscInt i;
151751dcc849SKerry Stevens   int status;
15180ca81413SKerry Stevens 
1519f09cb4aaSKerry Stevens   pVal = (int*)malloc(N*sizeof(int));
1520ba61063dSBarry Smith   /* allocate memory in the heap for the thread structure */
152151dcc849SKerry Stevens   PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1522ba61063dSBarry Smith   BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
152351d315f7SKerry Stevens   job_true.pdata = (void**)malloc(N*sizeof(void*));
152451dcc849SKerry Stevens   for(i=0; i<N; i++) {
1525f09cb4aaSKerry Stevens     pVal[i] = i;
1526f09cb4aaSKerry Stevens     status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1527ba61063dSBarry Smith     /* error check to ensure proper thread creation */
152851dcc849SKerry Stevens     status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1529ba61063dSBarry Smith     /* should check error */
153051dcc849SKerry Stevens   }
15316c9b208dSKerry Stevens   //printf("Finished True Thread Pool Initialization\n");
153251dcc849SKerry Stevens   return NULL;
153351dcc849SKerry Stevens }
153451dcc849SKerry Stevens 
1535f09cb4aaSKerry Stevens 
1536f09cb4aaSKerry Stevens #undef __FUNCT__
153751d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
153851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
153951dcc849SKerry Stevens   int i,ierr;
154051dcc849SKerry Stevens   void* jstatus;
154151dcc849SKerry Stevens 
154251dcc849SKerry Stevens   PetscFunctionBegin;
1543cfcfc605SKerry Stevens 
1544ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1545ba61063dSBarry Smith   /* join the threads */
154651dcc849SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
154751dcc849SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
154851dcc849SKerry Stevens   }
154951dcc849SKerry Stevens   free(BarrPoint);
155051dcc849SKerry Stevens   free(PetscThreadPoint);
1551cfcfc605SKerry Stevens 
155251dcc849SKerry Stevens   PetscFunctionReturn(0);
155351dcc849SKerry Stevens }
155451dcc849SKerry Stevens 
1555f09cb4aaSKerry Stevens #undef __FUNCT__
155651d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
155751d315f7SKerry Stevens void MainWait_True() {
155851dcc849SKerry Stevens   int ierr;
15596c9b208dSKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
156051d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
156151d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
156251dcc849SKerry Stevens   }
156351d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
156451dcc849SKerry Stevens }
156551dcc849SKerry Stevens 
1566f09cb4aaSKerry Stevens #undef __FUNCT__
156751d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
156851d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
156951dcc849SKerry Stevens   int ierr;
15700ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
157136d20dc5SKerry Stevens 
15720ca81413SKerry Stevens   MainWait();
157351d315f7SKerry Stevens   job_true.pfunc = pFunc;
157451d315f7SKerry Stevens   job_true.pdata = data;
157551d315f7SKerry Stevens   job_true.pbarr = &BarrPoint[n];
157651d315f7SKerry Stevens   job_true.iNumJobThreads = n;
157751d315f7SKerry Stevens   job_true.startJob = PETSC_TRUE;
157851d315f7SKerry Stevens   ierr = pthread_cond_broadcast(&job_true.cond);
15790ca81413SKerry Stevens   if(pFunc!=FuncFinish) {
1580ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done */
15810ca81413SKerry Stevens   }
158236d20dc5SKerry Stevens 
15830ca81413SKerry Stevens   if(ithreaderr) {
15840ca81413SKerry Stevens     ijoberr = ithreaderr;
15850ca81413SKerry Stevens   }
15860ca81413SKerry Stevens   return ijoberr;
158751dcc849SKerry Stevens }
1588683509dcSKerry Stevens /**** NO THREAD POOL FUNCTION ****/
1589683509dcSKerry Stevens #undef __FUNCT__
1590683509dcSKerry Stevens #define __FUNCT__ "MainJob_Spawn"
1591683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void** data,PetscInt n) {
1592683509dcSKerry Stevens   PetscErrorCode ijoberr = 0;
1593683509dcSKerry Stevens 
1594683509dcSKerry Stevens   pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
1595cfcfc605SKerry Stevens   PetscThreadPoint = apThread; /* point to same place */
1596683509dcSKerry Stevens   PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1597683509dcSKerry Stevens   PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
1598683509dcSKerry Stevens   free(apThread);
1599683509dcSKerry Stevens 
1600683509dcSKerry Stevens   return ijoberr;
1601683509dcSKerry Stevens }
160251d315f7SKerry Stevens /****  ****/
1603ba61063dSBarry Smith #endif
160451dcc849SKerry Stevens 
160551dcc849SKerry Stevens void* FuncFinish(void* arg) {
160651dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16070ca81413SKerry Stevens   return(0);
160851dcc849SKerry Stevens }
1609ba61063dSBarry Smith 
1610ba61063dSBarry Smith #endif
1611