xref: /petsc/src/sys/objects/init.c (revision fdfc40db95c1000e3753cc51e039c4c5f46aeef0)
1e5c89e4eSSatish Balay /*
2e5c89e4eSSatish Balay 
3e5c89e4eSSatish Balay    This file defines part of the initialization of PETSc
4e5c89e4eSSatish Balay 
5e5c89e4eSSatish Balay   This file uses regular malloc and free because it cannot know
6e5c89e4eSSatish Balay   what malloc is being used until it has already processed the input.
7e5c89e4eSSatish Balay */
8ef386f4bSSatish Balay 
9ef386f4bSSatish Balay #include <petscsys.h>        /*I  "petscsys.h"   I*/
10ef386f4bSSatish Balay 
118f54b378SBarry Smith #ifndef _GNU_SOURCE
128f54b378SBarry Smith #define _GNU_SOURCE
138f54b378SBarry Smith #endif
14cd723cd1SBarry Smith #if defined(PETSC_HAVE_SCHED_H)
1524edbb7eSSatish Balay #ifndef __USE_GNU
1624edbb7eSSatish Balay #define __USE_GNU
1724edbb7eSSatish Balay #endif
188f54b378SBarry Smith #include <sched.h>
19cd723cd1SBarry Smith #endif
20bda8bf91SBarry Smith #if defined(PETSC_HAVE_PTHREAD_H)
2151dcc849SKerry Stevens #include <pthread.h>
22ba61063dSBarry Smith #endif
23ef386f4bSSatish Balay 
24ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
2551d315f7SKerry Stevens #include <sys/sysinfo.h>
26ba61063dSBarry Smith #endif
27121deb67SSatish Balay #if defined(PETSC_HAVE_UNISTD_H)
2851d315f7SKerry Stevens #include <unistd.h>
29121deb67SSatish Balay #endif
30e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
31e5c89e4eSSatish Balay #include <stdlib.h>
32e5c89e4eSSatish Balay #endif
33e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
34e5c89e4eSSatish Balay #include <malloc.h>
35e5c89e4eSSatish Balay #endif
36555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
37555d055bSBarry Smith #include <valgrind/valgrind.h>
38555d055bSBarry Smith #endif
39555d055bSBarry Smith 
40e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
41e5c89e4eSSatish Balay /*
42e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
43e5c89e4eSSatish Balay    already started before PETSc was initialized.
44e5c89e4eSSatish Balay */
457087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
467087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
477087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
4851dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
4951dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
507087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
517087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
52ba61063dSBarry Smith 
53ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
5451dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
5551dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
56ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
57ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
58ba61063dSBarry Smith #endif
5951d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
60f09cb4aaSKerry Stevens int*         pVal;
6151dcc849SKerry Stevens 
62ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
6351d315f7SKerry Stevens int* ThreadCoreAffinity;
6451d315f7SKerry Stevens 
65ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
6651d315f7SKerry Stevens 
6751d315f7SKerry Stevens typedef struct {
6851d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
6951d315f7SKerry Stevens   pthread_cond_t**  cond1array;
7051d315f7SKerry Stevens   pthread_cond_t** cond2array;
7151d315f7SKerry Stevens   void* (*pfunc)(void*);
7251d315f7SKerry Stevens   void** pdata;
7351d315f7SKerry Stevens   PetscBool startJob;
7451d315f7SKerry Stevens   estat eJobStat;
7551d315f7SKerry Stevens   PetscBool** arrThreadStarted;
7651d315f7SKerry Stevens   PetscBool** arrThreadReady;
7751d315f7SKerry Stevens } sjob_tree;
7851d315f7SKerry Stevens sjob_tree job_tree;
7951d315f7SKerry Stevens typedef struct {
8051d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8151d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8251d315f7SKerry Stevens   pthread_cond_t** cond2array;
8351d315f7SKerry Stevens   void* (*pfunc)(void*);
8451d315f7SKerry Stevens   void** pdata;
8551d315f7SKerry Stevens   PetscBool** arrThreadReady;
8651d315f7SKerry Stevens } sjob_main;
8751d315f7SKerry Stevens sjob_main job_main;
8851d315f7SKerry Stevens typedef struct {
8951d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
9051d315f7SKerry Stevens   pthread_cond_t**  cond1array;
9151d315f7SKerry Stevens   pthread_cond_t** cond2array;
9251d315f7SKerry Stevens   void* (*pfunc)(void*);
9351d315f7SKerry Stevens   void** pdata;
9451d315f7SKerry Stevens   PetscBool startJob;
9551d315f7SKerry Stevens   estat eJobStat;
9651d315f7SKerry Stevens   PetscBool** arrThreadStarted;
9751d315f7SKerry Stevens   PetscBool** arrThreadReady;
9851d315f7SKerry Stevens } sjob_chain;
9951d315f7SKerry Stevens sjob_chain job_chain;
100ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
10151dcc849SKerry Stevens typedef struct {
10251dcc849SKerry Stevens   pthread_mutex_t mutex;
10351dcc849SKerry Stevens   pthread_cond_t cond;
10451dcc849SKerry Stevens   void* (*pfunc)(void*);
10551dcc849SKerry Stevens   void** pdata;
10651dcc849SKerry Stevens   pthread_barrier_t* pbarr;
10751dcc849SKerry Stevens   int iNumJobThreads;
10851dcc849SKerry Stevens   int iNumReadyThreads;
10951dcc849SKerry Stevens   PetscBool startJob;
11051d315f7SKerry Stevens } sjob_true;
11151d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
112ba61063dSBarry Smith #endif
11351dcc849SKerry Stevens 
114ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
115ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
116ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
117ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
118ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
119ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
12051dcc849SKerry Stevens 
12151d315f7SKerry Stevens /* Function Pointers */
12251d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
12351d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
12451d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
12551d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
12651d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
12736d20dc5SKerry Stevens /**** Tree Thread Pool Functions ****/
12851d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
12951d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
13051d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
13151d315f7SKerry Stevens void           MainWait_Tree(void);
13251d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
13336d20dc5SKerry Stevens /**** Main Thread Pool Functions ****/
13451d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
13551d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
13651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
13751d315f7SKerry Stevens void           MainWait_Main(void);
13851d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
13936d20dc5SKerry Stevens /**** Chain Thread Pool Functions ****/
14051d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
14151d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
14251d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
14351d315f7SKerry Stevens void           MainWait_Chain(void);
14451d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
14536d20dc5SKerry Stevens /**** True Thread Pool Functions ****/
14651d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
14751d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
14851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
14951d315f7SKerry Stevens void           MainWait_True(void);
15051d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
15136d20dc5SKerry Stevens /**** NO Thread Pool Function  ****/
152683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void**,PetscInt);
15336d20dc5SKerry Stevens /****  ****/
15451dcc849SKerry Stevens void* FuncFinish(void*);
1550ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1560ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
157ba61063dSBarry Smith #endif
158e5c89e4eSSatish Balay 
159e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
160e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
161e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
162e5c89e4eSSatish Balay #endif
1632c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
164500d8756SSatish Balay MPI_Datatype   MPIU_C_DOUBLE_COMPLEX;
165500d8756SSatish Balay MPI_Datatype   MPIU_C_COMPLEX;
1662c876bd9SBarry Smith #endif
1677087cfbeSBarry Smith PetscScalar    PETSC_i;
168e5c89e4eSSatish Balay #else
1697087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
170e5c89e4eSSatish Balay #endif
171ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
172c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
173c90a1750SBarry Smith #endif
1747087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1757087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
17675567043SBarry Smith 
177e5c89e4eSSatish Balay /*
178e5c89e4eSSatish Balay      These are needed by petscbt.h
179e5c89e4eSSatish Balay */
180c6db04a5SJed Brown #include <petscbt.h>
1817087cfbeSBarry Smith char      _BT_mask = ' ';
1827087cfbeSBarry Smith char      _BT_c = ' ';
1837087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
184e5c89e4eSSatish Balay 
185e5c89e4eSSatish Balay /*
186e5c89e4eSSatish Balay        Function that is called to display all error messages
187e5c89e4eSSatish Balay */
1887087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1897087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
190238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1917087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
192238ccf28SShri Abhyankar #else
1937087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
194238ccf28SShri Abhyankar #endif
195bab1f7e6SVictor Minden /*
1968154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1978154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
198bab1f7e6SVictor Minden 
199e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
200e5c89e4eSSatish Balay /*
201e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
202e5c89e4eSSatish Balay */
203e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
204e5c89e4eSSatish Balay 
205e5c89e4eSSatish Balay #undef __FUNCT__
206f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
2077087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
208e5c89e4eSSatish Balay {
209e5c89e4eSSatish Balay   PetscErrorCode ierr;
210e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
211e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
212e5c89e4eSSatish Balay   char           version[256];
213e5c89e4eSSatish Balay 
214e5c89e4eSSatish Balay   PetscFunctionBegin;
215e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
216e5c89e4eSSatish Balay   if (!rank) {
217e5c89e4eSSatish Balay     char        arch[10];
218f56c2debSBarry Smith     int         err;
21988c29154SBarry Smith     PetscViewer viewer;
220f56c2debSBarry Smith 
221e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
222e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
223a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
224e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
225e5c89e4eSSatish Balay     if (filename) {
226e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
227e5c89e4eSSatish Balay     } else {
228e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
229e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
230e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
231e5c89e4eSSatish Balay     }
232e5c89e4eSSatish Balay 
233e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
234e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
235e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
236e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
237e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
23888c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
23988c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2406bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
241e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
242f56c2debSBarry Smith     err = fflush(*fd);
243e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
244e5c89e4eSSatish Balay   }
245e5c89e4eSSatish Balay   PetscFunctionReturn(0);
246e5c89e4eSSatish Balay }
247e5c89e4eSSatish Balay 
248e5c89e4eSSatish Balay #undef __FUNCT__
249f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2507087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
251e5c89e4eSSatish Balay {
252e5c89e4eSSatish Balay   PetscErrorCode ierr;
253e5c89e4eSSatish Balay   PetscMPIInt    rank;
254e5c89e4eSSatish Balay   char           date[64];
255f56c2debSBarry Smith   int            err;
256e5c89e4eSSatish Balay 
257e5c89e4eSSatish Balay   PetscFunctionBegin;
258e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
259e5c89e4eSSatish Balay   if (!rank) {
260e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
261e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
262e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
263e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
264f56c2debSBarry Smith     err = fflush(*fd);
265e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
266f56c2debSBarry Smith     err = fclose(*fd);
267e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
268e5c89e4eSSatish Balay   }
269e5c89e4eSSatish Balay   PetscFunctionReturn(0);
270e5c89e4eSSatish Balay }
271e5c89e4eSSatish Balay 
272e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
273e5c89e4eSSatish Balay 
274e5c89e4eSSatish Balay /*
275e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
276e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
277e5c89e4eSSatish Balay 
278e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2793c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2803c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
281e5c89e4eSSatish Balay */
282e5c89e4eSSatish Balay 
283e5c89e4eSSatish Balay #undef __FUNCT__
284e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
285e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
286e5c89e4eSSatish Balay {
287e5c89e4eSSatish Balay   PetscFunctionBegin;
2883c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
289e5c89e4eSSatish Balay   abort();
290e5c89e4eSSatish Balay }
291e5c89e4eSSatish Balay 
292e5c89e4eSSatish Balay #undef __FUNCT__
293e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
294e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
295e5c89e4eSSatish Balay {
296e5c89e4eSSatish Balay   PetscErrorCode ierr;
297e5c89e4eSSatish Balay 
298e5c89e4eSSatish Balay   PetscFunctionBegin;
2993c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
300e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
301e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
3023c311c98SBarry Smith     MPI_Abort(*comm,*flag);
303e5c89e4eSSatish Balay   }
304e5c89e4eSSatish Balay }
305e5c89e4eSSatish Balay 
306e5c89e4eSSatish Balay #undef __FUNCT__
307e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
308e5c89e4eSSatish Balay /*@C
309e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
310e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
311e5c89e4eSSatish Balay 
312e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
313e5c89e4eSSatish Balay 
314e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
315e5c89e4eSSatish Balay 
316e5c89e4eSSatish Balay    Level: advanced
317e5c89e4eSSatish Balay 
318e5c89e4eSSatish Balay    Note:
319e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
320e5c89e4eSSatish Balay 
32188c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
322e5c89e4eSSatish Balay @*/
3237087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
324e5c89e4eSSatish Balay {
325e5c89e4eSSatish Balay   PetscFunctionBegin;
326e5c89e4eSSatish Balay   PetscFinalize();
327e5c89e4eSSatish Balay   exit(0);
328e5c89e4eSSatish Balay   return 0;
329e5c89e4eSSatish Balay }
330e5c89e4eSSatish Balay 
331ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
33209573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
333ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
334e5c89e4eSSatish Balay static char       emacsmachinename[256];
335e5c89e4eSSatish Balay 
336e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
337e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
338e5c89e4eSSatish Balay 
339e5c89e4eSSatish Balay #undef __FUNCT__
340e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
341e5c89e4eSSatish Balay /*@C
342e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
343e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
344e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
345e5c89e4eSSatish Balay 
346e5c89e4eSSatish Balay    Input Parameter:
347e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
348da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
349e5c89e4eSSatish Balay 
350e5c89e4eSSatish Balay    Level: developer
351e5c89e4eSSatish Balay 
352e5c89e4eSSatish Balay    Concepts: package help message
353e5c89e4eSSatish Balay 
354e5c89e4eSSatish Balay @*/
3557087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
356e5c89e4eSSatish Balay {
357e5c89e4eSSatish Balay   PetscFunctionBegin;
358e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
359e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
360e5c89e4eSSatish Balay   PetscFunctionReturn(0);
361e5c89e4eSSatish Balay }
362e5c89e4eSSatish Balay 
363*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_CPU_SET_T)
364*fdfc40dbSShri Abhyankar PETSC_STATIC_INLINE PetscErrorCode PetscPthreadSetAffinity(PetscInt icorr)
365*fdfc40dbSShri Abhyankar {
366*fdfc40dbSShri Abhyankar   cpu_set_t mset;
367*fdfc40dbSShri Abhyankar   int ncorr = get_nprocs();
368*fdfc40dbSShri Abhyankar 
369*fdfc40dbSShri Abhyankar   CPU_ZERO(&mset);
370*fdfc40dbSShri Abhyankar   CPU_SET(icorr%ncorr,&mset);
371*fdfc40dbSShri Abhyankar   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
372*fdfc40dbSShri Abhyankar   return 0;
373*fdfc40dbSShri Abhyankar }
374*fdfc40dbSShri Abhyankar #endif
375*fdfc40dbSShri Abhyankar 
376*fdfc40dbSShri Abhyankar 
377e5c89e4eSSatish Balay #undef __FUNCT__
378e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3797087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
380e5c89e4eSSatish Balay {
381e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
382e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
383ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
384e5c89e4eSSatish Balay   PetscErrorCode ierr;
385a6d0e24fSJed Brown   PetscReal      si;
386e5c89e4eSSatish Balay   int            i;
387e5c89e4eSSatish Balay   PetscMPIInt    rank;
388e5c89e4eSSatish Balay   char           version[256];
389e5c89e4eSSatish Balay 
390e5c89e4eSSatish Balay   PetscFunctionBegin;
391e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
392e5c89e4eSSatish Balay 
393e5c89e4eSSatish Balay   /*
394e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
395e5c89e4eSSatish Balay   */
3968bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
39781b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
398acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
399e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
400555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
401555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
402555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
403555d055bSBarry Smith #endif
404e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
405555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
406555d055bSBarry Smith     }
407555d055bSBarry Smith #endif
408e5c89e4eSSatish Balay   }
409e5c89e4eSSatish Balay #else
410acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
411acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
412e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
413e5c89e4eSSatish Balay #endif
414e5c89e4eSSatish Balay   if (flg3) {
415e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
416e5c89e4eSSatish Balay   }
41790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
418acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
419e5c89e4eSSatish Balay   if (flg1) {
420e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
421e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
422e5c89e4eSSatish Balay   }
423e5c89e4eSSatish Balay 
42490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
425acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4267783f70dSSatish Balay   if (!flg1) {
42790d69ab7SBarry Smith     flg1 = PETSC_FALSE;
428acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4297783f70dSSatish Balay   }
430e5c89e4eSSatish Balay   if (flg1) {
431e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
432e5c89e4eSSatish Balay   }
433e5c89e4eSSatish Balay 
434e5c89e4eSSatish Balay   /*
435e5c89e4eSSatish Balay       Set the display variable for graphics
436e5c89e4eSSatish Balay   */
437e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
438e5c89e4eSSatish Balay 
43951dcc849SKerry Stevens 
44051dcc849SKerry Stevens   /*
441e5c89e4eSSatish Balay       Print the PETSc version information
442e5c89e4eSSatish Balay   */
443e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
444e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
445e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
446e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
447e5c89e4eSSatish Balay 
448e5c89e4eSSatish Balay     /*
449e5c89e4eSSatish Balay        Print "higher-level" package version message
450e5c89e4eSSatish Balay     */
451e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
452e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
453e5c89e4eSSatish Balay     }
454e5c89e4eSSatish Balay 
455a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
456e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
457e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
458e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
459e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
460e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
46184e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
462e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
463e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
464e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
465e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
466e5c89e4eSSatish Balay   }
467e5c89e4eSSatish Balay 
468e5c89e4eSSatish Balay   /*
469e5c89e4eSSatish Balay        Print "higher-level" package help message
470e5c89e4eSSatish Balay   */
471e5c89e4eSSatish Balay   if (flg3){
472e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
473e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
474e5c89e4eSSatish Balay     }
475e5c89e4eSSatish Balay   }
476e5c89e4eSSatish Balay 
477e5c89e4eSSatish Balay   /*
478e5c89e4eSSatish Balay       Setup the error handling
479e5c89e4eSSatish Balay   */
48090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
481acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
482cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
48390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
484acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
485cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
48690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
487acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
488e5c89e4eSSatish Balay   if (flg1) {
489e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
490e5c89e4eSSatish Balay   }
49190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
492acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
493cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
49496cc47afSJed Brown   flg1 = PETSC_FALSE;
495acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
49696cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
497e5c89e4eSSatish Balay 
498e5c89e4eSSatish Balay   /*
499e5c89e4eSSatish Balay       Setup debugger information
500e5c89e4eSSatish Balay   */
501e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
502e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
503e5c89e4eSSatish Balay   if (flg1) {
504e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
505e5c89e4eSSatish Balay 
506e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
507e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
508e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
509e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
510e5c89e4eSSatish Balay   }
5115e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
5125e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
513e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
514e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
515e5c89e4eSSatish Balay   if (flg1 || flg2) {
516e5c89e4eSSatish Balay     PetscMPIInt    size;
517e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
518e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
519e5c89e4eSSatish Balay     /*
520e5c89e4eSSatish Balay        we have to make sure that all processors have opened
521e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
522e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
523e5c89e4eSSatish Balay        and kill the program.
524e5c89e4eSSatish Balay     */
525e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
526e5c89e4eSSatish Balay     if (size > 2) {
527533163c2SBarry Smith       PetscMPIInt dummy = 0;
528e5c89e4eSSatish Balay       MPI_Status  status;
529e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
530e5c89e4eSSatish Balay         if (rank != i) {
531e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
532e5c89e4eSSatish Balay         }
533e5c89e4eSSatish Balay       }
534e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
535e5c89e4eSSatish Balay         if (rank != i) {
536e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
537e5c89e4eSSatish Balay         }
538e5c89e4eSSatish Balay       }
539e5c89e4eSSatish Balay     }
540e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
541e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
542e5c89e4eSSatish Balay     lsize = size;
543e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
544e5c89e4eSSatish Balay     if (flag) {
545e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
546e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
547e5c89e4eSSatish Balay       }
548e5c89e4eSSatish Balay     }
549e5c89e4eSSatish Balay     if (!flag) {
550e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
551e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
552e5c89e4eSSatish Balay       if (flg1) {
553e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
554e5c89e4eSSatish Balay       } else {
555e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
556e5c89e4eSSatish Balay       }
557e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
558e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
559e5c89e4eSSatish Balay     }
560e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
561e5c89e4eSSatish Balay   }
562e5c89e4eSSatish Balay 
563e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
564cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
565e5c89e4eSSatish Balay 
56693ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
56722b84c2fSbcordonn   /*
56822b84c2fSbcordonn     Activates new sockets for zope if needed
56922b84c2fSbcordonn   */
57084ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
571d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5726dc8fec2Sbcordonn   if (flgz){
57322b84c2fSbcordonn     int  sockfd;
574f1384234SBarry Smith     char hostname[256];
57522b84c2fSbcordonn     char username[256];
5766dc8fec2Sbcordonn     int  remoteport = 9999;
5779c4c166aSBarry Smith 
57884ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
57984ab5442Sbcordonn     if (!hostname[0]){
5809c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5819c4c166aSBarry Smith     }
58222b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5839c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
58422b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
58522b84c2fSbcordonn     if (flgzout){
58622b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
587606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5886dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5899c4c166aSBarry Smith     } else {
590d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
591d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5929c4c166aSBarry Smith     }
5939c4c166aSBarry Smith   }
59493ba235fSBarry Smith #endif
595ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
596ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
597ffc871a5SBarry Smith   if (flgz){
598ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
599ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
600ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
601ffc871a5SBarry Smith   }
602ffc871a5SBarry Smith #endif
6036dc8fec2Sbcordonn 
604e5c89e4eSSatish Balay   /*
605e5c89e4eSSatish Balay         Setup profiling and logging
606e5c89e4eSSatish Balay   */
6076cf91177SBarry Smith #if defined (PETSC_USE_INFO)
6088bb29257SSatish Balay   {
609e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
6106cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
6118bb29257SSatish Balay     if (flg1 && logname[0]) {
612fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
6138bb29257SSatish Balay     } else if (flg1) {
614fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
615e5c89e4eSSatish Balay     }
616e5c89e4eSSatish Balay   }
617865f6aa8SSatish Balay #endif
618865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
619865f6aa8SSatish Balay   mname[0] = 0;
620f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
621865f6aa8SSatish Balay   if (flg1) {
622865f6aa8SSatish Balay     if (mname[0]) {
623f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
624865f6aa8SSatish Balay     } else {
625f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
626865f6aa8SSatish Balay     }
627865f6aa8SSatish Balay   }
628e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
62990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
630fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
631e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
632e5c89e4eSSatish Balay #endif
63390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
63490d69ab7SBarry Smith   flg2 = PETSC_FALSE;
63590d69ab7SBarry Smith   flg3 = PETSC_FALSE;
636acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
637acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
638d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6399f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
640e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6419f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
642e5c89e4eSSatish Balay 
643e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
644e5c89e4eSSatish Balay   if (flg1) {
645e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
646e5c89e4eSSatish Balay     FILE *file;
647e5c89e4eSSatish Balay     if (mname[0]) {
648e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
649e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
650e5c89e4eSSatish Balay       file = fopen(fname,"w");
651f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
652e5c89e4eSSatish Balay     } else {
653da9f1d6bSBarry Smith       file = PETSC_STDOUT;
654e5c89e4eSSatish Balay     }
655e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
656e5c89e4eSSatish Balay   }
657e5c89e4eSSatish Balay #endif
658e5c89e4eSSatish Balay 
659e5c89e4eSSatish Balay   /*
660e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
661e5c89e4eSSatish Balay   */
66263d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
663e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
664e5c89e4eSSatish Balay #endif
665e5c89e4eSSatish Balay 
666acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
667e5c89e4eSSatish Balay 
668ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
669af359df3SBarry Smith   /*
670af359df3SBarry Smith       Determine whether user specified maximum number of threads
671af359df3SBarry Smith    */
672af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
673af359df3SBarry Smith 
674b154f58aSKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-main",&flg1);CHKERRQ(ierr);
675b154f58aSKerry Stevens   if(flg1) {
676*fdfc40dbSShri Abhyankar     PetscInt icorr;
677b154f58aSKerry Stevens     ierr = PetscOptionsGetInt(PETSC_NULL,"-main",&icorr,PETSC_NULL);CHKERRQ(ierr);
678*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_CPU_SET_T)
679*fdfc40dbSShri Abhyankar     ierr = PetscPthreadSetAffinity(icorr);CHKERRQ(ierr);
680*fdfc40dbSShri Abhyankar #endif
681b154f58aSKerry Stevens   }
682b154f58aSKerry Stevens 
683*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_CPU_SET_T)
684af359df3SBarry Smith   PetscInt N_CORES = get_nprocs();
685af359df3SBarry Smith   ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
686af359df3SBarry Smith   char tstr[9];
687af359df3SBarry Smith   char tbuf[2];
688af359df3SBarry Smith   strcpy(tstr,"-thread");
689af359df3SBarry Smith   for(i=0;i<PetscMaxThreads;i++) {
690af359df3SBarry Smith     ThreadCoreAffinity[i] = i;
691af359df3SBarry Smith     sprintf(tbuf,"%d",i);
692af359df3SBarry Smith     strcat(tstr,tbuf);
693af359df3SBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
694af359df3SBarry Smith     if(flg1) {
695af359df3SBarry Smith       ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
696af359df3SBarry Smith       ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
697af359df3SBarry Smith     }
698af359df3SBarry Smith     tstr[7] = '\0';
699af359df3SBarry Smith   }
700*fdfc40dbSShri Abhyankar #endif
701cfcfc605SKerry Stevens 
702cfcfc605SKerry Stevens   /*
703cfcfc605SKerry Stevens       Determine whether to use thread pool
704cfcfc605SKerry Stevens    */
705cfcfc605SKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
706cfcfc605SKerry Stevens   if (flg1) {
707cfcfc605SKerry Stevens     PetscUseThreadPool = PETSC_TRUE;
708af359df3SBarry Smith     /* get the thread pool type */
709af359df3SBarry Smith     PetscInt ipool = 0;
710af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
711af359df3SBarry Smith 
712af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
713af359df3SBarry Smith     switch(ipool) {
714af359df3SBarry Smith     case 1:
715af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
716af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
717af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
718af359df3SBarry Smith       MainWait              = &MainWait_Tree;
719af359df3SBarry Smith       MainJob               = &MainJob_Tree;
720af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
721af359df3SBarry Smith       break;
722af359df3SBarry Smith     case 2:
723af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
724af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
725af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
726af359df3SBarry Smith       MainWait              = &MainWait_Main;
727af359df3SBarry Smith       MainJob               = &MainJob_Main;
728af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
729af359df3SBarry Smith       break;
730af359df3SBarry Smith     case 3:
731af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
732af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
733af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
734af359df3SBarry Smith       MainWait              = &MainWait_Chain;
735af359df3SBarry Smith       MainJob               = &MainJob_Chain;
736af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
737af359df3SBarry Smith       break;
738af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
739af359df3SBarry Smith     default:
740af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
741af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
742af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
743af359df3SBarry Smith       MainWait              = &MainWait_True;
744af359df3SBarry Smith       MainJob               = &MainJob_True;
745af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
746af359df3SBarry Smith       break;
747*fdfc40dbSShri Abhyankar # else
748*fdfc40dbSShri Abhyankar     default:
749*fdfc40dbSShri Abhyankar       PetscThreadFunc       = &PetscThreadFunc_Chain;
750*fdfc40dbSShri Abhyankar       PetscThreadInitialize = &PetscThreadInitialize_Chain;
751*fdfc40dbSShri Abhyankar       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
752*fdfc40dbSShri Abhyankar       MainWait              = &MainWait_Chain;
753*fdfc40dbSShri Abhyankar       MainJob               = &MainJob_Chain;
754*fdfc40dbSShri Abhyankar       PetscInfo(PETSC_NULL,"Cannot use true thread pool since pthread_barrier_t is not available,using chain thread pool instead\n");
755*fdfc40dbSShri Abhyankar       break;
756af359df3SBarry Smith #endif
757af359df3SBarry Smith     }
758af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
759683509dcSKerry Stevens   } else {
760683509dcSKerry Stevens     //need to define these in the case on 'no threads' or 'thread create/destroy'
761683509dcSKerry Stevens     //could take any of the above versions
762683509dcSKerry Stevens     MainJob               = &MainJob_Spawn;
763af359df3SBarry Smith   }
764af359df3SBarry Smith #endif
765e5c89e4eSSatish Balay   /*
766e5c89e4eSSatish Balay        Print basic help message
767e5c89e4eSSatish Balay   */
768e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
769e5c89e4eSSatish Balay   if (flg1) {
770e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
771301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
772301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
773301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
774e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
775e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
776e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
777e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
778e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
779e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
780e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
781e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
782e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
783e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
784e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
785e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
786e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
787e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
788e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
789e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
790e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
791e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
792e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7934161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7944161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
795e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
796e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
797e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
798e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
799e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
800e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
801a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
802e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
80340ab9619SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -server <port>: Run PETSc webserver (default port is 8080) see PetscWebServe()\n");CHKERRQ(ierr);
804e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
805e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
8067ef452c0SMatthew G Knepley     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary _summary_python]: logging objects and events\n");CHKERRQ(ierr);
807e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
808e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
809e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
810e5c89e4eSSatish Balay #endif
8116cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
812e5c89e4eSSatish Balay #endif
813e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
814e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
815e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
816e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
817e5c89e4eSSatish Balay   }
818e5c89e4eSSatish Balay 
819a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
820e5c89e4eSSatish Balay   if (flg1) {
821e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
822e5c89e4eSSatish Balay   }
823e5c89e4eSSatish Balay 
8246cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
825e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
826e5c89e4eSSatish Balay   if (f) {
8276cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
828e5c89e4eSSatish Balay   }
829827f890bSBarry Smith 
8308154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
831c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
83273113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
83373113deaSBarry Smith   else flg1 = PETSC_FALSE;
8348154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
8358154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
836bab1f7e6SVictor Minden #endif
837192daf7cSBarry Smith 
838e5c89e4eSSatish Balay   PetscFunctionReturn(0);
839e5c89e4eSSatish Balay }
840df413903SBarry Smith 
841ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
842ba61063dSBarry Smith 
84351d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
84451d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
84551d315f7SKerry Stevens   PetscErrorCode iterr;
846*fdfc40dbSShri Abhyankar   int ierr;
84751d315f7SKerry Stevens   int* pId = (int*)arg;
84851d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
84951d315f7SKerry Stevens   PetscBool PeeOn;
850*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_CPU_SET_T)
851*fdfc40dbSShri Abhyankar   iterr = PetscPthreadSetAffinity(ThreadCoreAffinity[ThreadId]);CHKERRQ(iterr);
852*fdfc40dbSShri Abhyankar #endif
8539e800a48SKerry Stevens   //printf("Thread %d In Tree Thread Function\n",ThreadId);
85451d315f7SKerry Stevens 
85551d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
85651d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
85751d315f7SKerry Stevens   }
85851d315f7SKerry Stevens   else {
85951d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
86051d315f7SKerry Stevens   }
86151d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
862ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
86351d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
86451d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
86551d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
86651d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
86751d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
868ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
869ba61063dSBarry Smith            upon return, has the lock */
87051d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
87151d315f7SKerry Stevens         }
87251d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
87351d315f7SKerry Stevens       }
87451d315f7SKerry Stevens     }
875ba61063dSBarry Smith     /* your subordinates are now ready */
87651d315f7SKerry Stevens   }
87751d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
878ba61063dSBarry Smith   /* update your ready status */
87951d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
88051d315f7SKerry Stevens   if(ThreadId==0) {
88151d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
882ba61063dSBarry Smith     /* ignal main */
88351d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
88451d315f7SKerry Stevens   }
88551d315f7SKerry Stevens   else {
886ba61063dSBarry Smith     /* tell your boss that you're ready to work */
88751d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
88851d315f7SKerry Stevens   }
889ba61063dSBarry Smith   /* the while loop needs to have an exit
890ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
891ba61063dSBarry Smith    and calling FuncFinish */
89251d315f7SKerry Stevens   while(PetscThreadGo) {
893ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
894ba61063dSBarry Smith       waiting when you don't have to causes problems
895ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
89651d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
897ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
898ba61063dSBarry Smith        upon return, has the lock */
89951d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
90051d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
90151d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
90251d315f7SKerry Stevens     }
90351d315f7SKerry Stevens     if(ThreadId==0) {
90451d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
90551d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
90651d315f7SKerry Stevens     }
90751d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
90851d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
909ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
91051d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
91151d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
91251d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
91351d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
91451d315f7SKerry Stevens         }
91551d315f7SKerry Stevens       }
91651d315f7SKerry Stevens     }
917ba61063dSBarry Smith     /* do your job */
91851d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
91951d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
92051d315f7SKerry Stevens     }
92151d315f7SKerry Stevens     else {
92251d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
92351d315f7SKerry Stevens     }
92451d315f7SKerry Stevens     if(iterr!=0) {
92551d315f7SKerry Stevens       ithreaderr = 1;
92651d315f7SKerry Stevens     }
92751d315f7SKerry Stevens     if(PetscThreadGo) {
928ba61063dSBarry Smith       /* reset job, get ready for more */
92951d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
930ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
931ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
93251d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
93351d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
93451d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
93551d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
93651d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
937ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
938ba61063dSBarry Smith                upon return, has the lock */
93951d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
94051d315f7SKerry Stevens             }
94151d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
94251d315f7SKerry Stevens           }
94351d315f7SKerry Stevens 	}
944ba61063dSBarry Smith         /* your subordinates are now ready */
94551d315f7SKerry Stevens       }
94651d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
94751d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
94851d315f7SKerry Stevens       if(ThreadId==0) {
949ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
950ba61063dSBarry Smith         /* root thread signals 'main' */
95151d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
95251d315f7SKerry Stevens       }
95351d315f7SKerry Stevens       else {
954ba61063dSBarry Smith         /* signal your boss before you go to sleep */
95551d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
95651d315f7SKerry Stevens       }
95751d315f7SKerry Stevens     }
95851d315f7SKerry Stevens   }
95951d315f7SKerry Stevens   return NULL;
96051d315f7SKerry Stevens }
96151d315f7SKerry Stevens 
96251d315f7SKerry Stevens #undef __FUNCT__
96351d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
96451d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
96551d315f7SKerry Stevens   PetscInt i,ierr;
96651d315f7SKerry Stevens   int status;
96751d315f7SKerry Stevens 
96851d315f7SKerry Stevens   if(PetscUseThreadPool) {
969*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_MEMALIGN)
97051d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
97151d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
97251d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
97351d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
97451d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
97551d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
97651d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
977*fdfc40dbSShri Abhyankar #else
978*fdfc40dbSShri Abhyankar     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
979*fdfc40dbSShri Abhyankar     arrmutex = (char*)malloc(Val2);
980*fdfc40dbSShri Abhyankar     arrcond1 = (char*)malloc(Val2);
981*fdfc40dbSShri Abhyankar     arrcond2 = (char*)malloc(Val2);
982*fdfc40dbSShri Abhyankar     arrstart = (char*)malloc(Val2);
983*fdfc40dbSShri Abhyankar     arrready = (char*)malloc(Val2);
984*fdfc40dbSShri Abhyankar #endif
98551d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
98651d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
98751d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
98851d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
98951d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
990ba61063dSBarry Smith     /* initialize job structure */
99151d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
99251d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
99351d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
99451d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
99551d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
99651d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
99751d315f7SKerry Stevens     }
99851d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
99951d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
100051d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
100151d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
100251d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
100351d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
100451d315f7SKerry Stevens     }
100551d315f7SKerry Stevens     job_tree.pfunc = NULL;
100651d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
100751d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
100851d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
100951d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1010ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
101151d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1012ba61063dSBarry Smith     /* create threads */
101351d315f7SKerry Stevens     for(i=0; i<N; i++) {
101451d315f7SKerry Stevens       pVal[i] = i;
101551d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1016ba61063dSBarry Smith       /* should check status */
101751d315f7SKerry Stevens     }
101851d315f7SKerry Stevens   }
101951d315f7SKerry Stevens   return NULL;
102051d315f7SKerry Stevens }
102151d315f7SKerry Stevens 
102251d315f7SKerry Stevens #undef __FUNCT__
102351d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
102451d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
102551d315f7SKerry Stevens   int i,ierr;
102651d315f7SKerry Stevens   void* jstatus;
102751d315f7SKerry Stevens 
102851d315f7SKerry Stevens   PetscFunctionBegin;
102951d315f7SKerry Stevens 
1030ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1031ba61063dSBarry Smith   /* join the threads */
103251d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
103351d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1034ba61063dSBarry Smith     /* do error checking*/
103551d315f7SKerry Stevens   }
103651d315f7SKerry Stevens   free(PetscThreadPoint);
103751d315f7SKerry Stevens   free(arrmutex);
103851d315f7SKerry Stevens   free(arrcond1);
103951d315f7SKerry Stevens   free(arrcond2);
104051d315f7SKerry Stevens   free(arrstart);
104151d315f7SKerry Stevens   free(arrready);
104251d315f7SKerry Stevens   free(job_tree.pdata);
104351d315f7SKerry Stevens   free(pVal);
1044cfcfc605SKerry Stevens 
104551d315f7SKerry Stevens   PetscFunctionReturn(0);
104651d315f7SKerry Stevens }
104751d315f7SKerry Stevens 
104851d315f7SKerry Stevens #undef __FUNCT__
104951d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
105051d315f7SKerry Stevens void MainWait_Tree() {
105151d315f7SKerry Stevens   int ierr;
105251d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
105351d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
105451d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
105551d315f7SKerry Stevens   }
105651d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
105751d315f7SKerry Stevens }
105851d315f7SKerry Stevens 
105951d315f7SKerry Stevens #undef __FUNCT__
106051d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
106151d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
106251d315f7SKerry Stevens   int i,ierr;
106351d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
106436d20dc5SKerry Stevens 
106551d315f7SKerry Stevens   MainWait();
106651d315f7SKerry Stevens   job_tree.pfunc = pFunc;
106751d315f7SKerry Stevens   job_tree.pdata = data;
106851d315f7SKerry Stevens   job_tree.startJob = PETSC_TRUE;
106951d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
107051d315f7SKerry Stevens     *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
107151d315f7SKerry Stevens   }
107251d315f7SKerry Stevens   job_tree.eJobStat = JobInitiated;
107351d315f7SKerry Stevens   ierr = pthread_cond_signal(job_tree.cond2array[0]);
107451d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1075ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
107651d315f7SKerry Stevens   }
107736d20dc5SKerry Stevens 
107851d315f7SKerry Stevens   if(ithreaderr) {
107951d315f7SKerry Stevens     ijoberr = ithreaderr;
108051d315f7SKerry Stevens   }
108151d315f7SKerry Stevens   return ijoberr;
108251d315f7SKerry Stevens }
108351d315f7SKerry Stevens /****  ****/
108451d315f7SKerry Stevens 
108551d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
108651d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
108751d315f7SKerry Stevens   PetscErrorCode iterr;
1088*fdfc40dbSShri Abhyankar   int ierr;
108951d315f7SKerry Stevens   int* pId = (int*)arg;
109051d315f7SKerry Stevens   int ThreadId = *pId;
10919e800a48SKerry Stevens   //printf("Thread %d In Main Thread Function\n",ThreadId);
1092*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_CPU_SET_T)
1093*fdfc40dbSShri Abhyankar     iterr = PetscPthreadSetAffinity(ThreadCoreAffinity[ThreadId]);CHKERRQ(iterr);
1094*fdfc40dbSShri Abhyankar #endif
109551d315f7SKerry Stevens 
109651d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1097ba61063dSBarry Smith   /* update your ready status */
109851d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1099ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
110051d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
110151d315f7SKerry Stevens 
1102ba61063dSBarry Smith   /* the while loop needs to have an exit
1103ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1104ba61063dSBarry Smith      and calling FuncFinish */
110551d315f7SKerry Stevens   while(PetscThreadGo) {
1106ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1107ba61063dSBarry Smith        waiting when you don't have to causes problems
1108ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
110951d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1110ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1111ba61063dSBarry Smith        upon return, has the lock */
111251d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1113ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
111451d315f7SKerry Stevens     }
111551d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
111651d315f7SKerry Stevens     if(job_main.pdata==NULL) {
111751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
111851d315f7SKerry Stevens     }
111951d315f7SKerry Stevens     else {
112051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
112151d315f7SKerry Stevens     }
112251d315f7SKerry Stevens     if(iterr!=0) {
112351d315f7SKerry Stevens       ithreaderr = 1;
112451d315f7SKerry Stevens     }
112551d315f7SKerry Stevens     if(PetscThreadGo) {
1126ba61063dSBarry Smith       /* reset job, get ready for more */
112751d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
112851d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1129ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
113051d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
113151d315f7SKerry Stevens     }
113251d315f7SKerry Stevens   }
113351d315f7SKerry Stevens   return NULL;
113451d315f7SKerry Stevens }
113551d315f7SKerry Stevens 
113651d315f7SKerry Stevens #undef __FUNCT__
113751d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
113851d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
113951d315f7SKerry Stevens   PetscInt i,ierr;
114051d315f7SKerry Stevens   int status;
114151d315f7SKerry Stevens 
114251d315f7SKerry Stevens   if(PetscUseThreadPool) {
1143*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_MEMALIGN)
114451d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
114551d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
114651d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
114751d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
114851d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
114951d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
115051d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
1151*fdfc40dbSShri Abhyankar #else
1152*fdfc40dbSShri Abhyankar     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
1153*fdfc40dbSShri Abhyankar     arrmutex = (char*)malloc(Val2);
1154*fdfc40dbSShri Abhyankar     arrcond1 = (char*)malloc(Val2);
1155*fdfc40dbSShri Abhyankar     arrcond2 = (char*)malloc(Val2);
1156*fdfc40dbSShri Abhyankar     arrstart = (char*)malloc(Val2);
1157*fdfc40dbSShri Abhyankar     arrready = (char*)malloc(Val2);
1158*fdfc40dbSShri Abhyankar #endif
1159*fdfc40dbSShri Abhyankar 
116051d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
116151d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
116251d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
116351d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1164ba61063dSBarry Smith     /* initialize job structure */
116551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
116651d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
116751d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
116851d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
116951d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
117051d315f7SKerry Stevens     }
117151d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
117251d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
117351d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
117451d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
117551d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
117651d315f7SKerry Stevens     }
117751d315f7SKerry Stevens     job_main.pfunc = NULL;
117851d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
117951d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1180ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
118151d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1182ba61063dSBarry Smith     /* create threads */
118351d315f7SKerry Stevens     for(i=0; i<N; i++) {
118451d315f7SKerry Stevens       pVal[i] = i;
118551d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1186ba61063dSBarry Smith       /* error check */
118751d315f7SKerry Stevens     }
118851d315f7SKerry Stevens   }
118951d315f7SKerry Stevens   else {
119051d315f7SKerry Stevens   }
119151d315f7SKerry Stevens   return NULL;
119251d315f7SKerry Stevens }
119351d315f7SKerry Stevens 
119451d315f7SKerry Stevens #undef __FUNCT__
119551d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
119651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
119751d315f7SKerry Stevens   int i,ierr;
119851d315f7SKerry Stevens   void* jstatus;
119951d315f7SKerry Stevens 
120051d315f7SKerry Stevens   PetscFunctionBegin;
120151d315f7SKerry Stevens 
1202ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1203ba61063dSBarry Smith   /* join the threads */
120451d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1205ba61063dSBarry Smith     ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
120651d315f7SKerry Stevens   }
120751d315f7SKerry Stevens   free(PetscThreadPoint);
120851d315f7SKerry Stevens   free(arrmutex);
120951d315f7SKerry Stevens   free(arrcond1);
121051d315f7SKerry Stevens   free(arrcond2);
121151d315f7SKerry Stevens   free(arrstart);
121251d315f7SKerry Stevens   free(arrready);
121351d315f7SKerry Stevens   free(job_main.pdata);
121451d315f7SKerry Stevens   free(pVal);
1215cfcfc605SKerry Stevens 
121651d315f7SKerry Stevens   PetscFunctionReturn(0);
121751d315f7SKerry Stevens }
121851d315f7SKerry Stevens 
121951d315f7SKerry Stevens #undef __FUNCT__
122051d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
122151d315f7SKerry Stevens void MainWait_Main() {
122251d315f7SKerry Stevens   int i,ierr;
122351d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
122451d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
122551d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
122651d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
122751d315f7SKerry Stevens     }
122851d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
122951d315f7SKerry Stevens   }
123051d315f7SKerry Stevens }
123151d315f7SKerry Stevens 
123251d315f7SKerry Stevens #undef __FUNCT__
123351d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
123451d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
123551d315f7SKerry Stevens   int i,ierr;
123651d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
123736d20dc5SKerry Stevens 
1238ba61063dSBarry Smith   MainWait(); /* you know everyone is waiting to be signalled! */
123951d315f7SKerry Stevens   job_main.pfunc = pFunc;
124051d315f7SKerry Stevens   job_main.pdata = data;
124151d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1242ba61063dSBarry Smith     *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
124351d315f7SKerry Stevens   }
1244ba61063dSBarry Smith   /* tell the threads to go to work */
124551d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
124651d315f7SKerry Stevens     ierr = pthread_cond_signal(job_main.cond2array[i]);
124751d315f7SKerry Stevens   }
124851d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1249ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
125051d315f7SKerry Stevens   }
125136d20dc5SKerry Stevens 
125251d315f7SKerry Stevens   if(ithreaderr) {
125351d315f7SKerry Stevens     ijoberr = ithreaderr;
125451d315f7SKerry Stevens   }
125551d315f7SKerry Stevens   return ijoberr;
125651d315f7SKerry Stevens }
125751d315f7SKerry Stevens /****  ****/
125851d315f7SKerry Stevens 
125951d315f7SKerry Stevens /**** Chain Thread Functions ****/
126051d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
126151d315f7SKerry Stevens   PetscErrorCode iterr;
1262*fdfc40dbSShri Abhyankar   int ierr;
126351d315f7SKerry Stevens   int* pId = (int*)arg;
126451d315f7SKerry Stevens   int ThreadId = *pId;
126551d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
126651d315f7SKerry Stevens   PetscBool PeeOn;
12679e800a48SKerry Stevens   //printf("Thread %d In Chain Thread Function\n",ThreadId);
1268*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_CPU_SET_T)
1269*fdfc40dbSShri Abhyankar   iterr = PetscPthreadSetAffinity(ThreadCoreAffinity[ThreadId]);CHKERRQ(iterr);
1270*fdfc40dbSShri Abhyankar #endif
127151d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
127251d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
127351d315f7SKerry Stevens   }
127451d315f7SKerry Stevens   else {
127551d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
127651d315f7SKerry Stevens   }
127751d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1278ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
127951d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
128051d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1281ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1282ba61063dSBarry Smith        upon return, has the lock */
128351d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
128451d315f7SKerry Stevens     }
128551d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1286ba61063dSBarry Smith     /* your subordinate is now ready*/
128751d315f7SKerry Stevens   }
128851d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1289ba61063dSBarry Smith   /* update your ready status */
129051d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
129151d315f7SKerry Stevens   if(ThreadId==0) {
129251d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1293ba61063dSBarry Smith     /* signal main */
129451d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
129551d315f7SKerry Stevens   }
129651d315f7SKerry Stevens   else {
1297ba61063dSBarry Smith     /* tell your boss that you're ready to work */
129851d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
129951d315f7SKerry Stevens   }
1300ba61063dSBarry Smith   /*  the while loop needs to have an exit
1301ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1302ba61063dSBarry Smith    and calling FuncFinish */
130351d315f7SKerry Stevens   while(PetscThreadGo) {
1304ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1305ba61063dSBarry Smith        waiting when you don't have to causes problems
1306ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
130751d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1308ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1309ba61063dSBarry Smith        upon return, has the lock */
131051d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
131151d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
131251d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
131351d315f7SKerry Stevens     }
131451d315f7SKerry Stevens     if(ThreadId==0) {
131551d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
131651d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
131751d315f7SKerry Stevens     }
131851d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
131951d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1320ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
132151d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
132251d315f7SKerry Stevens     }
1323ba61063dSBarry Smith     /* do your job */
132451d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
132551d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
132651d315f7SKerry Stevens     }
132751d315f7SKerry Stevens     else {
132851d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
132951d315f7SKerry Stevens     }
133051d315f7SKerry Stevens     if(iterr!=0) {
133151d315f7SKerry Stevens       ithreaderr = 1;
133251d315f7SKerry Stevens     }
133351d315f7SKerry Stevens     if(PetscThreadGo) {
1334ba61063dSBarry Smith       /* reset job, get ready for more */
133551d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1336ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1337ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
133851d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
133951d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1340ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1341ba61063dSBarry Smith            upon return, has the lock */
134251d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
134351d315f7SKerry Stevens         }
134451d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1345ba61063dSBarry Smith         /* your subordinate is now ready */
134651d315f7SKerry Stevens       }
134751d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
134851d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
134951d315f7SKerry Stevens       if(ThreadId==0) {
1350ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1351ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
135251d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
135351d315f7SKerry Stevens       }
135451d315f7SKerry Stevens       else {
1355ba61063dSBarry Smith         /* signal your boss before you go to sleep */
135651d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
135751d315f7SKerry Stevens       }
135851d315f7SKerry Stevens     }
135951d315f7SKerry Stevens   }
136051d315f7SKerry Stevens   return NULL;
136151d315f7SKerry Stevens }
136251d315f7SKerry Stevens 
136351d315f7SKerry Stevens #undef __FUNCT__
136451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
136551d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
136651d315f7SKerry Stevens   PetscInt i,ierr;
136751d315f7SKerry Stevens   int status;
136851d315f7SKerry Stevens 
136951d315f7SKerry Stevens   if(PetscUseThreadPool) {
1370*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_MEMALIGN)
137151d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
137251d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
137351d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
137451d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
137551d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
137651d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
137751d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
1378*fdfc40dbSShri Abhyankar #else
1379*fdfc40dbSShri Abhyankar     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
1380*fdfc40dbSShri Abhyankar     arrmutex = (char*)malloc(Val2);
1381*fdfc40dbSShri Abhyankar     arrcond1 = (char*)malloc(Val2);
1382*fdfc40dbSShri Abhyankar     arrcond2 = (char*)malloc(Val2);
1383*fdfc40dbSShri Abhyankar     arrstart = (char*)malloc(Val2);
1384*fdfc40dbSShri Abhyankar     arrready = (char*)malloc(Val2);
1385*fdfc40dbSShri Abhyankar #endif
1386*fdfc40dbSShri Abhyankar 
138751d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
138851d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
138951d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
139051d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
139151d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1392ba61063dSBarry Smith     /* initialize job structure */
139351d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
139451d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
139551d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
139651d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
139751d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
139851d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
139951d315f7SKerry Stevens     }
140051d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
140151d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
140251d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
140351d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
140451d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
140551d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
140651d315f7SKerry Stevens     }
140751d315f7SKerry Stevens     job_chain.pfunc = NULL;
140851d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
140951d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
141051d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
141151d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1412ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
141351d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1414ba61063dSBarry Smith     /* create threads */
141551d315f7SKerry Stevens     for(i=0; i<N; i++) {
141651d315f7SKerry Stevens       pVal[i] = i;
141751d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1418ba61063dSBarry Smith       /* should check error */
141951d315f7SKerry Stevens     }
142051d315f7SKerry Stevens   }
142151d315f7SKerry Stevens   else {
142251d315f7SKerry Stevens   }
142351d315f7SKerry Stevens   return NULL;
142451d315f7SKerry Stevens }
142551d315f7SKerry Stevens 
142651d315f7SKerry Stevens 
142751d315f7SKerry Stevens #undef __FUNCT__
142851d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
142951d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
143051d315f7SKerry Stevens   int i,ierr;
143151d315f7SKerry Stevens   void* jstatus;
143251d315f7SKerry Stevens 
143351d315f7SKerry Stevens   PetscFunctionBegin;
143451d315f7SKerry Stevens 
1435ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1436ba61063dSBarry Smith   /* join the threads */
143751d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
143851d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1439ba61063dSBarry Smith     /* should check error */
144051d315f7SKerry Stevens   }
144151d315f7SKerry Stevens   free(PetscThreadPoint);
144251d315f7SKerry Stevens   free(arrmutex);
144351d315f7SKerry Stevens   free(arrcond1);
144451d315f7SKerry Stevens   free(arrcond2);
144551d315f7SKerry Stevens   free(arrstart);
144651d315f7SKerry Stevens   free(arrready);
144751d315f7SKerry Stevens   free(job_chain.pdata);
144851d315f7SKerry Stevens   free(pVal);
1449cfcfc605SKerry Stevens 
145051d315f7SKerry Stevens   PetscFunctionReturn(0);
145151d315f7SKerry Stevens }
145251d315f7SKerry Stevens 
145351d315f7SKerry Stevens #undef __FUNCT__
145451d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
145551d315f7SKerry Stevens void MainWait_Chain() {
145651d315f7SKerry Stevens   int ierr;
145751d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
145851d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
145951d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
146051d315f7SKerry Stevens   }
146151d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
146251d315f7SKerry Stevens }
146351d315f7SKerry Stevens 
146451d315f7SKerry Stevens #undef __FUNCT__
146551d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
146651d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
146751d315f7SKerry Stevens   int i,ierr;
146851d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
146936d20dc5SKerry Stevens 
147051d315f7SKerry Stevens   MainWait();
147151d315f7SKerry Stevens   job_chain.pfunc = pFunc;
147251d315f7SKerry Stevens   job_chain.pdata = data;
147351d315f7SKerry Stevens   job_chain.startJob = PETSC_TRUE;
147451d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
147551d315f7SKerry Stevens     *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
147651d315f7SKerry Stevens   }
147751d315f7SKerry Stevens   job_chain.eJobStat = JobInitiated;
147851d315f7SKerry Stevens   ierr = pthread_cond_signal(job_chain.cond2array[0]);
147951d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1480ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
148151d315f7SKerry Stevens   }
148236d20dc5SKerry Stevens 
148351d315f7SKerry Stevens   if(ithreaderr) {
148451d315f7SKerry Stevens     ijoberr = ithreaderr;
148551d315f7SKerry Stevens   }
148651d315f7SKerry Stevens   return ijoberr;
148751d315f7SKerry Stevens }
148851d315f7SKerry Stevens /****  ****/
148951d315f7SKerry Stevens 
1490ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
149151d315f7SKerry Stevens /**** True Thread Functions ****/
149251d315f7SKerry Stevens void* PetscThreadFunc_True(void* arg) {
1493*fdfc40dbSShri Abhyankar   int ierr,iVal;
149451dcc849SKerry Stevens   int* pId = (int*)arg;
149551dcc849SKerry Stevens   int ThreadId = *pId;
14960ca81413SKerry Stevens   PetscErrorCode iterr;
14979e800a48SKerry Stevens   //printf("Thread %d In True Pool Thread Function\n",ThreadId);
1498*fdfc40dbSShri Abhyankar #if defined(PETSC_HAVE_CPU_SET_T)
1499*fdfc40dbSShri Abhyankar   iterr = PetscPthreadSetAffinity(ThreadCoreAffinity[ThreadId]);CHKERRQ(iterr);
1500*fdfc40dbSShri Abhyankar #endif
150151d315f7SKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
150251d315f7SKerry Stevens   job_true.iNumReadyThreads++;
150351d315f7SKerry Stevens   if(job_true.iNumReadyThreads==PetscMaxThreads) {
150451dcc849SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
150551dcc849SKerry Stevens   }
1506ba61063dSBarry Smith   /*the while loop needs to have an exit
1507ba61063dSBarry Smith     the 'main' thread can terminate all the threads by performing a broadcast
1508ba61063dSBarry Smith    and calling FuncFinish */
150951dcc849SKerry Stevens   while(PetscThreadGo) {
1510ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
1511ba61063dSBarry Smith       waiting when you don't have to causes problems
1512ba61063dSBarry Smith      also need to wait if another thread sneaks in and messes with the predicate */
151351d315f7SKerry Stevens     while(job_true.startJob==PETSC_FALSE&&job_true.iNumJobThreads==0) {
1514ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1515ba61063dSBarry Smith        upon return, has the lock */
15166c9b208dSKerry Stevens       //printf("Thread %d Going to Sleep!\n",ThreadId);
151751d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
151851dcc849SKerry Stevens     }
151951d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
152051d315f7SKerry Stevens     job_true.iNumJobThreads--;
152151d315f7SKerry Stevens     job_true.iNumReadyThreads--;
152251d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
152351d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
152451d315f7SKerry Stevens     if(job_true.pdata==NULL) {
152551d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
152651dcc849SKerry Stevens     }
152751dcc849SKerry Stevens     else {
152851d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
152951dcc849SKerry Stevens     }
15300ca81413SKerry Stevens     if(iterr!=0) {
15310ca81413SKerry Stevens       ithreaderr = 1;
15320ca81413SKerry Stevens     }
15336c9b208dSKerry Stevens     //printf("Thread %d Finished Job\n",ThreadId);
1534ba61063dSBarry Smith     /* the barrier is necessary BECAUSE: look at job_true.iNumReadyThreads
1535ba61063dSBarry Smith       what happens if a thread finishes before they all start? BAD!
1536ba61063dSBarry Smith      what happens if a thread finishes before any else start? BAD! */
1537ba61063dSBarry Smith     pthread_barrier_wait(job_true.pbarr); /* ensures all threads are finished */
1538ba61063dSBarry Smith     /* reset job */
153951dcc849SKerry Stevens     if(PetscThreadGo) {
154051d315f7SKerry Stevens       pthread_mutex_lock(&job_true.mutex);
154151d315f7SKerry Stevens       job_true.iNumReadyThreads++;
154251d315f7SKerry Stevens       if(job_true.iNumReadyThreads==PetscMaxThreads) {
1543ba61063dSBarry Smith 	/* signal the 'main' thread that the job is done! (only done once) */
154451dcc849SKerry Stevens 	ierr = pthread_cond_signal(&main_cond);
154551dcc849SKerry Stevens       }
154651dcc849SKerry Stevens     }
154751dcc849SKerry Stevens   }
154851dcc849SKerry Stevens   return NULL;
154951dcc849SKerry Stevens }
155051dcc849SKerry Stevens 
1551f09cb4aaSKerry Stevens #undef __FUNCT__
155251d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_True"
155351d315f7SKerry Stevens void* PetscThreadInitialize_True(PetscInt N) {
155451dcc849SKerry Stevens   PetscInt i;
155551dcc849SKerry Stevens   int status;
15560ca81413SKerry Stevens 
1557f09cb4aaSKerry Stevens   pVal = (int*)malloc(N*sizeof(int));
1558ba61063dSBarry Smith   /* allocate memory in the heap for the thread structure */
155951dcc849SKerry Stevens   PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1560ba61063dSBarry Smith   BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
156151d315f7SKerry Stevens   job_true.pdata = (void**)malloc(N*sizeof(void*));
156251dcc849SKerry Stevens   for(i=0; i<N; i++) {
1563f09cb4aaSKerry Stevens     pVal[i] = i;
1564f09cb4aaSKerry Stevens     status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1565ba61063dSBarry Smith     /* error check to ensure proper thread creation */
156651dcc849SKerry Stevens     status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1567ba61063dSBarry Smith     /* should check error */
156851dcc849SKerry Stevens   }
15696c9b208dSKerry Stevens   //printf("Finished True Thread Pool Initialization\n");
157051dcc849SKerry Stevens   return NULL;
157151dcc849SKerry Stevens }
157251dcc849SKerry Stevens 
1573f09cb4aaSKerry Stevens 
1574f09cb4aaSKerry Stevens #undef __FUNCT__
157551d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
157651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
157751dcc849SKerry Stevens   int i,ierr;
157851dcc849SKerry Stevens   void* jstatus;
157951dcc849SKerry Stevens 
158051dcc849SKerry Stevens   PetscFunctionBegin;
1581cfcfc605SKerry Stevens 
1582ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1583ba61063dSBarry Smith   /* join the threads */
158451dcc849SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
158551dcc849SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
158651dcc849SKerry Stevens   }
158751dcc849SKerry Stevens   free(BarrPoint);
158851dcc849SKerry Stevens   free(PetscThreadPoint);
1589cfcfc605SKerry Stevens 
159051dcc849SKerry Stevens   PetscFunctionReturn(0);
159151dcc849SKerry Stevens }
159251dcc849SKerry Stevens 
1593f09cb4aaSKerry Stevens #undef __FUNCT__
159451d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
159551d315f7SKerry Stevens void MainWait_True() {
159651dcc849SKerry Stevens   int ierr;
15976c9b208dSKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
159851d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
159951d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
160051dcc849SKerry Stevens   }
160151d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
160251dcc849SKerry Stevens }
160351dcc849SKerry Stevens 
1604f09cb4aaSKerry Stevens #undef __FUNCT__
160551d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
160651d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
160751dcc849SKerry Stevens   int ierr;
16080ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
160936d20dc5SKerry Stevens 
16100ca81413SKerry Stevens   MainWait();
161151d315f7SKerry Stevens   job_true.pfunc = pFunc;
161251d315f7SKerry Stevens   job_true.pdata = data;
161351d315f7SKerry Stevens   job_true.pbarr = &BarrPoint[n];
161451d315f7SKerry Stevens   job_true.iNumJobThreads = n;
161551d315f7SKerry Stevens   job_true.startJob = PETSC_TRUE;
161651d315f7SKerry Stevens   ierr = pthread_cond_broadcast(&job_true.cond);
16170ca81413SKerry Stevens   if(pFunc!=FuncFinish) {
1618ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done */
16190ca81413SKerry Stevens   }
162036d20dc5SKerry Stevens 
16210ca81413SKerry Stevens   if(ithreaderr) {
16220ca81413SKerry Stevens     ijoberr = ithreaderr;
16230ca81413SKerry Stevens   }
16240ca81413SKerry Stevens   return ijoberr;
162551dcc849SKerry Stevens }
1626*fdfc40dbSShri Abhyankar #endif
1627*fdfc40dbSShri Abhyankar 
1628683509dcSKerry Stevens /**** NO THREAD POOL FUNCTION ****/
1629683509dcSKerry Stevens #undef __FUNCT__
1630683509dcSKerry Stevens #define __FUNCT__ "MainJob_Spawn"
1631683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void** data,PetscInt n) {
1632683509dcSKerry Stevens   PetscErrorCode ijoberr = 0;
1633683509dcSKerry Stevens 
1634683509dcSKerry Stevens   pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
1635cfcfc605SKerry Stevens   PetscThreadPoint = apThread; /* point to same place */
1636683509dcSKerry Stevens   PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1637683509dcSKerry Stevens   PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
1638683509dcSKerry Stevens   free(apThread);
1639683509dcSKerry Stevens 
1640683509dcSKerry Stevens   return ijoberr;
1641683509dcSKerry Stevens }
164251d315f7SKerry Stevens /****  ****/
1643ba61063dSBarry Smith #endif
164451dcc849SKerry Stevens 
164551dcc849SKerry Stevens void* FuncFinish(void* arg) {
164651dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16470ca81413SKerry Stevens   return(0);
164851dcc849SKerry Stevens }
1649