xref: /petsc/src/sys/objects/init.c (revision ef386f4b9e9a5f76b306de1e9a9e07c881c3254b)
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 */
8*ef386f4bSSatish Balay 
9*ef386f4bSSatish Balay #include <petscsys.h>        /*I  "petscsys.h"   I*/
10*ef386f4bSSatish Balay 
11*ef386f4bSSatish Balay #if defined(PETSC_USE_PTHREAD)
128f54b378SBarry Smith #ifndef _GNU_SOURCE
138f54b378SBarry Smith #define _GNU_SOURCE
148f54b378SBarry Smith #endif
15*ef386f4bSSatish Balay #if defined(PETSC_HAVE_SCHED_H)
168f54b378SBarry Smith #include <sched.h>
17*ef386f4bSSatish Balay #endif
1851dcc849SKerry Stevens #include <pthread.h>
19ba61063dSBarry Smith #endif
20*ef386f4bSSatish Balay 
21ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
2251d315f7SKerry Stevens #include <sys/sysinfo.h>
23ba61063dSBarry Smith #endif
2451d315f7SKerry Stevens #include <unistd.h>
25e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
26e5c89e4eSSatish Balay #include <stdlib.h>
27e5c89e4eSSatish Balay #endif
28e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
29e5c89e4eSSatish Balay #include <malloc.h>
30e5c89e4eSSatish Balay #endif
31555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
32555d055bSBarry Smith #include <valgrind/valgrind.h>
33555d055bSBarry Smith #endif
34555d055bSBarry Smith 
35e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
36e5c89e4eSSatish Balay /*
37e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
38e5c89e4eSSatish Balay    already started before PETSc was initialized.
39e5c89e4eSSatish Balay */
407087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
417087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
427087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
4351dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
4451dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
457087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
467087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
47ba61063dSBarry Smith 
48ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
4951dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
5051dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
51af359df3SBarry Smith #define PETSC_HAVE_PTHREAD_BARRIER
52ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
53ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
54ba61063dSBarry Smith #endif
5551d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
56f09cb4aaSKerry Stevens int*         pVal;
5751dcc849SKerry Stevens 
58ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
5951d315f7SKerry Stevens int* ThreadCoreAffinity;
6051d315f7SKerry Stevens 
61ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
6251d315f7SKerry Stevens 
6351d315f7SKerry Stevens typedef struct {
6451d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
6551d315f7SKerry Stevens   pthread_cond_t**  cond1array;
6651d315f7SKerry Stevens   pthread_cond_t** cond2array;
6751d315f7SKerry Stevens   void* (*pfunc)(void*);
6851d315f7SKerry Stevens   void** pdata;
6951d315f7SKerry Stevens   PetscBool startJob;
7051d315f7SKerry Stevens   estat eJobStat;
7151d315f7SKerry Stevens   PetscBool** arrThreadStarted;
7251d315f7SKerry Stevens   PetscBool** arrThreadReady;
7351d315f7SKerry Stevens } sjob_tree;
7451d315f7SKerry Stevens sjob_tree job_tree;
7551d315f7SKerry Stevens typedef struct {
7651d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
7751d315f7SKerry Stevens   pthread_cond_t**  cond1array;
7851d315f7SKerry Stevens   pthread_cond_t** cond2array;
7951d315f7SKerry Stevens   void* (*pfunc)(void*);
8051d315f7SKerry Stevens   void** pdata;
8151d315f7SKerry Stevens   PetscBool** arrThreadReady;
8251d315f7SKerry Stevens } sjob_main;
8351d315f7SKerry Stevens sjob_main job_main;
8451d315f7SKerry Stevens typedef struct {
8551d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8651d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8751d315f7SKerry Stevens   pthread_cond_t** cond2array;
8851d315f7SKerry Stevens   void* (*pfunc)(void*);
8951d315f7SKerry Stevens   void** pdata;
9051d315f7SKerry Stevens   PetscBool startJob;
9151d315f7SKerry Stevens   estat eJobStat;
9251d315f7SKerry Stevens   PetscBool** arrThreadStarted;
9351d315f7SKerry Stevens   PetscBool** arrThreadReady;
9451d315f7SKerry Stevens } sjob_chain;
9551d315f7SKerry Stevens sjob_chain job_chain;
96ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
9751dcc849SKerry Stevens typedef struct {
9851dcc849SKerry Stevens   pthread_mutex_t mutex;
9951dcc849SKerry Stevens   pthread_cond_t cond;
10051dcc849SKerry Stevens   void* (*pfunc)(void*);
10151dcc849SKerry Stevens   void** pdata;
10251dcc849SKerry Stevens   pthread_barrier_t* pbarr;
10351dcc849SKerry Stevens   int iNumJobThreads;
10451dcc849SKerry Stevens   int iNumReadyThreads;
10551dcc849SKerry Stevens   PetscBool startJob;
10651d315f7SKerry Stevens } sjob_true;
10751d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
108ba61063dSBarry Smith #endif
10951dcc849SKerry Stevens 
110ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
111ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
112ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
113ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
114ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
115ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
11651dcc849SKerry Stevens 
11751d315f7SKerry Stevens /* Function Pointers */
11851d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
11951d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
12051d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
12151d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
12251d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
12336d20dc5SKerry Stevens /**** Tree Thread Pool Functions ****/
12451d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
12551d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
12651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
12751d315f7SKerry Stevens void           MainWait_Tree(void);
12851d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
12936d20dc5SKerry Stevens /**** Main Thread Pool Functions ****/
13051d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
13151d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
13251d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
13351d315f7SKerry Stevens void           MainWait_Main(void);
13451d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
13536d20dc5SKerry Stevens /**** Chain Thread Pool Functions ****/
13651d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
13751d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
13851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
13951d315f7SKerry Stevens void           MainWait_Chain(void);
14051d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
14136d20dc5SKerry Stevens /**** True Thread Pool Functions ****/
14251d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
14351d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
14451d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
14551d315f7SKerry Stevens void           MainWait_True(void);
14651d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
14736d20dc5SKerry Stevens /**** NO Thread Pool Function  ****/
148683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void**,PetscInt);
14936d20dc5SKerry Stevens /****  ****/
15051dcc849SKerry Stevens void* FuncFinish(void*);
1510ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1520ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
153ba61063dSBarry Smith #endif
154e5c89e4eSSatish Balay 
155e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
156e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
157e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
158e5c89e4eSSatish Balay #endif
1592c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1607087cfbeSBarry Smith MPI_Datatype   MPI_C_DOUBLE_COMPLEX;
1617087cfbeSBarry Smith MPI_Datatype   MPI_C_COMPLEX;
1622c876bd9SBarry Smith #endif
1637087cfbeSBarry Smith PetscScalar    PETSC_i;
164e5c89e4eSSatish Balay #else
1657087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
166e5c89e4eSSatish Balay #endif
167ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
168c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
169c90a1750SBarry Smith #endif
1707087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1717087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
17275567043SBarry Smith 
173e5c89e4eSSatish Balay /*
174e5c89e4eSSatish Balay      These are needed by petscbt.h
175e5c89e4eSSatish Balay */
176c6db04a5SJed Brown #include <petscbt.h>
1777087cfbeSBarry Smith char      _BT_mask = ' ';
1787087cfbeSBarry Smith char      _BT_c = ' ';
1797087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
180e5c89e4eSSatish Balay 
181e5c89e4eSSatish Balay /*
182e5c89e4eSSatish Balay        Function that is called to display all error messages
183e5c89e4eSSatish Balay */
1847087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1857087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
186238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1877087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
188238ccf28SShri Abhyankar #else
1897087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
190238ccf28SShri Abhyankar #endif
191bab1f7e6SVictor Minden /*
1928154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1938154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
194bab1f7e6SVictor Minden 
195e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
196e5c89e4eSSatish Balay /*
197e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
198e5c89e4eSSatish Balay */
199e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
200e5c89e4eSSatish Balay 
201e5c89e4eSSatish Balay #undef __FUNCT__
202f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
2037087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
204e5c89e4eSSatish Balay {
205e5c89e4eSSatish Balay   PetscErrorCode ierr;
206e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
207e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
208e5c89e4eSSatish Balay   char           version[256];
209e5c89e4eSSatish Balay 
210e5c89e4eSSatish Balay   PetscFunctionBegin;
211e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
212e5c89e4eSSatish Balay   if (!rank) {
213e5c89e4eSSatish Balay     char        arch[10];
214f56c2debSBarry Smith     int         err;
21588c29154SBarry Smith     PetscViewer viewer;
216f56c2debSBarry Smith 
217e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
218e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
219a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
220e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
221e5c89e4eSSatish Balay     if (filename) {
222e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
223e5c89e4eSSatish Balay     } else {
224e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
225e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
226e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
227e5c89e4eSSatish Balay     }
228e5c89e4eSSatish Balay 
229e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
230e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
231e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
232e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
233e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
23488c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
23588c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2366bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
237e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
238f56c2debSBarry Smith     err = fflush(*fd);
239e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
240e5c89e4eSSatish Balay   }
241e5c89e4eSSatish Balay   PetscFunctionReturn(0);
242e5c89e4eSSatish Balay }
243e5c89e4eSSatish Balay 
244e5c89e4eSSatish Balay #undef __FUNCT__
245f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2467087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
247e5c89e4eSSatish Balay {
248e5c89e4eSSatish Balay   PetscErrorCode ierr;
249e5c89e4eSSatish Balay   PetscMPIInt    rank;
250e5c89e4eSSatish Balay   char           date[64];
251f56c2debSBarry Smith   int            err;
252e5c89e4eSSatish Balay 
253e5c89e4eSSatish Balay   PetscFunctionBegin;
254e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
255e5c89e4eSSatish Balay   if (!rank) {
256e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
257e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
258e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
259e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
260f56c2debSBarry Smith     err = fflush(*fd);
261e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
262f56c2debSBarry Smith     err = fclose(*fd);
263e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
264e5c89e4eSSatish Balay   }
265e5c89e4eSSatish Balay   PetscFunctionReturn(0);
266e5c89e4eSSatish Balay }
267e5c89e4eSSatish Balay 
268e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
269e5c89e4eSSatish Balay 
270e5c89e4eSSatish Balay /*
271e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
272e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
273e5c89e4eSSatish Balay 
274e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2753c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2763c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
277e5c89e4eSSatish Balay */
278e5c89e4eSSatish Balay 
279e5c89e4eSSatish Balay #undef __FUNCT__
280e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
281e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
282e5c89e4eSSatish Balay {
283e5c89e4eSSatish Balay   PetscFunctionBegin;
2843c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
285e5c89e4eSSatish Balay   abort();
286e5c89e4eSSatish Balay }
287e5c89e4eSSatish Balay 
288e5c89e4eSSatish Balay #undef __FUNCT__
289e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
290e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
291e5c89e4eSSatish Balay {
292e5c89e4eSSatish Balay   PetscErrorCode ierr;
293e5c89e4eSSatish Balay 
294e5c89e4eSSatish Balay   PetscFunctionBegin;
2953c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
296e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
297e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
2983c311c98SBarry Smith     MPI_Abort(*comm,*flag);
299e5c89e4eSSatish Balay   }
300e5c89e4eSSatish Balay }
301e5c89e4eSSatish Balay 
302e5c89e4eSSatish Balay #undef __FUNCT__
303e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
304e5c89e4eSSatish Balay /*@C
305e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
306e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
307e5c89e4eSSatish Balay 
308e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
309e5c89e4eSSatish Balay 
310e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
311e5c89e4eSSatish Balay 
312e5c89e4eSSatish Balay    Level: advanced
313e5c89e4eSSatish Balay 
314e5c89e4eSSatish Balay    Note:
315e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
316e5c89e4eSSatish Balay 
31788c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
318e5c89e4eSSatish Balay @*/
3197087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
320e5c89e4eSSatish Balay {
321e5c89e4eSSatish Balay   PetscFunctionBegin;
322e5c89e4eSSatish Balay   PetscFinalize();
323e5c89e4eSSatish Balay   exit(0);
324e5c89e4eSSatish Balay   return 0;
325e5c89e4eSSatish Balay }
326e5c89e4eSSatish Balay 
327ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
32809573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
329ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
330e5c89e4eSSatish Balay static char       emacsmachinename[256];
331e5c89e4eSSatish Balay 
332e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
333e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
334e5c89e4eSSatish Balay 
335e5c89e4eSSatish Balay #undef __FUNCT__
336e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
337e5c89e4eSSatish Balay /*@C
338e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
339e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
340e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
341e5c89e4eSSatish Balay 
342e5c89e4eSSatish Balay    Input Parameter:
343e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
344da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
345e5c89e4eSSatish Balay 
346e5c89e4eSSatish Balay    Level: developer
347e5c89e4eSSatish Balay 
348e5c89e4eSSatish Balay    Concepts: package help message
349e5c89e4eSSatish Balay 
350e5c89e4eSSatish Balay @*/
3517087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
352e5c89e4eSSatish Balay {
353e5c89e4eSSatish Balay   PetscFunctionBegin;
354e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
355e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
356e5c89e4eSSatish Balay   PetscFunctionReturn(0);
357e5c89e4eSSatish Balay }
358e5c89e4eSSatish Balay 
359e5c89e4eSSatish Balay #undef __FUNCT__
360e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3617087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
362e5c89e4eSSatish Balay {
363e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
364e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
365ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
366e5c89e4eSSatish Balay   PetscErrorCode ierr;
367a6d0e24fSJed Brown   PetscReal      si;
368e5c89e4eSSatish Balay   int            i;
369e5c89e4eSSatish Balay   PetscMPIInt    rank;
370e5c89e4eSSatish Balay   char           version[256];
371e5c89e4eSSatish Balay 
372e5c89e4eSSatish Balay   PetscFunctionBegin;
373e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
374e5c89e4eSSatish Balay 
375e5c89e4eSSatish Balay   /*
376e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
377e5c89e4eSSatish Balay   */
3788bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
37981b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
380acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
381e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
382555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
383555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
384555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
385555d055bSBarry Smith #endif
386e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
387555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
388555d055bSBarry Smith     }
389555d055bSBarry Smith #endif
390e5c89e4eSSatish Balay   }
391e5c89e4eSSatish Balay #else
392acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
393acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
394e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
395e5c89e4eSSatish Balay #endif
396e5c89e4eSSatish Balay   if (flg3) {
397e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
398e5c89e4eSSatish Balay   }
39990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
400acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
401e5c89e4eSSatish Balay   if (flg1) {
402e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
403e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
404e5c89e4eSSatish Balay   }
405e5c89e4eSSatish Balay 
40690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
407acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4087783f70dSSatish Balay   if (!flg1) {
40990d69ab7SBarry Smith     flg1 = PETSC_FALSE;
410acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4117783f70dSSatish Balay   }
412e5c89e4eSSatish Balay   if (flg1) {
413e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
414e5c89e4eSSatish Balay   }
415e5c89e4eSSatish Balay 
416e5c89e4eSSatish Balay   /*
417e5c89e4eSSatish Balay       Set the display variable for graphics
418e5c89e4eSSatish Balay   */
419e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
420e5c89e4eSSatish Balay 
42151dcc849SKerry Stevens 
42251dcc849SKerry Stevens   /*
423e5c89e4eSSatish Balay       Print the PETSc version information
424e5c89e4eSSatish Balay   */
425e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
426e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
427e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
428e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
429e5c89e4eSSatish Balay 
430e5c89e4eSSatish Balay     /*
431e5c89e4eSSatish Balay        Print "higher-level" package version message
432e5c89e4eSSatish Balay     */
433e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
434e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
435e5c89e4eSSatish Balay     }
436e5c89e4eSSatish Balay 
437a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
438e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
439e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
440e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
441e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
442e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
44384e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
444e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
445e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
446e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
447e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
448e5c89e4eSSatish Balay   }
449e5c89e4eSSatish Balay 
450e5c89e4eSSatish Balay   /*
451e5c89e4eSSatish Balay        Print "higher-level" package help message
452e5c89e4eSSatish Balay   */
453e5c89e4eSSatish Balay   if (flg3){
454e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
455e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
456e5c89e4eSSatish Balay     }
457e5c89e4eSSatish Balay   }
458e5c89e4eSSatish Balay 
459e5c89e4eSSatish Balay   /*
460e5c89e4eSSatish Balay       Setup the error handling
461e5c89e4eSSatish Balay   */
46290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
463acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
464cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
46590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
466acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
467cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
46890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
469acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
470e5c89e4eSSatish Balay   if (flg1) {
471e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
472e5c89e4eSSatish Balay   }
47390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
474acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
475cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
47696cc47afSJed Brown   flg1 = PETSC_FALSE;
477acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
47896cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
479e5c89e4eSSatish Balay 
480e5c89e4eSSatish Balay   /*
481e5c89e4eSSatish Balay       Setup debugger information
482e5c89e4eSSatish Balay   */
483e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
484e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
485e5c89e4eSSatish Balay   if (flg1) {
486e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
487e5c89e4eSSatish Balay 
488e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
489e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
490e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
491e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
492e5c89e4eSSatish Balay   }
4935e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
4945e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
495e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
496e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
497e5c89e4eSSatish Balay   if (flg1 || flg2) {
498e5c89e4eSSatish Balay     PetscMPIInt    size;
499e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
500e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
501e5c89e4eSSatish Balay     /*
502e5c89e4eSSatish Balay        we have to make sure that all processors have opened
503e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
504e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
505e5c89e4eSSatish Balay        and kill the program.
506e5c89e4eSSatish Balay     */
507e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
508e5c89e4eSSatish Balay     if (size > 2) {
509533163c2SBarry Smith       PetscMPIInt dummy = 0;
510e5c89e4eSSatish Balay       MPI_Status  status;
511e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
512e5c89e4eSSatish Balay         if (rank != i) {
513e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
514e5c89e4eSSatish Balay         }
515e5c89e4eSSatish Balay       }
516e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
517e5c89e4eSSatish Balay         if (rank != i) {
518e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
519e5c89e4eSSatish Balay         }
520e5c89e4eSSatish Balay       }
521e5c89e4eSSatish Balay     }
522e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
523e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
524e5c89e4eSSatish Balay     lsize = size;
525e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
526e5c89e4eSSatish Balay     if (flag) {
527e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
528e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
529e5c89e4eSSatish Balay       }
530e5c89e4eSSatish Balay     }
531e5c89e4eSSatish Balay     if (!flag) {
532e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
533e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
534e5c89e4eSSatish Balay       if (flg1) {
535e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
536e5c89e4eSSatish Balay       } else {
537e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
538e5c89e4eSSatish Balay       }
539e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
540e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
541e5c89e4eSSatish Balay     }
542e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
543e5c89e4eSSatish Balay   }
544e5c89e4eSSatish Balay 
545e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
546cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
547e5c89e4eSSatish Balay 
54893ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
54922b84c2fSbcordonn   /*
55022b84c2fSbcordonn     Activates new sockets for zope if needed
55122b84c2fSbcordonn   */
55284ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
553d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5546dc8fec2Sbcordonn   if (flgz){
55522b84c2fSbcordonn     int  sockfd;
556f1384234SBarry Smith     char hostname[256];
55722b84c2fSbcordonn     char username[256];
5586dc8fec2Sbcordonn     int  remoteport = 9999;
5599c4c166aSBarry Smith 
56084ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
56184ab5442Sbcordonn     if (!hostname[0]){
5629c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5639c4c166aSBarry Smith     }
56422b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5659c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
56622b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
56722b84c2fSbcordonn     if (flgzout){
56822b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
569606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5706dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5719c4c166aSBarry Smith     } else {
572d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
573d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5749c4c166aSBarry Smith     }
5759c4c166aSBarry Smith   }
57693ba235fSBarry Smith #endif
577ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
578ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
579ffc871a5SBarry Smith   if (flgz){
580ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
581ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
582ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
583ffc871a5SBarry Smith   }
584ffc871a5SBarry Smith #endif
5856dc8fec2Sbcordonn 
586e5c89e4eSSatish Balay   /*
587e5c89e4eSSatish Balay         Setup profiling and logging
588e5c89e4eSSatish Balay   */
5896cf91177SBarry Smith #if defined (PETSC_USE_INFO)
5908bb29257SSatish Balay   {
591e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
5926cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
5938bb29257SSatish Balay     if (flg1 && logname[0]) {
594fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
5958bb29257SSatish Balay     } else if (flg1) {
596fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
597e5c89e4eSSatish Balay     }
598e5c89e4eSSatish Balay   }
599865f6aa8SSatish Balay #endif
600865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
601865f6aa8SSatish Balay   mname[0] = 0;
602f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
603865f6aa8SSatish Balay   if (flg1) {
604865f6aa8SSatish Balay     if (mname[0]) {
605f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
606865f6aa8SSatish Balay     } else {
607f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
608865f6aa8SSatish Balay     }
609865f6aa8SSatish Balay   }
610e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
61190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
612fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
613e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
614e5c89e4eSSatish Balay #endif
61590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
61690d69ab7SBarry Smith   flg2 = PETSC_FALSE;
61790d69ab7SBarry Smith   flg3 = PETSC_FALSE;
618acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
619acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
620d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6219f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
622e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6239f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
624e5c89e4eSSatish Balay 
625e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
626e5c89e4eSSatish Balay   if (flg1) {
627e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
628e5c89e4eSSatish Balay     FILE *file;
629e5c89e4eSSatish Balay     if (mname[0]) {
630e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
631e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
632e5c89e4eSSatish Balay       file = fopen(fname,"w");
633f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
634e5c89e4eSSatish Balay     } else {
635da9f1d6bSBarry Smith       file = PETSC_STDOUT;
636e5c89e4eSSatish Balay     }
637e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
638e5c89e4eSSatish Balay   }
639e5c89e4eSSatish Balay #endif
640e5c89e4eSSatish Balay 
641e5c89e4eSSatish Balay   /*
642e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
643e5c89e4eSSatish Balay   */
64463d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
645e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
646e5c89e4eSSatish Balay #endif
647e5c89e4eSSatish Balay 
648acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
649e5c89e4eSSatish Balay 
650ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
651af359df3SBarry Smith   /*
652af359df3SBarry Smith       Determine whether user specified maximum number of threads
653af359df3SBarry Smith    */
654af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
655af359df3SBarry Smith 
656b154f58aSKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-main",&flg1);CHKERRQ(ierr);
657b154f58aSKerry Stevens   if(flg1) {
658b154f58aSKerry Stevens     cpu_set_t mset;
659b154f58aSKerry Stevens     int icorr,ncorr = get_nprocs();
660b154f58aSKerry Stevens     ierr = PetscOptionsGetInt(PETSC_NULL,"-main",&icorr,PETSC_NULL);CHKERRQ(ierr);
661b154f58aSKerry Stevens     CPU_ZERO(&mset);
662b154f58aSKerry Stevens     CPU_SET(icorr%ncorr,&mset);
663b154f58aSKerry Stevens     sched_setaffinity(0,sizeof(cpu_set_t),&mset);
664b154f58aSKerry Stevens   }
665b154f58aSKerry Stevens 
666af359df3SBarry Smith   PetscInt N_CORES = get_nprocs();
667af359df3SBarry Smith   ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
668af359df3SBarry Smith   char tstr[9];
669af359df3SBarry Smith   char tbuf[2];
670af359df3SBarry Smith   strcpy(tstr,"-thread");
671af359df3SBarry Smith   for(i=0;i<PetscMaxThreads;i++) {
672af359df3SBarry Smith     ThreadCoreAffinity[i] = i;
673af359df3SBarry Smith     sprintf(tbuf,"%d",i);
674af359df3SBarry Smith     strcat(tstr,tbuf);
675af359df3SBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
676af359df3SBarry Smith     if(flg1) {
677af359df3SBarry Smith       ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
678af359df3SBarry Smith       ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
679af359df3SBarry Smith     }
680af359df3SBarry Smith     tstr[7] = '\0';
681af359df3SBarry Smith   }
682cfcfc605SKerry Stevens 
683cfcfc605SKerry Stevens   /*
684cfcfc605SKerry Stevens       Determine whether to use thread pool
685cfcfc605SKerry Stevens    */
686cfcfc605SKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
687cfcfc605SKerry Stevens   if (flg1) {
688cfcfc605SKerry Stevens     PetscUseThreadPool = PETSC_TRUE;
689af359df3SBarry Smith     /* get the thread pool type */
690af359df3SBarry Smith     PetscInt ipool = 0;
691af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
692af359df3SBarry Smith 
693af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
694af359df3SBarry Smith     switch(ipool) {
695af359df3SBarry Smith     case 1:
696af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
697af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
698af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
699af359df3SBarry Smith       MainWait              = &MainWait_Tree;
700af359df3SBarry Smith       MainJob               = &MainJob_Tree;
701af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
702af359df3SBarry Smith       break;
703af359df3SBarry Smith     case 2:
704af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
705af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
706af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
707af359df3SBarry Smith       MainWait              = &MainWait_Main;
708af359df3SBarry Smith       MainJob               = &MainJob_Main;
709af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
710af359df3SBarry Smith       break;
711af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
712af359df3SBarry Smith     case 3:
713af359df3SBarry Smith #else
714af359df3SBarry Smith     default:
715af359df3SBarry Smith #endif
716af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
717af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
718af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
719af359df3SBarry Smith       MainWait              = &MainWait_Chain;
720af359df3SBarry Smith       MainJob               = &MainJob_Chain;
721af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
722af359df3SBarry Smith       break;
723af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
724af359df3SBarry Smith     default:
725af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
726af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
727af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
728af359df3SBarry Smith       MainWait              = &MainWait_True;
729af359df3SBarry Smith       MainJob               = &MainJob_True;
730af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
731af359df3SBarry Smith       break;
732af359df3SBarry Smith #endif
733af359df3SBarry Smith     }
734af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
735683509dcSKerry Stevens   } else {
736683509dcSKerry Stevens     //need to define these in the case on 'no threads' or 'thread create/destroy'
737683509dcSKerry Stevens     //could take any of the above versions
738683509dcSKerry Stevens     MainJob               = &MainJob_Spawn;
739af359df3SBarry Smith   }
740af359df3SBarry Smith #endif
741e5c89e4eSSatish Balay   /*
742e5c89e4eSSatish Balay        Print basic help message
743e5c89e4eSSatish Balay   */
744e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
745e5c89e4eSSatish Balay   if (flg1) {
746e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
747301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
748301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
749301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
750e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
751e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
752e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
753e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
754e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
755e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
756e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
757e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
758e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
759e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
760e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
761e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
762e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
763e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
764e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
765e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
766e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
767e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
768e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7694161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7704161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
771e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
772e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
773e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
774e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
775e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
776e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
777a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
778e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
77940ab9619SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -server <port>: Run PETSc webserver (default port is 8080) see PetscWebServe()\n");CHKERRQ(ierr);
780e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
781e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
7827ef452c0SMatthew G Knepley     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary _summary_python]: logging objects and events\n");CHKERRQ(ierr);
783e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
784e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
785e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
786e5c89e4eSSatish Balay #endif
7876cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
788e5c89e4eSSatish Balay #endif
789e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
790e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
791e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
792e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
793e5c89e4eSSatish Balay   }
794e5c89e4eSSatish Balay 
795a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
796e5c89e4eSSatish Balay   if (flg1) {
797e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
798e5c89e4eSSatish Balay   }
799e5c89e4eSSatish Balay 
8006cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
801e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
802e5c89e4eSSatish Balay   if (f) {
8036cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
804e5c89e4eSSatish Balay   }
805827f890bSBarry Smith 
8068154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
807c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
80873113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
80973113deaSBarry Smith   else flg1 = PETSC_FALSE;
8108154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
8118154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
812bab1f7e6SVictor Minden #endif
813192daf7cSBarry Smith 
814e5c89e4eSSatish Balay   PetscFunctionReturn(0);
815e5c89e4eSSatish Balay }
816df413903SBarry Smith 
817ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
818ba61063dSBarry Smith 
81951d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
82051d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
82151d315f7SKerry Stevens   PetscErrorCode iterr;
82251d315f7SKerry Stevens   int icorr,ierr;
82351d315f7SKerry Stevens   int* pId = (int*)arg;
82451d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
82551d315f7SKerry Stevens   PetscBool PeeOn;
82651d315f7SKerry Stevens   cpu_set_t mset;
8279e800a48SKerry Stevens   //printf("Thread %d In Tree Thread Function\n",ThreadId);
82851d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
82951d315f7SKerry Stevens   CPU_ZERO(&mset);
83051d315f7SKerry Stevens   CPU_SET(icorr,&mset);
83151d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
83251d315f7SKerry Stevens 
83351d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
83451d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
83551d315f7SKerry Stevens   }
83651d315f7SKerry Stevens   else {
83751d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
83851d315f7SKerry Stevens   }
83951d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
840ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
84151d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
84251d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
84351d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
84451d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
84551d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
846ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
847ba61063dSBarry Smith            upon return, has the lock */
84851d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
84951d315f7SKerry Stevens         }
85051d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
85151d315f7SKerry Stevens       }
85251d315f7SKerry Stevens     }
853ba61063dSBarry Smith     /* your subordinates are now ready */
85451d315f7SKerry Stevens   }
85551d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
856ba61063dSBarry Smith   /* update your ready status */
85751d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
85851d315f7SKerry Stevens   if(ThreadId==0) {
85951d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
860ba61063dSBarry Smith     /* ignal main */
86151d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
86251d315f7SKerry Stevens   }
86351d315f7SKerry Stevens   else {
864ba61063dSBarry Smith     /* tell your boss that you're ready to work */
86551d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
86651d315f7SKerry Stevens   }
867ba61063dSBarry Smith   /* the while loop needs to have an exit
868ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
869ba61063dSBarry Smith    and calling FuncFinish */
87051d315f7SKerry Stevens   while(PetscThreadGo) {
871ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
872ba61063dSBarry Smith       waiting when you don't have to causes problems
873ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
87451d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
875ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
876ba61063dSBarry Smith        upon return, has the lock */
87751d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
87851d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
87951d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
88051d315f7SKerry Stevens     }
88151d315f7SKerry Stevens     if(ThreadId==0) {
88251d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
88351d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
88451d315f7SKerry Stevens     }
88551d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
88651d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
887ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
88851d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
88951d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
89051d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
89151d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
89251d315f7SKerry Stevens         }
89351d315f7SKerry Stevens       }
89451d315f7SKerry Stevens     }
895ba61063dSBarry Smith     /* do your job */
89651d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
89751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
89851d315f7SKerry Stevens     }
89951d315f7SKerry Stevens     else {
90051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
90151d315f7SKerry Stevens     }
90251d315f7SKerry Stevens     if(iterr!=0) {
90351d315f7SKerry Stevens       ithreaderr = 1;
90451d315f7SKerry Stevens     }
90551d315f7SKerry Stevens     if(PetscThreadGo) {
906ba61063dSBarry Smith       /* reset job, get ready for more */
90751d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
908ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
909ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
91051d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
91151d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
91251d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
91351d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
91451d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
915ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
916ba61063dSBarry Smith                upon return, has the lock */
91751d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
91851d315f7SKerry Stevens             }
91951d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
92051d315f7SKerry Stevens           }
92151d315f7SKerry Stevens 	}
922ba61063dSBarry Smith         /* your subordinates are now ready */
92351d315f7SKerry Stevens       }
92451d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
92551d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
92651d315f7SKerry Stevens       if(ThreadId==0) {
927ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
928ba61063dSBarry Smith         /* root thread signals 'main' */
92951d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
93051d315f7SKerry Stevens       }
93151d315f7SKerry Stevens       else {
932ba61063dSBarry Smith         /* signal your boss before you go to sleep */
93351d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
93451d315f7SKerry Stevens       }
93551d315f7SKerry Stevens     }
93651d315f7SKerry Stevens   }
93751d315f7SKerry Stevens   return NULL;
93851d315f7SKerry Stevens }
93951d315f7SKerry Stevens 
94051d315f7SKerry Stevens #undef __FUNCT__
94151d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
94251d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
94351d315f7SKerry Stevens   PetscInt i,ierr;
94451d315f7SKerry Stevens   int status;
94551d315f7SKerry Stevens 
94651d315f7SKerry Stevens   if(PetscUseThreadPool) {
94751d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
94851d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
94951d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
95051d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
95151d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
95251d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
95351d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
95451d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
95551d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95651d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95751d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
95851d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
959ba61063dSBarry Smith     /* initialize job structure */
96051d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
96151d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
96251d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
96351d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
96451d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
96551d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
96651d315f7SKerry Stevens     }
96751d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
96851d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
96951d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
97051d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
97151d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
97251d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
97351d315f7SKerry Stevens     }
97451d315f7SKerry Stevens     job_tree.pfunc = NULL;
97551d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
97651d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
97751d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
97851d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
979ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
98051d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
981ba61063dSBarry Smith     /* create threads */
98251d315f7SKerry Stevens     for(i=0; i<N; i++) {
98351d315f7SKerry Stevens       pVal[i] = i;
98451d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
985ba61063dSBarry Smith       /* should check status */
98651d315f7SKerry Stevens     }
98751d315f7SKerry Stevens   }
98851d315f7SKerry Stevens   return NULL;
98951d315f7SKerry Stevens }
99051d315f7SKerry Stevens 
99151d315f7SKerry Stevens #undef __FUNCT__
99251d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
99351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
99451d315f7SKerry Stevens   int i,ierr;
99551d315f7SKerry Stevens   void* jstatus;
99651d315f7SKerry Stevens 
99751d315f7SKerry Stevens   PetscFunctionBegin;
99851d315f7SKerry Stevens 
999ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1000ba61063dSBarry Smith   /* join the threads */
100151d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
100251d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1003ba61063dSBarry Smith     /* do error checking*/
100451d315f7SKerry Stevens   }
100551d315f7SKerry Stevens   free(PetscThreadPoint);
100651d315f7SKerry Stevens   free(arrmutex);
100751d315f7SKerry Stevens   free(arrcond1);
100851d315f7SKerry Stevens   free(arrcond2);
100951d315f7SKerry Stevens   free(arrstart);
101051d315f7SKerry Stevens   free(arrready);
101151d315f7SKerry Stevens   free(job_tree.pdata);
101251d315f7SKerry Stevens   free(pVal);
1013cfcfc605SKerry Stevens 
101451d315f7SKerry Stevens   PetscFunctionReturn(0);
101551d315f7SKerry Stevens }
101651d315f7SKerry Stevens 
101751d315f7SKerry Stevens #undef __FUNCT__
101851d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
101951d315f7SKerry Stevens void MainWait_Tree() {
102051d315f7SKerry Stevens   int ierr;
102151d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
102251d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
102351d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
102451d315f7SKerry Stevens   }
102551d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
102651d315f7SKerry Stevens }
102751d315f7SKerry Stevens 
102851d315f7SKerry Stevens #undef __FUNCT__
102951d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
103051d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
103151d315f7SKerry Stevens   int i,ierr;
103251d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
103336d20dc5SKerry Stevens 
103451d315f7SKerry Stevens   MainWait();
103551d315f7SKerry Stevens   job_tree.pfunc = pFunc;
103651d315f7SKerry Stevens   job_tree.pdata = data;
103751d315f7SKerry Stevens   job_tree.startJob = PETSC_TRUE;
103851d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
103951d315f7SKerry Stevens     *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
104051d315f7SKerry Stevens   }
104151d315f7SKerry Stevens   job_tree.eJobStat = JobInitiated;
104251d315f7SKerry Stevens   ierr = pthread_cond_signal(job_tree.cond2array[0]);
104351d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1044ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
104551d315f7SKerry Stevens   }
104636d20dc5SKerry Stevens 
104751d315f7SKerry Stevens   if(ithreaderr) {
104851d315f7SKerry Stevens     ijoberr = ithreaderr;
104951d315f7SKerry Stevens   }
105051d315f7SKerry Stevens   return ijoberr;
105151d315f7SKerry Stevens }
105251d315f7SKerry Stevens /****  ****/
105351d315f7SKerry Stevens 
105451d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
105551d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
105651d315f7SKerry Stevens   PetscErrorCode iterr;
105751d315f7SKerry Stevens   int icorr,ierr;
105851d315f7SKerry Stevens   int* pId = (int*)arg;
105951d315f7SKerry Stevens   int ThreadId = *pId;
106051d315f7SKerry Stevens   cpu_set_t mset;
10619e800a48SKerry Stevens   //printf("Thread %d In Main Thread Function\n",ThreadId);
106251d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
106351d315f7SKerry Stevens   CPU_ZERO(&mset);
106451d315f7SKerry Stevens   CPU_SET(icorr,&mset);
106551d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
106651d315f7SKerry Stevens 
106751d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1068ba61063dSBarry Smith   /* update your ready status */
106951d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1070ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
107151d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
107251d315f7SKerry Stevens 
1073ba61063dSBarry Smith   /* the while loop needs to have an exit
1074ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1075ba61063dSBarry Smith      and calling FuncFinish */
107651d315f7SKerry Stevens   while(PetscThreadGo) {
1077ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1078ba61063dSBarry Smith        waiting when you don't have to causes problems
1079ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
108051d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1081ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1082ba61063dSBarry Smith        upon return, has the lock */
108351d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1084ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
108551d315f7SKerry Stevens     }
108651d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
108751d315f7SKerry Stevens     if(job_main.pdata==NULL) {
108851d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
108951d315f7SKerry Stevens     }
109051d315f7SKerry Stevens     else {
109151d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
109251d315f7SKerry Stevens     }
109351d315f7SKerry Stevens     if(iterr!=0) {
109451d315f7SKerry Stevens       ithreaderr = 1;
109551d315f7SKerry Stevens     }
109651d315f7SKerry Stevens     if(PetscThreadGo) {
1097ba61063dSBarry Smith       /* reset job, get ready for more */
109851d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
109951d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1100ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
110151d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
110251d315f7SKerry Stevens     }
110351d315f7SKerry Stevens   }
110451d315f7SKerry Stevens   return NULL;
110551d315f7SKerry Stevens }
110651d315f7SKerry Stevens 
110751d315f7SKerry Stevens #undef __FUNCT__
110851d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
110951d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
111051d315f7SKerry Stevens   PetscInt i,ierr;
111151d315f7SKerry Stevens   int status;
111251d315f7SKerry Stevens 
111351d315f7SKerry Stevens   if(PetscUseThreadPool) {
111451d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
111551d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
111651d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
111751d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
111851d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
111951d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
112051d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
112151d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
112251d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112351d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112451d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1125ba61063dSBarry Smith     /* initialize job structure */
112651d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
112751d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
112851d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
112951d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
113051d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
113151d315f7SKerry Stevens     }
113251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
113351d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
113451d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
113551d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
113651d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
113751d315f7SKerry Stevens     }
113851d315f7SKerry Stevens     job_main.pfunc = NULL;
113951d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
114051d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1141ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
114251d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1143ba61063dSBarry Smith     /* create threads */
114451d315f7SKerry Stevens     for(i=0; i<N; i++) {
114551d315f7SKerry Stevens       pVal[i] = i;
114651d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1147ba61063dSBarry Smith       /* error check */
114851d315f7SKerry Stevens     }
114951d315f7SKerry Stevens   }
115051d315f7SKerry Stevens   else {
115151d315f7SKerry Stevens   }
115251d315f7SKerry Stevens   return NULL;
115351d315f7SKerry Stevens }
115451d315f7SKerry Stevens 
115551d315f7SKerry Stevens #undef __FUNCT__
115651d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
115751d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
115851d315f7SKerry Stevens   int i,ierr;
115951d315f7SKerry Stevens   void* jstatus;
116051d315f7SKerry Stevens 
116151d315f7SKerry Stevens   PetscFunctionBegin;
116251d315f7SKerry Stevens 
1163ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1164ba61063dSBarry Smith   /* join the threads */
116551d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1166ba61063dSBarry Smith     ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
116751d315f7SKerry Stevens   }
116851d315f7SKerry Stevens   free(PetscThreadPoint);
116951d315f7SKerry Stevens   free(arrmutex);
117051d315f7SKerry Stevens   free(arrcond1);
117151d315f7SKerry Stevens   free(arrcond2);
117251d315f7SKerry Stevens   free(arrstart);
117351d315f7SKerry Stevens   free(arrready);
117451d315f7SKerry Stevens   free(job_main.pdata);
117551d315f7SKerry Stevens   free(pVal);
1176cfcfc605SKerry Stevens 
117751d315f7SKerry Stevens   PetscFunctionReturn(0);
117851d315f7SKerry Stevens }
117951d315f7SKerry Stevens 
118051d315f7SKerry Stevens #undef __FUNCT__
118151d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
118251d315f7SKerry Stevens void MainWait_Main() {
118351d315f7SKerry Stevens   int i,ierr;
118451d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
118551d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
118651d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
118751d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
118851d315f7SKerry Stevens     }
118951d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
119051d315f7SKerry Stevens   }
119151d315f7SKerry Stevens }
119251d315f7SKerry Stevens 
119351d315f7SKerry Stevens #undef __FUNCT__
119451d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
119551d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
119651d315f7SKerry Stevens   int i,ierr;
119751d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
119836d20dc5SKerry Stevens 
1199ba61063dSBarry Smith   MainWait(); /* you know everyone is waiting to be signalled! */
120051d315f7SKerry Stevens   job_main.pfunc = pFunc;
120151d315f7SKerry Stevens   job_main.pdata = data;
120251d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1203ba61063dSBarry Smith     *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
120451d315f7SKerry Stevens   }
1205ba61063dSBarry Smith   /* tell the threads to go to work */
120651d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
120751d315f7SKerry Stevens     ierr = pthread_cond_signal(job_main.cond2array[i]);
120851d315f7SKerry Stevens   }
120951d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1210ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
121151d315f7SKerry Stevens   }
121236d20dc5SKerry Stevens 
121351d315f7SKerry Stevens   if(ithreaderr) {
121451d315f7SKerry Stevens     ijoberr = ithreaderr;
121551d315f7SKerry Stevens   }
121651d315f7SKerry Stevens   return ijoberr;
121751d315f7SKerry Stevens }
121851d315f7SKerry Stevens /****  ****/
121951d315f7SKerry Stevens 
122051d315f7SKerry Stevens /**** Chain Thread Functions ****/
122151d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
122251d315f7SKerry Stevens   PetscErrorCode iterr;
122351d315f7SKerry Stevens   int icorr,ierr;
122451d315f7SKerry Stevens   int* pId = (int*)arg;
122551d315f7SKerry Stevens   int ThreadId = *pId;
122651d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
122751d315f7SKerry Stevens   PetscBool PeeOn;
122851d315f7SKerry Stevens   cpu_set_t mset;
12299e800a48SKerry Stevens   //printf("Thread %d In Chain Thread Function\n",ThreadId);
123051d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
123151d315f7SKerry Stevens   CPU_ZERO(&mset);
123251d315f7SKerry Stevens   CPU_SET(icorr,&mset);
123351d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
123451d315f7SKerry Stevens 
123551d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
123651d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
123751d315f7SKerry Stevens   }
123851d315f7SKerry Stevens   else {
123951d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
124051d315f7SKerry Stevens   }
124151d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1242ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
124351d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
124451d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1245ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1246ba61063dSBarry Smith        upon return, has the lock */
124751d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
124851d315f7SKerry Stevens     }
124951d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1250ba61063dSBarry Smith     /* your subordinate is now ready*/
125151d315f7SKerry Stevens   }
125251d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1253ba61063dSBarry Smith   /* update your ready status */
125451d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
125551d315f7SKerry Stevens   if(ThreadId==0) {
125651d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1257ba61063dSBarry Smith     /* signal main */
125851d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
125951d315f7SKerry Stevens   }
126051d315f7SKerry Stevens   else {
1261ba61063dSBarry Smith     /* tell your boss that you're ready to work */
126251d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
126351d315f7SKerry Stevens   }
1264ba61063dSBarry Smith   /*  the while loop needs to have an exit
1265ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1266ba61063dSBarry Smith    and calling FuncFinish */
126751d315f7SKerry Stevens   while(PetscThreadGo) {
1268ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1269ba61063dSBarry Smith        waiting when you don't have to causes problems
1270ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
127151d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1272ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1273ba61063dSBarry Smith        upon return, has the lock */
127451d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
127551d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
127651d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
127751d315f7SKerry Stevens     }
127851d315f7SKerry Stevens     if(ThreadId==0) {
127951d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
128051d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
128151d315f7SKerry Stevens     }
128251d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
128351d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1284ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
128551d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
128651d315f7SKerry Stevens     }
1287ba61063dSBarry Smith     /* do your job */
128851d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
128951d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
129051d315f7SKerry Stevens     }
129151d315f7SKerry Stevens     else {
129251d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
129351d315f7SKerry Stevens     }
129451d315f7SKerry Stevens     if(iterr!=0) {
129551d315f7SKerry Stevens       ithreaderr = 1;
129651d315f7SKerry Stevens     }
129751d315f7SKerry Stevens     if(PetscThreadGo) {
1298ba61063dSBarry Smith       /* reset job, get ready for more */
129951d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1300ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1301ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
130251d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
130351d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1304ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1305ba61063dSBarry Smith            upon return, has the lock */
130651d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
130751d315f7SKerry Stevens         }
130851d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1309ba61063dSBarry Smith         /* your subordinate is now ready */
131051d315f7SKerry Stevens       }
131151d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
131251d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
131351d315f7SKerry Stevens       if(ThreadId==0) {
1314ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1315ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
131651d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
131751d315f7SKerry Stevens       }
131851d315f7SKerry Stevens       else {
1319ba61063dSBarry Smith         /* signal your boss before you go to sleep */
132051d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
132151d315f7SKerry Stevens       }
132251d315f7SKerry Stevens     }
132351d315f7SKerry Stevens   }
132451d315f7SKerry Stevens   return NULL;
132551d315f7SKerry Stevens }
132651d315f7SKerry Stevens 
132751d315f7SKerry Stevens #undef __FUNCT__
132851d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
132951d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
133051d315f7SKerry Stevens   PetscInt i,ierr;
133151d315f7SKerry Stevens   int status;
133251d315f7SKerry Stevens 
133351d315f7SKerry Stevens   if(PetscUseThreadPool) {
133451d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
133551d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
133651d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
133751d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
133851d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
133951d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
134051d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
134151d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
134251d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134351d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134451d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
134551d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1346ba61063dSBarry Smith     /* initialize job structure */
134751d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
134851d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
134951d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
135051d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
135151d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
135251d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
135351d315f7SKerry Stevens     }
135451d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
135551d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
135651d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
135751d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
135851d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
135951d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
136051d315f7SKerry Stevens     }
136151d315f7SKerry Stevens     job_chain.pfunc = NULL;
136251d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
136351d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
136451d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
136551d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1366ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
136751d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1368ba61063dSBarry Smith     /* create threads */
136951d315f7SKerry Stevens     for(i=0; i<N; i++) {
137051d315f7SKerry Stevens       pVal[i] = i;
137151d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1372ba61063dSBarry Smith       /* should check error */
137351d315f7SKerry Stevens     }
137451d315f7SKerry Stevens   }
137551d315f7SKerry Stevens   else {
137651d315f7SKerry Stevens   }
137751d315f7SKerry Stevens   return NULL;
137851d315f7SKerry Stevens }
137951d315f7SKerry Stevens 
138051d315f7SKerry Stevens 
138151d315f7SKerry Stevens #undef __FUNCT__
138251d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
138351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
138451d315f7SKerry Stevens   int i,ierr;
138551d315f7SKerry Stevens   void* jstatus;
138651d315f7SKerry Stevens 
138751d315f7SKerry Stevens   PetscFunctionBegin;
138851d315f7SKerry Stevens 
1389ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1390ba61063dSBarry Smith   /* join the threads */
139151d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
139251d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1393ba61063dSBarry Smith     /* should check error */
139451d315f7SKerry Stevens   }
139551d315f7SKerry Stevens   free(PetscThreadPoint);
139651d315f7SKerry Stevens   free(arrmutex);
139751d315f7SKerry Stevens   free(arrcond1);
139851d315f7SKerry Stevens   free(arrcond2);
139951d315f7SKerry Stevens   free(arrstart);
140051d315f7SKerry Stevens   free(arrready);
140151d315f7SKerry Stevens   free(job_chain.pdata);
140251d315f7SKerry Stevens   free(pVal);
1403cfcfc605SKerry Stevens 
140451d315f7SKerry Stevens   PetscFunctionReturn(0);
140551d315f7SKerry Stevens }
140651d315f7SKerry Stevens 
140751d315f7SKerry Stevens #undef __FUNCT__
140851d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
140951d315f7SKerry Stevens void MainWait_Chain() {
141051d315f7SKerry Stevens   int ierr;
141151d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
141251d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
141351d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
141451d315f7SKerry Stevens   }
141551d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
141651d315f7SKerry Stevens }
141751d315f7SKerry Stevens 
141851d315f7SKerry Stevens #undef __FUNCT__
141951d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
142051d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
142151d315f7SKerry Stevens   int i,ierr;
142251d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
142336d20dc5SKerry Stevens 
142451d315f7SKerry Stevens   MainWait();
142551d315f7SKerry Stevens   job_chain.pfunc = pFunc;
142651d315f7SKerry Stevens   job_chain.pdata = data;
142751d315f7SKerry Stevens   job_chain.startJob = PETSC_TRUE;
142851d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
142951d315f7SKerry Stevens     *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
143051d315f7SKerry Stevens   }
143151d315f7SKerry Stevens   job_chain.eJobStat = JobInitiated;
143251d315f7SKerry Stevens   ierr = pthread_cond_signal(job_chain.cond2array[0]);
143351d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1434ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
143551d315f7SKerry Stevens   }
143636d20dc5SKerry Stevens 
143751d315f7SKerry Stevens   if(ithreaderr) {
143851d315f7SKerry Stevens     ijoberr = ithreaderr;
143951d315f7SKerry Stevens   }
144051d315f7SKerry Stevens   return ijoberr;
144151d315f7SKerry Stevens }
144251d315f7SKerry Stevens /****  ****/
144351d315f7SKerry Stevens 
1444ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
144551d315f7SKerry Stevens /**** True Thread Functions ****/
144651d315f7SKerry Stevens void* PetscThreadFunc_True(void* arg) {
144751d315f7SKerry Stevens   int icorr,ierr,iVal;
144851dcc849SKerry Stevens   int* pId = (int*)arg;
144951dcc849SKerry Stevens   int ThreadId = *pId;
14500ca81413SKerry Stevens   PetscErrorCode iterr;
145151d315f7SKerry Stevens   cpu_set_t mset;
14529e800a48SKerry Stevens   //printf("Thread %d In True Pool Thread Function\n",ThreadId);
145351d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
145451d315f7SKerry Stevens   CPU_ZERO(&mset);
145551d315f7SKerry Stevens   CPU_SET(icorr,&mset);
145651d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
145751d315f7SKerry Stevens 
145851d315f7SKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
145951d315f7SKerry Stevens   job_true.iNumReadyThreads++;
146051d315f7SKerry Stevens   if(job_true.iNumReadyThreads==PetscMaxThreads) {
146151dcc849SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
146251dcc849SKerry Stevens   }
1463ba61063dSBarry Smith   /*the while loop needs to have an exit
1464ba61063dSBarry Smith     the 'main' thread can terminate all the threads by performing a broadcast
1465ba61063dSBarry Smith    and calling FuncFinish */
146651dcc849SKerry Stevens   while(PetscThreadGo) {
1467ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
1468ba61063dSBarry Smith       waiting when you don't have to causes problems
1469ba61063dSBarry Smith      also need to wait if another thread sneaks in and messes with the predicate */
147051d315f7SKerry Stevens     while(job_true.startJob==PETSC_FALSE&&job_true.iNumJobThreads==0) {
1471ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1472ba61063dSBarry Smith        upon return, has the lock */
14736c9b208dSKerry Stevens       //printf("Thread %d Going to Sleep!\n",ThreadId);
147451d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
147551dcc849SKerry Stevens     }
147651d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
147751d315f7SKerry Stevens     job_true.iNumJobThreads--;
147851d315f7SKerry Stevens     job_true.iNumReadyThreads--;
147951d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
148051d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
148151d315f7SKerry Stevens     if(job_true.pdata==NULL) {
148251d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
148351dcc849SKerry Stevens     }
148451dcc849SKerry Stevens     else {
148551d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
148651dcc849SKerry Stevens     }
14870ca81413SKerry Stevens     if(iterr!=0) {
14880ca81413SKerry Stevens       ithreaderr = 1;
14890ca81413SKerry Stevens     }
14906c9b208dSKerry Stevens     //printf("Thread %d Finished Job\n",ThreadId);
1491ba61063dSBarry Smith     /* the barrier is necessary BECAUSE: look at job_true.iNumReadyThreads
1492ba61063dSBarry Smith       what happens if a thread finishes before they all start? BAD!
1493ba61063dSBarry Smith      what happens if a thread finishes before any else start? BAD! */
1494ba61063dSBarry Smith     pthread_barrier_wait(job_true.pbarr); /* ensures all threads are finished */
1495ba61063dSBarry Smith     /* reset job */
149651dcc849SKerry Stevens     if(PetscThreadGo) {
149751d315f7SKerry Stevens       pthread_mutex_lock(&job_true.mutex);
149851d315f7SKerry Stevens       job_true.iNumReadyThreads++;
149951d315f7SKerry Stevens       if(job_true.iNumReadyThreads==PetscMaxThreads) {
1500ba61063dSBarry Smith 	/* signal the 'main' thread that the job is done! (only done once) */
150151dcc849SKerry Stevens 	ierr = pthread_cond_signal(&main_cond);
150251dcc849SKerry Stevens       }
150351dcc849SKerry Stevens     }
150451dcc849SKerry Stevens   }
150551dcc849SKerry Stevens   return NULL;
150651dcc849SKerry Stevens }
150751dcc849SKerry Stevens 
1508f09cb4aaSKerry Stevens #undef __FUNCT__
150951d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_True"
151051d315f7SKerry Stevens void* PetscThreadInitialize_True(PetscInt N) {
151151dcc849SKerry Stevens   PetscInt i;
151251dcc849SKerry Stevens   int status;
15130ca81413SKerry Stevens 
1514f09cb4aaSKerry Stevens   pVal = (int*)malloc(N*sizeof(int));
1515ba61063dSBarry Smith   /* allocate memory in the heap for the thread structure */
151651dcc849SKerry Stevens   PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1517ba61063dSBarry Smith   BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
151851d315f7SKerry Stevens   job_true.pdata = (void**)malloc(N*sizeof(void*));
151951dcc849SKerry Stevens   for(i=0; i<N; i++) {
1520f09cb4aaSKerry Stevens     pVal[i] = i;
1521f09cb4aaSKerry Stevens     status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1522ba61063dSBarry Smith     /* error check to ensure proper thread creation */
152351dcc849SKerry Stevens     status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1524ba61063dSBarry Smith     /* should check error */
152551dcc849SKerry Stevens   }
15266c9b208dSKerry Stevens   //printf("Finished True Thread Pool Initialization\n");
152751dcc849SKerry Stevens   return NULL;
152851dcc849SKerry Stevens }
152951dcc849SKerry Stevens 
1530f09cb4aaSKerry Stevens 
1531f09cb4aaSKerry Stevens #undef __FUNCT__
153251d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
153351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
153451dcc849SKerry Stevens   int i,ierr;
153551dcc849SKerry Stevens   void* jstatus;
153651dcc849SKerry Stevens 
153751dcc849SKerry Stevens   PetscFunctionBegin;
1538cfcfc605SKerry Stevens 
1539ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1540ba61063dSBarry Smith   /* join the threads */
154151dcc849SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
154251dcc849SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
154351dcc849SKerry Stevens   }
154451dcc849SKerry Stevens   free(BarrPoint);
154551dcc849SKerry Stevens   free(PetscThreadPoint);
1546cfcfc605SKerry Stevens 
154751dcc849SKerry Stevens   PetscFunctionReturn(0);
154851dcc849SKerry Stevens }
154951dcc849SKerry Stevens 
1550f09cb4aaSKerry Stevens #undef __FUNCT__
155151d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
155251d315f7SKerry Stevens void MainWait_True() {
155351dcc849SKerry Stevens   int ierr;
15546c9b208dSKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
155551d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
155651d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
155751dcc849SKerry Stevens   }
155851d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
155951dcc849SKerry Stevens }
156051dcc849SKerry Stevens 
1561f09cb4aaSKerry Stevens #undef __FUNCT__
156251d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
156351d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
156451dcc849SKerry Stevens   int ierr;
15650ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
156636d20dc5SKerry Stevens 
15670ca81413SKerry Stevens   MainWait();
156851d315f7SKerry Stevens   job_true.pfunc = pFunc;
156951d315f7SKerry Stevens   job_true.pdata = data;
157051d315f7SKerry Stevens   job_true.pbarr = &BarrPoint[n];
157151d315f7SKerry Stevens   job_true.iNumJobThreads = n;
157251d315f7SKerry Stevens   job_true.startJob = PETSC_TRUE;
157351d315f7SKerry Stevens   ierr = pthread_cond_broadcast(&job_true.cond);
15740ca81413SKerry Stevens   if(pFunc!=FuncFinish) {
1575ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done */
15760ca81413SKerry Stevens   }
157736d20dc5SKerry Stevens 
15780ca81413SKerry Stevens   if(ithreaderr) {
15790ca81413SKerry Stevens     ijoberr = ithreaderr;
15800ca81413SKerry Stevens   }
15810ca81413SKerry Stevens   return ijoberr;
158251dcc849SKerry Stevens }
1583683509dcSKerry Stevens /**** NO THREAD POOL FUNCTION ****/
1584683509dcSKerry Stevens #undef __FUNCT__
1585683509dcSKerry Stevens #define __FUNCT__ "MainJob_Spawn"
1586683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void** data,PetscInt n) {
1587683509dcSKerry Stevens   PetscErrorCode ijoberr = 0;
1588683509dcSKerry Stevens 
1589683509dcSKerry Stevens   pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
1590cfcfc605SKerry Stevens   PetscThreadPoint = apThread; /* point to same place */
1591683509dcSKerry Stevens   PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1592683509dcSKerry Stevens   PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
1593683509dcSKerry Stevens   free(apThread);
1594683509dcSKerry Stevens 
1595683509dcSKerry Stevens   return ijoberr;
1596683509dcSKerry Stevens }
159751d315f7SKerry Stevens /****  ****/
1598ba61063dSBarry Smith #endif
159951dcc849SKerry Stevens 
160051dcc849SKerry Stevens void* FuncFinish(void* arg) {
160151dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16020ca81413SKerry Stevens   return(0);
160351dcc849SKerry Stevens }
1604ba61063dSBarry Smith 
1605ba61063dSBarry Smith #endif
1606