xref: /petsc/src/sys/objects/init.c (revision cd723cd1e71baf63298996d3aa19d25195c37816)
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 */
84d4e02f7SKerry Stevens //#if defined(PETSC_HAVE_SCHED_H) && defined(PETSC_USE_PTHREAD)
98f54b378SBarry Smith #ifndef _GNU_SOURCE
108f54b378SBarry Smith #define _GNU_SOURCE
118f54b378SBarry Smith #endif
12*cd723cd1SBarry Smith #if defined(PETSC_HAVE_SCHED_H)
138f54b378SBarry Smith #include <sched.h>
14*cd723cd1SBarry Smith #endif
154d4e02f7SKerry Stevens #include <petscsys.h>        /*I  "petscsys.h"   I*/
16ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD)
1751dcc849SKerry Stevens #include <pthread.h>
18ba61063dSBarry Smith #endif
19ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
2051d315f7SKerry Stevens #include <sys/sysinfo.h>
21ba61063dSBarry Smith #endif
2251d315f7SKerry Stevens #include <unistd.h>
23e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
24e5c89e4eSSatish Balay #include <stdlib.h>
25e5c89e4eSSatish Balay #endif
26e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
27e5c89e4eSSatish Balay #include <malloc.h>
28e5c89e4eSSatish Balay #endif
29555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
30555d055bSBarry Smith #include <valgrind/valgrind.h>
31555d055bSBarry Smith #endif
32555d055bSBarry Smith 
33e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
34e5c89e4eSSatish Balay /*
35e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
36e5c89e4eSSatish Balay    already started before PETSc was initialized.
37e5c89e4eSSatish Balay */
387087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
397087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
407087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
4151dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
4251dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
437087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
447087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
45ba61063dSBarry Smith 
46ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
4751dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
4851dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
49af359df3SBarry Smith #define PETSC_HAVE_PTHREAD_BARRIER
50ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
51ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
52ba61063dSBarry Smith #endif
5351d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
54f09cb4aaSKerry Stevens int*         pVal;
5551dcc849SKerry Stevens 
56ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
5751d315f7SKerry Stevens int* ThreadCoreAffinity;
5851d315f7SKerry Stevens 
59ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
6051d315f7SKerry Stevens 
6151d315f7SKerry Stevens typedef struct {
6251d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
6351d315f7SKerry Stevens   pthread_cond_t**  cond1array;
6451d315f7SKerry Stevens   pthread_cond_t** cond2array;
6551d315f7SKerry Stevens   void* (*pfunc)(void*);
6651d315f7SKerry Stevens   void** pdata;
6751d315f7SKerry Stevens   PetscBool startJob;
6851d315f7SKerry Stevens   estat eJobStat;
6951d315f7SKerry Stevens   PetscBool** arrThreadStarted;
7051d315f7SKerry Stevens   PetscBool** arrThreadReady;
7151d315f7SKerry Stevens } sjob_tree;
7251d315f7SKerry Stevens sjob_tree job_tree;
7351d315f7SKerry Stevens typedef struct {
7451d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
7551d315f7SKerry Stevens   pthread_cond_t**  cond1array;
7651d315f7SKerry Stevens   pthread_cond_t** cond2array;
7751d315f7SKerry Stevens   void* (*pfunc)(void*);
7851d315f7SKerry Stevens   void** pdata;
7951d315f7SKerry Stevens   PetscBool** arrThreadReady;
8051d315f7SKerry Stevens } sjob_main;
8151d315f7SKerry Stevens sjob_main job_main;
8251d315f7SKerry Stevens typedef struct {
8351d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8451d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8551d315f7SKerry Stevens   pthread_cond_t** cond2array;
8651d315f7SKerry Stevens   void* (*pfunc)(void*);
8751d315f7SKerry Stevens   void** pdata;
8851d315f7SKerry Stevens   PetscBool startJob;
8951d315f7SKerry Stevens   estat eJobStat;
9051d315f7SKerry Stevens   PetscBool** arrThreadStarted;
9151d315f7SKerry Stevens   PetscBool** arrThreadReady;
9251d315f7SKerry Stevens } sjob_chain;
9351d315f7SKerry Stevens sjob_chain job_chain;
94ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
9551dcc849SKerry Stevens typedef struct {
9651dcc849SKerry Stevens   pthread_mutex_t mutex;
9751dcc849SKerry Stevens   pthread_cond_t cond;
9851dcc849SKerry Stevens   void* (*pfunc)(void*);
9951dcc849SKerry Stevens   void** pdata;
10051dcc849SKerry Stevens   pthread_barrier_t* pbarr;
10151dcc849SKerry Stevens   int iNumJobThreads;
10251dcc849SKerry Stevens   int iNumReadyThreads;
10351dcc849SKerry Stevens   PetscBool startJob;
10451d315f7SKerry Stevens } sjob_true;
10551d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
106ba61063dSBarry Smith #endif
10751dcc849SKerry Stevens 
108ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
109ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
110ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
111ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
112ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
113ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
11451dcc849SKerry Stevens 
11551d315f7SKerry Stevens /* Function Pointers */
11651d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
11751d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
11851d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
11951d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
12051d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
12136d20dc5SKerry Stevens /**** Tree Thread Pool Functions ****/
12251d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
12351d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
12451d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
12551d315f7SKerry Stevens void           MainWait_Tree(void);
12651d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
12736d20dc5SKerry Stevens /**** Main Thread Pool Functions ****/
12851d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
12951d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
13051d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
13151d315f7SKerry Stevens void           MainWait_Main(void);
13251d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
13336d20dc5SKerry Stevens /**** Chain Thread Pool Functions ****/
13451d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
13551d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
13651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
13751d315f7SKerry Stevens void           MainWait_Chain(void);
13851d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
13936d20dc5SKerry Stevens /**** True Thread Pool Functions ****/
14051d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
14151d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
14251d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
14351d315f7SKerry Stevens void           MainWait_True(void);
14451d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
14536d20dc5SKerry Stevens /**** NO Thread Pool Function  ****/
146683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void**,PetscInt);
14736d20dc5SKerry Stevens /****  ****/
14851dcc849SKerry Stevens void* FuncFinish(void*);
1490ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1500ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
151ba61063dSBarry Smith #endif
152e5c89e4eSSatish Balay 
153e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
154e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
155e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
156e5c89e4eSSatish Balay #endif
1572c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1587087cfbeSBarry Smith MPI_Datatype   MPI_C_DOUBLE_COMPLEX;
1597087cfbeSBarry Smith MPI_Datatype   MPI_C_COMPLEX;
1602c876bd9SBarry Smith #endif
1617087cfbeSBarry Smith PetscScalar    PETSC_i;
162e5c89e4eSSatish Balay #else
1637087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
164e5c89e4eSSatish Balay #endif
165ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
166c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
167c90a1750SBarry Smith #endif
1687087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1697087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
17075567043SBarry Smith 
171e5c89e4eSSatish Balay /*
172e5c89e4eSSatish Balay      These are needed by petscbt.h
173e5c89e4eSSatish Balay */
174c6db04a5SJed Brown #include <petscbt.h>
1757087cfbeSBarry Smith char      _BT_mask = ' ';
1767087cfbeSBarry Smith char      _BT_c = ' ';
1777087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
178e5c89e4eSSatish Balay 
179e5c89e4eSSatish Balay /*
180e5c89e4eSSatish Balay        Function that is called to display all error messages
181e5c89e4eSSatish Balay */
1827087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1837087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
184238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1857087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
186238ccf28SShri Abhyankar #else
1877087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
188238ccf28SShri Abhyankar #endif
189bab1f7e6SVictor Minden /*
1908154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1918154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
192bab1f7e6SVictor Minden 
193e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
194e5c89e4eSSatish Balay /*
195e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
196e5c89e4eSSatish Balay */
197e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
198e5c89e4eSSatish Balay 
199e5c89e4eSSatish Balay #undef __FUNCT__
200f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
2017087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
202e5c89e4eSSatish Balay {
203e5c89e4eSSatish Balay   PetscErrorCode ierr;
204e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
205e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
206e5c89e4eSSatish Balay   char           version[256];
207e5c89e4eSSatish Balay 
208e5c89e4eSSatish Balay   PetscFunctionBegin;
209e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
210e5c89e4eSSatish Balay   if (!rank) {
211e5c89e4eSSatish Balay     char        arch[10];
212f56c2debSBarry Smith     int         err;
21388c29154SBarry Smith     PetscViewer viewer;
214f56c2debSBarry Smith 
215e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
216e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
217a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
218e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
219e5c89e4eSSatish Balay     if (filename) {
220e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
221e5c89e4eSSatish Balay     } else {
222e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
223e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
224e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
225e5c89e4eSSatish Balay     }
226e5c89e4eSSatish Balay 
227e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
228e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
229e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
230e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
231e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
23288c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
23388c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2346bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
235e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
236f56c2debSBarry Smith     err = fflush(*fd);
237e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
238e5c89e4eSSatish Balay   }
239e5c89e4eSSatish Balay   PetscFunctionReturn(0);
240e5c89e4eSSatish Balay }
241e5c89e4eSSatish Balay 
242e5c89e4eSSatish Balay #undef __FUNCT__
243f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2447087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
245e5c89e4eSSatish Balay {
246e5c89e4eSSatish Balay   PetscErrorCode ierr;
247e5c89e4eSSatish Balay   PetscMPIInt    rank;
248e5c89e4eSSatish Balay   char           date[64];
249f56c2debSBarry Smith   int            err;
250e5c89e4eSSatish Balay 
251e5c89e4eSSatish Balay   PetscFunctionBegin;
252e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
253e5c89e4eSSatish Balay   if (!rank) {
254e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
255e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
256e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
257e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
258f56c2debSBarry Smith     err = fflush(*fd);
259e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
260f56c2debSBarry Smith     err = fclose(*fd);
261e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
262e5c89e4eSSatish Balay   }
263e5c89e4eSSatish Balay   PetscFunctionReturn(0);
264e5c89e4eSSatish Balay }
265e5c89e4eSSatish Balay 
266e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
267e5c89e4eSSatish Balay 
268e5c89e4eSSatish Balay /*
269e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
270e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
271e5c89e4eSSatish Balay 
272e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2733c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2743c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
275e5c89e4eSSatish Balay */
276e5c89e4eSSatish Balay 
277e5c89e4eSSatish Balay #undef __FUNCT__
278e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
279e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
280e5c89e4eSSatish Balay {
281e5c89e4eSSatish Balay   PetscFunctionBegin;
2823c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
283e5c89e4eSSatish Balay   abort();
284e5c89e4eSSatish Balay }
285e5c89e4eSSatish Balay 
286e5c89e4eSSatish Balay #undef __FUNCT__
287e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
288e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
289e5c89e4eSSatish Balay {
290e5c89e4eSSatish Balay   PetscErrorCode ierr;
291e5c89e4eSSatish Balay 
292e5c89e4eSSatish Balay   PetscFunctionBegin;
2933c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
294e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
295e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
2963c311c98SBarry Smith     MPI_Abort(*comm,*flag);
297e5c89e4eSSatish Balay   }
298e5c89e4eSSatish Balay }
299e5c89e4eSSatish Balay 
300e5c89e4eSSatish Balay #undef __FUNCT__
301e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
302e5c89e4eSSatish Balay /*@C
303e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
304e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
305e5c89e4eSSatish Balay 
306e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
307e5c89e4eSSatish Balay 
308e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
309e5c89e4eSSatish Balay 
310e5c89e4eSSatish Balay    Level: advanced
311e5c89e4eSSatish Balay 
312e5c89e4eSSatish Balay    Note:
313e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
314e5c89e4eSSatish Balay 
31588c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
316e5c89e4eSSatish Balay @*/
3177087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
318e5c89e4eSSatish Balay {
319e5c89e4eSSatish Balay   PetscFunctionBegin;
320e5c89e4eSSatish Balay   PetscFinalize();
321e5c89e4eSSatish Balay   exit(0);
322e5c89e4eSSatish Balay   return 0;
323e5c89e4eSSatish Balay }
324e5c89e4eSSatish Balay 
325ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
32609573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
327ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
328e5c89e4eSSatish Balay static char       emacsmachinename[256];
329e5c89e4eSSatish Balay 
330e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
331e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
332e5c89e4eSSatish Balay 
333e5c89e4eSSatish Balay #undef __FUNCT__
334e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
335e5c89e4eSSatish Balay /*@C
336e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
337e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
338e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
339e5c89e4eSSatish Balay 
340e5c89e4eSSatish Balay    Input Parameter:
341e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
342da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
343e5c89e4eSSatish Balay 
344e5c89e4eSSatish Balay    Level: developer
345e5c89e4eSSatish Balay 
346e5c89e4eSSatish Balay    Concepts: package help message
347e5c89e4eSSatish Balay 
348e5c89e4eSSatish Balay @*/
3497087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
350e5c89e4eSSatish Balay {
351e5c89e4eSSatish Balay   PetscFunctionBegin;
352e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
353e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
354e5c89e4eSSatish Balay   PetscFunctionReturn(0);
355e5c89e4eSSatish Balay }
356e5c89e4eSSatish Balay 
357e5c89e4eSSatish Balay #undef __FUNCT__
358e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3597087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
360e5c89e4eSSatish Balay {
361e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
362e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
363ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
364e5c89e4eSSatish Balay   PetscErrorCode ierr;
365a6d0e24fSJed Brown   PetscReal      si;
366e5c89e4eSSatish Balay   int            i;
367e5c89e4eSSatish Balay   PetscMPIInt    rank;
368e5c89e4eSSatish Balay   char           version[256];
369e5c89e4eSSatish Balay 
370e5c89e4eSSatish Balay   PetscFunctionBegin;
371e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
372e5c89e4eSSatish Balay 
373e5c89e4eSSatish Balay   /*
374e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
375e5c89e4eSSatish Balay   */
3768bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
37781b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
378acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
379e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
380555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
381555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
382555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
383555d055bSBarry Smith #endif
384e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
385555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
386555d055bSBarry Smith     }
387555d055bSBarry Smith #endif
388e5c89e4eSSatish Balay   }
389e5c89e4eSSatish Balay #else
390acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
391acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
392e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
393e5c89e4eSSatish Balay #endif
394e5c89e4eSSatish Balay   if (flg3) {
395e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
396e5c89e4eSSatish Balay   }
39790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
398acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
399e5c89e4eSSatish Balay   if (flg1) {
400e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
401e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
402e5c89e4eSSatish Balay   }
403e5c89e4eSSatish Balay 
40490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
405acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4067783f70dSSatish Balay   if (!flg1) {
40790d69ab7SBarry Smith     flg1 = PETSC_FALSE;
408acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4097783f70dSSatish Balay   }
410e5c89e4eSSatish Balay   if (flg1) {
411e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
412e5c89e4eSSatish Balay   }
413e5c89e4eSSatish Balay 
414e5c89e4eSSatish Balay   /*
415e5c89e4eSSatish Balay       Set the display variable for graphics
416e5c89e4eSSatish Balay   */
417e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
418e5c89e4eSSatish Balay 
41951dcc849SKerry Stevens 
42051dcc849SKerry Stevens   /*
421e5c89e4eSSatish Balay       Print the PETSc version information
422e5c89e4eSSatish Balay   */
423e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
424e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
425e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
426e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
427e5c89e4eSSatish Balay 
428e5c89e4eSSatish Balay     /*
429e5c89e4eSSatish Balay        Print "higher-level" package version message
430e5c89e4eSSatish Balay     */
431e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
432e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
433e5c89e4eSSatish Balay     }
434e5c89e4eSSatish Balay 
435a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
436e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
437e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
438e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
439e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
440e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
44184e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
442e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
443e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
444e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
445e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
446e5c89e4eSSatish Balay   }
447e5c89e4eSSatish Balay 
448e5c89e4eSSatish Balay   /*
449e5c89e4eSSatish Balay        Print "higher-level" package help message
450e5c89e4eSSatish Balay   */
451e5c89e4eSSatish Balay   if (flg3){
452e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
453e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
454e5c89e4eSSatish Balay     }
455e5c89e4eSSatish Balay   }
456e5c89e4eSSatish Balay 
457e5c89e4eSSatish Balay   /*
458e5c89e4eSSatish Balay       Setup the error handling
459e5c89e4eSSatish Balay   */
46090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
461acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
462cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
46390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
464acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
465cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
46690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
467acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
468e5c89e4eSSatish Balay   if (flg1) {
469e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
470e5c89e4eSSatish Balay   }
47190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
472acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
473cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
47496cc47afSJed Brown   flg1 = PETSC_FALSE;
475acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
47696cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
477e5c89e4eSSatish Balay 
478e5c89e4eSSatish Balay   /*
479e5c89e4eSSatish Balay       Setup debugger information
480e5c89e4eSSatish Balay   */
481e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
482e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
483e5c89e4eSSatish Balay   if (flg1) {
484e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
485e5c89e4eSSatish Balay 
486e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
487e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
488e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
489e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
490e5c89e4eSSatish Balay   }
4915e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
4925e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
493e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
494e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
495e5c89e4eSSatish Balay   if (flg1 || flg2) {
496e5c89e4eSSatish Balay     PetscMPIInt    size;
497e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
498e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
499e5c89e4eSSatish Balay     /*
500e5c89e4eSSatish Balay        we have to make sure that all processors have opened
501e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
502e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
503e5c89e4eSSatish Balay        and kill the program.
504e5c89e4eSSatish Balay     */
505e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
506e5c89e4eSSatish Balay     if (size > 2) {
507533163c2SBarry Smith       PetscMPIInt dummy = 0;
508e5c89e4eSSatish Balay       MPI_Status  status;
509e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
510e5c89e4eSSatish Balay         if (rank != i) {
511e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
512e5c89e4eSSatish Balay         }
513e5c89e4eSSatish Balay       }
514e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
515e5c89e4eSSatish Balay         if (rank != i) {
516e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
517e5c89e4eSSatish Balay         }
518e5c89e4eSSatish Balay       }
519e5c89e4eSSatish Balay     }
520e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
521e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
522e5c89e4eSSatish Balay     lsize = size;
523e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
524e5c89e4eSSatish Balay     if (flag) {
525e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
526e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
527e5c89e4eSSatish Balay       }
528e5c89e4eSSatish Balay     }
529e5c89e4eSSatish Balay     if (!flag) {
530e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
531e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
532e5c89e4eSSatish Balay       if (flg1) {
533e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
534e5c89e4eSSatish Balay       } else {
535e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
536e5c89e4eSSatish Balay       }
537e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
538e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
539e5c89e4eSSatish Balay     }
540e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
541e5c89e4eSSatish Balay   }
542e5c89e4eSSatish Balay 
543e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
544cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
545e5c89e4eSSatish Balay 
54693ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
54722b84c2fSbcordonn   /*
54822b84c2fSbcordonn     Activates new sockets for zope if needed
54922b84c2fSbcordonn   */
55084ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
551d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5526dc8fec2Sbcordonn   if (flgz){
55322b84c2fSbcordonn     int  sockfd;
554f1384234SBarry Smith     char hostname[256];
55522b84c2fSbcordonn     char username[256];
5566dc8fec2Sbcordonn     int  remoteport = 9999;
5579c4c166aSBarry Smith 
55884ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
55984ab5442Sbcordonn     if (!hostname[0]){
5609c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5619c4c166aSBarry Smith     }
56222b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5639c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
56422b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
56522b84c2fSbcordonn     if (flgzout){
56622b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
567606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5686dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5699c4c166aSBarry Smith     } else {
570d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
571d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5729c4c166aSBarry Smith     }
5739c4c166aSBarry Smith   }
57493ba235fSBarry Smith #endif
575ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
576ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
577ffc871a5SBarry Smith   if (flgz){
578ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
579ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
580ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
581ffc871a5SBarry Smith   }
582ffc871a5SBarry Smith #endif
5836dc8fec2Sbcordonn 
584e5c89e4eSSatish Balay   /*
585e5c89e4eSSatish Balay         Setup profiling and logging
586e5c89e4eSSatish Balay   */
5876cf91177SBarry Smith #if defined (PETSC_USE_INFO)
5888bb29257SSatish Balay   {
589e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
5906cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
5918bb29257SSatish Balay     if (flg1 && logname[0]) {
592fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
5938bb29257SSatish Balay     } else if (flg1) {
594fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
595e5c89e4eSSatish Balay     }
596e5c89e4eSSatish Balay   }
597865f6aa8SSatish Balay #endif
598865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
599865f6aa8SSatish Balay   mname[0] = 0;
600f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
601865f6aa8SSatish Balay   if (flg1) {
602865f6aa8SSatish Balay     if (mname[0]) {
603f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
604865f6aa8SSatish Balay     } else {
605f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
606865f6aa8SSatish Balay     }
607865f6aa8SSatish Balay   }
608e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
60990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
610fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
611e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
612e5c89e4eSSatish Balay #endif
61390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
61490d69ab7SBarry Smith   flg2 = PETSC_FALSE;
61590d69ab7SBarry Smith   flg3 = PETSC_FALSE;
616acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
617acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
618d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6199f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
620e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6219f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
622e5c89e4eSSatish Balay 
623e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
624e5c89e4eSSatish Balay   if (flg1) {
625e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
626e5c89e4eSSatish Balay     FILE *file;
627e5c89e4eSSatish Balay     if (mname[0]) {
628e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
629e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
630e5c89e4eSSatish Balay       file = fopen(fname,"w");
631f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
632e5c89e4eSSatish Balay     } else {
633da9f1d6bSBarry Smith       file = PETSC_STDOUT;
634e5c89e4eSSatish Balay     }
635e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
636e5c89e4eSSatish Balay   }
637e5c89e4eSSatish Balay #endif
638e5c89e4eSSatish Balay 
639e5c89e4eSSatish Balay   /*
640e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
641e5c89e4eSSatish Balay   */
64263d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
643e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
644e5c89e4eSSatish Balay #endif
645e5c89e4eSSatish Balay 
646acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
647e5c89e4eSSatish Balay 
648ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
649af359df3SBarry Smith   /*
650af359df3SBarry Smith       Determine whether user specified maximum number of threads
651af359df3SBarry Smith    */
652af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
653af359df3SBarry Smith 
654b154f58aSKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-main",&flg1);CHKERRQ(ierr);
655b154f58aSKerry Stevens   if(flg1) {
656b154f58aSKerry Stevens     cpu_set_t mset;
657b154f58aSKerry Stevens     int icorr,ncorr = get_nprocs();
658b154f58aSKerry Stevens     ierr = PetscOptionsGetInt(PETSC_NULL,"-main",&icorr,PETSC_NULL);CHKERRQ(ierr);
659b154f58aSKerry Stevens     CPU_ZERO(&mset);
660b154f58aSKerry Stevens     CPU_SET(icorr%ncorr,&mset);
661b154f58aSKerry Stevens     sched_setaffinity(0,sizeof(cpu_set_t),&mset);
662b154f58aSKerry Stevens   }
663b154f58aSKerry Stevens 
664af359df3SBarry Smith   PetscInt N_CORES = get_nprocs();
665af359df3SBarry Smith   ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
666af359df3SBarry Smith   char tstr[9];
667af359df3SBarry Smith   char tbuf[2];
668af359df3SBarry Smith   strcpy(tstr,"-thread");
669af359df3SBarry Smith   for(i=0;i<PetscMaxThreads;i++) {
670af359df3SBarry Smith     ThreadCoreAffinity[i] = i;
671af359df3SBarry Smith     sprintf(tbuf,"%d",i);
672af359df3SBarry Smith     strcat(tstr,tbuf);
673af359df3SBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
674af359df3SBarry Smith     if(flg1) {
675af359df3SBarry Smith       ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
676af359df3SBarry Smith       ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
677af359df3SBarry Smith     }
678af359df3SBarry Smith     tstr[7] = '\0';
679af359df3SBarry Smith   }
680cfcfc605SKerry Stevens 
681cfcfc605SKerry Stevens   /*
682cfcfc605SKerry Stevens       Determine whether to use thread pool
683cfcfc605SKerry Stevens    */
684cfcfc605SKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
685cfcfc605SKerry Stevens   if (flg1) {
686cfcfc605SKerry Stevens     PetscUseThreadPool = PETSC_TRUE;
687af359df3SBarry Smith     /* get the thread pool type */
688af359df3SBarry Smith     PetscInt ipool = 0;
689af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
690af359df3SBarry Smith 
691af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
692af359df3SBarry Smith     switch(ipool) {
693af359df3SBarry Smith     case 1:
694af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
695af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
696af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
697af359df3SBarry Smith       MainWait              = &MainWait_Tree;
698af359df3SBarry Smith       MainJob               = &MainJob_Tree;
699af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
700af359df3SBarry Smith       break;
701af359df3SBarry Smith     case 2:
702af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
703af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
704af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
705af359df3SBarry Smith       MainWait              = &MainWait_Main;
706af359df3SBarry Smith       MainJob               = &MainJob_Main;
707af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
708af359df3SBarry Smith       break;
709af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
710af359df3SBarry Smith     case 3:
711af359df3SBarry Smith #else
712af359df3SBarry Smith     default:
713af359df3SBarry Smith #endif
714af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
715af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
716af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
717af359df3SBarry Smith       MainWait              = &MainWait_Chain;
718af359df3SBarry Smith       MainJob               = &MainJob_Chain;
719af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
720af359df3SBarry Smith       break;
721af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
722af359df3SBarry Smith     default:
723af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
724af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
725af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
726af359df3SBarry Smith       MainWait              = &MainWait_True;
727af359df3SBarry Smith       MainJob               = &MainJob_True;
728af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
729af359df3SBarry Smith       break;
730af359df3SBarry Smith #endif
731af359df3SBarry Smith     }
732af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
733683509dcSKerry Stevens   } else {
734683509dcSKerry Stevens     //need to define these in the case on 'no threads' or 'thread create/destroy'
735683509dcSKerry Stevens     //could take any of the above versions
736683509dcSKerry Stevens     MainJob               = &MainJob_Spawn;
737af359df3SBarry Smith   }
738af359df3SBarry Smith #endif
739e5c89e4eSSatish Balay   /*
740e5c89e4eSSatish Balay        Print basic help message
741e5c89e4eSSatish Balay   */
742e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
743e5c89e4eSSatish Balay   if (flg1) {
744e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
745301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
746301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
747301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
748e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
749e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
750e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
751e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
752e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
753e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
754e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
755e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
756e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
757e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
758e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
759e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
760e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
761e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
762e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
763e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
764e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
765e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
766e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7674161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7684161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
769e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
770e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
771e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
772e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
773e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
774e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
775a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
776e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
77740ab9619SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -server <port>: Run PETSc webserver (default port is 8080) see PetscWebServe()\n");CHKERRQ(ierr);
778e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
779e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
7807ef452c0SMatthew G Knepley     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary _summary_python]: logging objects and events\n");CHKERRQ(ierr);
781e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
782e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
783e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
784e5c89e4eSSatish Balay #endif
7856cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
786e5c89e4eSSatish Balay #endif
787e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
788e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
789e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
790e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
791e5c89e4eSSatish Balay   }
792e5c89e4eSSatish Balay 
793a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
794e5c89e4eSSatish Balay   if (flg1) {
795e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
796e5c89e4eSSatish Balay   }
797e5c89e4eSSatish Balay 
7986cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
799e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
800e5c89e4eSSatish Balay   if (f) {
8016cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
802e5c89e4eSSatish Balay   }
803827f890bSBarry Smith 
8048154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
805c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
80673113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
80773113deaSBarry Smith   else flg1 = PETSC_FALSE;
8088154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
8098154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
810bab1f7e6SVictor Minden #endif
811192daf7cSBarry Smith 
812e5c89e4eSSatish Balay   PetscFunctionReturn(0);
813e5c89e4eSSatish Balay }
814df413903SBarry Smith 
815ff34cdc8SBarry Smith #if defined(PETSC_HAVE_PTHREADCLASSES)
816ba61063dSBarry Smith 
81751d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
81851d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
81951d315f7SKerry Stevens   PetscErrorCode iterr;
82051d315f7SKerry Stevens   int icorr,ierr;
82151d315f7SKerry Stevens   int* pId = (int*)arg;
82251d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
82351d315f7SKerry Stevens   PetscBool PeeOn;
82451d315f7SKerry Stevens   cpu_set_t mset;
8259e800a48SKerry Stevens   //printf("Thread %d In Tree Thread Function\n",ThreadId);
82651d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
82751d315f7SKerry Stevens   CPU_ZERO(&mset);
82851d315f7SKerry Stevens   CPU_SET(icorr,&mset);
82951d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
83051d315f7SKerry Stevens 
83151d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
83251d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
83351d315f7SKerry Stevens   }
83451d315f7SKerry Stevens   else {
83551d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
83651d315f7SKerry Stevens   }
83751d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
838ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
83951d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
84051d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
84151d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
84251d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
84351d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
844ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
845ba61063dSBarry Smith            upon return, has the lock */
84651d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
84751d315f7SKerry Stevens         }
84851d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
84951d315f7SKerry Stevens       }
85051d315f7SKerry Stevens     }
851ba61063dSBarry Smith     /* your subordinates are now ready */
85251d315f7SKerry Stevens   }
85351d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
854ba61063dSBarry Smith   /* update your ready status */
85551d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
85651d315f7SKerry Stevens   if(ThreadId==0) {
85751d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
858ba61063dSBarry Smith     /* ignal main */
85951d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
86051d315f7SKerry Stevens   }
86151d315f7SKerry Stevens   else {
862ba61063dSBarry Smith     /* tell your boss that you're ready to work */
86351d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
86451d315f7SKerry Stevens   }
865ba61063dSBarry Smith   /* the while loop needs to have an exit
866ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
867ba61063dSBarry Smith    and calling FuncFinish */
86851d315f7SKerry Stevens   while(PetscThreadGo) {
869ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
870ba61063dSBarry Smith       waiting when you don't have to causes problems
871ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
87251d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
873ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
874ba61063dSBarry Smith        upon return, has the lock */
87551d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
87651d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
87751d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
87851d315f7SKerry Stevens     }
87951d315f7SKerry Stevens     if(ThreadId==0) {
88051d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
88151d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
88251d315f7SKerry Stevens     }
88351d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
88451d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
885ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
88651d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
88751d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
88851d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
88951d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
89051d315f7SKerry Stevens         }
89151d315f7SKerry Stevens       }
89251d315f7SKerry Stevens     }
893ba61063dSBarry Smith     /* do your job */
89451d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
89551d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
89651d315f7SKerry Stevens     }
89751d315f7SKerry Stevens     else {
89851d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
89951d315f7SKerry Stevens     }
90051d315f7SKerry Stevens     if(iterr!=0) {
90151d315f7SKerry Stevens       ithreaderr = 1;
90251d315f7SKerry Stevens     }
90351d315f7SKerry Stevens     if(PetscThreadGo) {
904ba61063dSBarry Smith       /* reset job, get ready for more */
90551d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
906ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
907ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
90851d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
90951d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
91051d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
91151d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
91251d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
913ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
914ba61063dSBarry Smith                upon return, has the lock */
91551d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
91651d315f7SKerry Stevens             }
91751d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
91851d315f7SKerry Stevens           }
91951d315f7SKerry Stevens 	}
920ba61063dSBarry Smith         /* your subordinates are now ready */
92151d315f7SKerry Stevens       }
92251d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
92351d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
92451d315f7SKerry Stevens       if(ThreadId==0) {
925ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
926ba61063dSBarry Smith         /* root thread signals 'main' */
92751d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
92851d315f7SKerry Stevens       }
92951d315f7SKerry Stevens       else {
930ba61063dSBarry Smith         /* signal your boss before you go to sleep */
93151d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
93251d315f7SKerry Stevens       }
93351d315f7SKerry Stevens     }
93451d315f7SKerry Stevens   }
93551d315f7SKerry Stevens   return NULL;
93651d315f7SKerry Stevens }
93751d315f7SKerry Stevens 
93851d315f7SKerry Stevens #undef __FUNCT__
93951d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
94051d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
94151d315f7SKerry Stevens   PetscInt i,ierr;
94251d315f7SKerry Stevens   int status;
94351d315f7SKerry Stevens 
94451d315f7SKerry Stevens   if(PetscUseThreadPool) {
94551d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
94651d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
94751d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
94851d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
94951d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
95051d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
95151d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
95251d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
95351d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95451d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95551d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
95651d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
957ba61063dSBarry Smith     /* initialize job structure */
95851d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
95951d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
96051d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
96151d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
96251d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
96351d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
96451d315f7SKerry Stevens     }
96551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
96651d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
96751d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
96851d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
96951d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
97051d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
97151d315f7SKerry Stevens     }
97251d315f7SKerry Stevens     job_tree.pfunc = NULL;
97351d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
97451d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
97551d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
97651d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
977ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
97851d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
979ba61063dSBarry Smith     /* create threads */
98051d315f7SKerry Stevens     for(i=0; i<N; i++) {
98151d315f7SKerry Stevens       pVal[i] = i;
98251d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
983ba61063dSBarry Smith       /* should check status */
98451d315f7SKerry Stevens     }
98551d315f7SKerry Stevens   }
98651d315f7SKerry Stevens   return NULL;
98751d315f7SKerry Stevens }
98851d315f7SKerry Stevens 
98951d315f7SKerry Stevens #undef __FUNCT__
99051d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
99151d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
99251d315f7SKerry Stevens   int i,ierr;
99351d315f7SKerry Stevens   void* jstatus;
99451d315f7SKerry Stevens 
99551d315f7SKerry Stevens   PetscFunctionBegin;
99651d315f7SKerry Stevens 
997ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
998ba61063dSBarry Smith   /* join the threads */
99951d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
100051d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1001ba61063dSBarry Smith     /* do error checking*/
100251d315f7SKerry Stevens   }
100351d315f7SKerry Stevens   free(PetscThreadPoint);
100451d315f7SKerry Stevens   free(arrmutex);
100551d315f7SKerry Stevens   free(arrcond1);
100651d315f7SKerry Stevens   free(arrcond2);
100751d315f7SKerry Stevens   free(arrstart);
100851d315f7SKerry Stevens   free(arrready);
100951d315f7SKerry Stevens   free(job_tree.pdata);
101051d315f7SKerry Stevens   free(pVal);
1011cfcfc605SKerry Stevens 
101251d315f7SKerry Stevens   PetscFunctionReturn(0);
101351d315f7SKerry Stevens }
101451d315f7SKerry Stevens 
101551d315f7SKerry Stevens #undef __FUNCT__
101651d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
101751d315f7SKerry Stevens void MainWait_Tree() {
101851d315f7SKerry Stevens   int ierr;
101951d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
102051d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
102151d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
102251d315f7SKerry Stevens   }
102351d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
102451d315f7SKerry Stevens }
102551d315f7SKerry Stevens 
102651d315f7SKerry Stevens #undef __FUNCT__
102751d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
102851d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
102951d315f7SKerry Stevens   int i,ierr;
103051d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
103136d20dc5SKerry Stevens 
103251d315f7SKerry Stevens   MainWait();
103351d315f7SKerry Stevens   job_tree.pfunc = pFunc;
103451d315f7SKerry Stevens   job_tree.pdata = data;
103551d315f7SKerry Stevens   job_tree.startJob = PETSC_TRUE;
103651d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
103751d315f7SKerry Stevens     *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
103851d315f7SKerry Stevens   }
103951d315f7SKerry Stevens   job_tree.eJobStat = JobInitiated;
104051d315f7SKerry Stevens   ierr = pthread_cond_signal(job_tree.cond2array[0]);
104151d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1042ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
104351d315f7SKerry Stevens   }
104436d20dc5SKerry Stevens 
104551d315f7SKerry Stevens   if(ithreaderr) {
104651d315f7SKerry Stevens     ijoberr = ithreaderr;
104751d315f7SKerry Stevens   }
104851d315f7SKerry Stevens   return ijoberr;
104951d315f7SKerry Stevens }
105051d315f7SKerry Stevens /****  ****/
105151d315f7SKerry Stevens 
105251d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
105351d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
105451d315f7SKerry Stevens   PetscErrorCode iterr;
105551d315f7SKerry Stevens   int icorr,ierr;
105651d315f7SKerry Stevens   int* pId = (int*)arg;
105751d315f7SKerry Stevens   int ThreadId = *pId;
105851d315f7SKerry Stevens   cpu_set_t mset;
10599e800a48SKerry Stevens   //printf("Thread %d In Main Thread Function\n",ThreadId);
106051d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
106151d315f7SKerry Stevens   CPU_ZERO(&mset);
106251d315f7SKerry Stevens   CPU_SET(icorr,&mset);
106351d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
106451d315f7SKerry Stevens 
106551d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1066ba61063dSBarry Smith   /* update your ready status */
106751d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1068ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
106951d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
107051d315f7SKerry Stevens 
1071ba61063dSBarry Smith   /* the while loop needs to have an exit
1072ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1073ba61063dSBarry Smith      and calling FuncFinish */
107451d315f7SKerry Stevens   while(PetscThreadGo) {
1075ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1076ba61063dSBarry Smith        waiting when you don't have to causes problems
1077ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
107851d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1079ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1080ba61063dSBarry Smith        upon return, has the lock */
108151d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1082ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
108351d315f7SKerry Stevens     }
108451d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
108551d315f7SKerry Stevens     if(job_main.pdata==NULL) {
108651d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
108751d315f7SKerry Stevens     }
108851d315f7SKerry Stevens     else {
108951d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
109051d315f7SKerry Stevens     }
109151d315f7SKerry Stevens     if(iterr!=0) {
109251d315f7SKerry Stevens       ithreaderr = 1;
109351d315f7SKerry Stevens     }
109451d315f7SKerry Stevens     if(PetscThreadGo) {
1095ba61063dSBarry Smith       /* reset job, get ready for more */
109651d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
109751d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1098ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
109951d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
110051d315f7SKerry Stevens     }
110151d315f7SKerry Stevens   }
110251d315f7SKerry Stevens   return NULL;
110351d315f7SKerry Stevens }
110451d315f7SKerry Stevens 
110551d315f7SKerry Stevens #undef __FUNCT__
110651d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
110751d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
110851d315f7SKerry Stevens   PetscInt i,ierr;
110951d315f7SKerry Stevens   int status;
111051d315f7SKerry Stevens 
111151d315f7SKerry Stevens   if(PetscUseThreadPool) {
111251d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
111351d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
111451d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
111551d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
111651d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
111751d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
111851d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
111951d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
112051d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112151d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112251d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1123ba61063dSBarry Smith     /* initialize job structure */
112451d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
112551d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
112651d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
112751d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
112851d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
112951d315f7SKerry Stevens     }
113051d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
113151d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
113251d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
113351d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
113451d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
113551d315f7SKerry Stevens     }
113651d315f7SKerry Stevens     job_main.pfunc = NULL;
113751d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
113851d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1139ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
114051d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1141ba61063dSBarry Smith     /* create threads */
114251d315f7SKerry Stevens     for(i=0; i<N; i++) {
114351d315f7SKerry Stevens       pVal[i] = i;
114451d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1145ba61063dSBarry Smith       /* error check */
114651d315f7SKerry Stevens     }
114751d315f7SKerry Stevens   }
114851d315f7SKerry Stevens   else {
114951d315f7SKerry Stevens   }
115051d315f7SKerry Stevens   return NULL;
115151d315f7SKerry Stevens }
115251d315f7SKerry Stevens 
115351d315f7SKerry Stevens #undef __FUNCT__
115451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
115551d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
115651d315f7SKerry Stevens   int i,ierr;
115751d315f7SKerry Stevens   void* jstatus;
115851d315f7SKerry Stevens 
115951d315f7SKerry Stevens   PetscFunctionBegin;
116051d315f7SKerry Stevens 
1161ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1162ba61063dSBarry Smith   /* join the threads */
116351d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1164ba61063dSBarry Smith     ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
116551d315f7SKerry Stevens   }
116651d315f7SKerry Stevens   free(PetscThreadPoint);
116751d315f7SKerry Stevens   free(arrmutex);
116851d315f7SKerry Stevens   free(arrcond1);
116951d315f7SKerry Stevens   free(arrcond2);
117051d315f7SKerry Stevens   free(arrstart);
117151d315f7SKerry Stevens   free(arrready);
117251d315f7SKerry Stevens   free(job_main.pdata);
117351d315f7SKerry Stevens   free(pVal);
1174cfcfc605SKerry Stevens 
117551d315f7SKerry Stevens   PetscFunctionReturn(0);
117651d315f7SKerry Stevens }
117751d315f7SKerry Stevens 
117851d315f7SKerry Stevens #undef __FUNCT__
117951d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
118051d315f7SKerry Stevens void MainWait_Main() {
118151d315f7SKerry Stevens   int i,ierr;
118251d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
118351d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
118451d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
118551d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
118651d315f7SKerry Stevens     }
118751d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
118851d315f7SKerry Stevens   }
118951d315f7SKerry Stevens }
119051d315f7SKerry Stevens 
119151d315f7SKerry Stevens #undef __FUNCT__
119251d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
119351d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
119451d315f7SKerry Stevens   int i,ierr;
119551d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
119636d20dc5SKerry Stevens 
1197ba61063dSBarry Smith   MainWait(); /* you know everyone is waiting to be signalled! */
119851d315f7SKerry Stevens   job_main.pfunc = pFunc;
119951d315f7SKerry Stevens   job_main.pdata = data;
120051d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1201ba61063dSBarry Smith     *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
120251d315f7SKerry Stevens   }
1203ba61063dSBarry Smith   /* tell the threads to go to work */
120451d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
120551d315f7SKerry Stevens     ierr = pthread_cond_signal(job_main.cond2array[i]);
120651d315f7SKerry Stevens   }
120751d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1208ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
120951d315f7SKerry Stevens   }
121036d20dc5SKerry Stevens 
121151d315f7SKerry Stevens   if(ithreaderr) {
121251d315f7SKerry Stevens     ijoberr = ithreaderr;
121351d315f7SKerry Stevens   }
121451d315f7SKerry Stevens   return ijoberr;
121551d315f7SKerry Stevens }
121651d315f7SKerry Stevens /****  ****/
121751d315f7SKerry Stevens 
121851d315f7SKerry Stevens /**** Chain Thread Functions ****/
121951d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
122051d315f7SKerry Stevens   PetscErrorCode iterr;
122151d315f7SKerry Stevens   int icorr,ierr;
122251d315f7SKerry Stevens   int* pId = (int*)arg;
122351d315f7SKerry Stevens   int ThreadId = *pId;
122451d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
122551d315f7SKerry Stevens   PetscBool PeeOn;
122651d315f7SKerry Stevens   cpu_set_t mset;
12279e800a48SKerry Stevens   //printf("Thread %d In Chain Thread Function\n",ThreadId);
122851d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
122951d315f7SKerry Stevens   CPU_ZERO(&mset);
123051d315f7SKerry Stevens   CPU_SET(icorr,&mset);
123151d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
123251d315f7SKerry Stevens 
123351d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
123451d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
123551d315f7SKerry Stevens   }
123651d315f7SKerry Stevens   else {
123751d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
123851d315f7SKerry Stevens   }
123951d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1240ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
124151d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
124251d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1243ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1244ba61063dSBarry Smith        upon return, has the lock */
124551d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
124651d315f7SKerry Stevens     }
124751d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1248ba61063dSBarry Smith     /* your subordinate is now ready*/
124951d315f7SKerry Stevens   }
125051d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1251ba61063dSBarry Smith   /* update your ready status */
125251d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
125351d315f7SKerry Stevens   if(ThreadId==0) {
125451d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1255ba61063dSBarry Smith     /* signal main */
125651d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
125751d315f7SKerry Stevens   }
125851d315f7SKerry Stevens   else {
1259ba61063dSBarry Smith     /* tell your boss that you're ready to work */
126051d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
126151d315f7SKerry Stevens   }
1262ba61063dSBarry Smith   /*  the while loop needs to have an exit
1263ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1264ba61063dSBarry Smith    and calling FuncFinish */
126551d315f7SKerry Stevens   while(PetscThreadGo) {
1266ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1267ba61063dSBarry Smith        waiting when you don't have to causes problems
1268ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
126951d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1270ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1271ba61063dSBarry Smith        upon return, has the lock */
127251d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
127351d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
127451d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
127551d315f7SKerry Stevens     }
127651d315f7SKerry Stevens     if(ThreadId==0) {
127751d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
127851d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
127951d315f7SKerry Stevens     }
128051d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
128151d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1282ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
128351d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
128451d315f7SKerry Stevens     }
1285ba61063dSBarry Smith     /* do your job */
128651d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
128751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
128851d315f7SKerry Stevens     }
128951d315f7SKerry Stevens     else {
129051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
129151d315f7SKerry Stevens     }
129251d315f7SKerry Stevens     if(iterr!=0) {
129351d315f7SKerry Stevens       ithreaderr = 1;
129451d315f7SKerry Stevens     }
129551d315f7SKerry Stevens     if(PetscThreadGo) {
1296ba61063dSBarry Smith       /* reset job, get ready for more */
129751d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1298ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1299ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
130051d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
130151d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1302ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1303ba61063dSBarry Smith            upon return, has the lock */
130451d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
130551d315f7SKerry Stevens         }
130651d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1307ba61063dSBarry Smith         /* your subordinate is now ready */
130851d315f7SKerry Stevens       }
130951d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
131051d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
131151d315f7SKerry Stevens       if(ThreadId==0) {
1312ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1313ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
131451d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
131551d315f7SKerry Stevens       }
131651d315f7SKerry Stevens       else {
1317ba61063dSBarry Smith         /* signal your boss before you go to sleep */
131851d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
131951d315f7SKerry Stevens       }
132051d315f7SKerry Stevens     }
132151d315f7SKerry Stevens   }
132251d315f7SKerry Stevens   return NULL;
132351d315f7SKerry Stevens }
132451d315f7SKerry Stevens 
132551d315f7SKerry Stevens #undef __FUNCT__
132651d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
132751d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
132851d315f7SKerry Stevens   PetscInt i,ierr;
132951d315f7SKerry Stevens   int status;
133051d315f7SKerry Stevens 
133151d315f7SKerry Stevens   if(PetscUseThreadPool) {
133251d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
133351d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
133451d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
133551d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
133651d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
133751d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
133851d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
133951d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
134051d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134151d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134251d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
134351d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1344ba61063dSBarry Smith     /* initialize job structure */
134551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
134651d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
134751d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
134851d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
134951d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
135051d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
135151d315f7SKerry Stevens     }
135251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
135351d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
135451d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
135551d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
135651d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
135751d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
135851d315f7SKerry Stevens     }
135951d315f7SKerry Stevens     job_chain.pfunc = NULL;
136051d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
136151d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
136251d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
136351d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1364ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
136551d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1366ba61063dSBarry Smith     /* create threads */
136751d315f7SKerry Stevens     for(i=0; i<N; i++) {
136851d315f7SKerry Stevens       pVal[i] = i;
136951d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1370ba61063dSBarry Smith       /* should check error */
137151d315f7SKerry Stevens     }
137251d315f7SKerry Stevens   }
137351d315f7SKerry Stevens   else {
137451d315f7SKerry Stevens   }
137551d315f7SKerry Stevens   return NULL;
137651d315f7SKerry Stevens }
137751d315f7SKerry Stevens 
137851d315f7SKerry Stevens 
137951d315f7SKerry Stevens #undef __FUNCT__
138051d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
138151d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
138251d315f7SKerry Stevens   int i,ierr;
138351d315f7SKerry Stevens   void* jstatus;
138451d315f7SKerry Stevens 
138551d315f7SKerry Stevens   PetscFunctionBegin;
138651d315f7SKerry Stevens 
1387ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1388ba61063dSBarry Smith   /* join the threads */
138951d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
139051d315f7SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1391ba61063dSBarry Smith     /* should check error */
139251d315f7SKerry Stevens   }
139351d315f7SKerry Stevens   free(PetscThreadPoint);
139451d315f7SKerry Stevens   free(arrmutex);
139551d315f7SKerry Stevens   free(arrcond1);
139651d315f7SKerry Stevens   free(arrcond2);
139751d315f7SKerry Stevens   free(arrstart);
139851d315f7SKerry Stevens   free(arrready);
139951d315f7SKerry Stevens   free(job_chain.pdata);
140051d315f7SKerry Stevens   free(pVal);
1401cfcfc605SKerry Stevens 
140251d315f7SKerry Stevens   PetscFunctionReturn(0);
140351d315f7SKerry Stevens }
140451d315f7SKerry Stevens 
140551d315f7SKerry Stevens #undef __FUNCT__
140651d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
140751d315f7SKerry Stevens void MainWait_Chain() {
140851d315f7SKerry Stevens   int ierr;
140951d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
141051d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
141151d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
141251d315f7SKerry Stevens   }
141351d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
141451d315f7SKerry Stevens }
141551d315f7SKerry Stevens 
141651d315f7SKerry Stevens #undef __FUNCT__
141751d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
141851d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
141951d315f7SKerry Stevens   int i,ierr;
142051d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
142136d20dc5SKerry Stevens 
142251d315f7SKerry Stevens   MainWait();
142351d315f7SKerry Stevens   job_chain.pfunc = pFunc;
142451d315f7SKerry Stevens   job_chain.pdata = data;
142551d315f7SKerry Stevens   job_chain.startJob = PETSC_TRUE;
142651d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
142751d315f7SKerry Stevens     *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
142851d315f7SKerry Stevens   }
142951d315f7SKerry Stevens   job_chain.eJobStat = JobInitiated;
143051d315f7SKerry Stevens   ierr = pthread_cond_signal(job_chain.cond2array[0]);
143151d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1432ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
143351d315f7SKerry Stevens   }
143436d20dc5SKerry Stevens 
143551d315f7SKerry Stevens   if(ithreaderr) {
143651d315f7SKerry Stevens     ijoberr = ithreaderr;
143751d315f7SKerry Stevens   }
143851d315f7SKerry Stevens   return ijoberr;
143951d315f7SKerry Stevens }
144051d315f7SKerry Stevens /****  ****/
144151d315f7SKerry Stevens 
1442ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
144351d315f7SKerry Stevens /**** True Thread Functions ****/
144451d315f7SKerry Stevens void* PetscThreadFunc_True(void* arg) {
144551d315f7SKerry Stevens   int icorr,ierr,iVal;
144651dcc849SKerry Stevens   int* pId = (int*)arg;
144751dcc849SKerry Stevens   int ThreadId = *pId;
14480ca81413SKerry Stevens   PetscErrorCode iterr;
144951d315f7SKerry Stevens   cpu_set_t mset;
14509e800a48SKerry Stevens   //printf("Thread %d In True Pool Thread Function\n",ThreadId);
145151d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
145251d315f7SKerry Stevens   CPU_ZERO(&mset);
145351d315f7SKerry Stevens   CPU_SET(icorr,&mset);
145451d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
145551d315f7SKerry Stevens 
145651d315f7SKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
145751d315f7SKerry Stevens   job_true.iNumReadyThreads++;
145851d315f7SKerry Stevens   if(job_true.iNumReadyThreads==PetscMaxThreads) {
145951dcc849SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
146051dcc849SKerry Stevens   }
1461ba61063dSBarry Smith   /*the while loop needs to have an exit
1462ba61063dSBarry Smith     the 'main' thread can terminate all the threads by performing a broadcast
1463ba61063dSBarry Smith    and calling FuncFinish */
146451dcc849SKerry Stevens   while(PetscThreadGo) {
1465ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
1466ba61063dSBarry Smith       waiting when you don't have to causes problems
1467ba61063dSBarry Smith      also need to wait if another thread sneaks in and messes with the predicate */
146851d315f7SKerry Stevens     while(job_true.startJob==PETSC_FALSE&&job_true.iNumJobThreads==0) {
1469ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1470ba61063dSBarry Smith        upon return, has the lock */
14716c9b208dSKerry Stevens       //printf("Thread %d Going to Sleep!\n",ThreadId);
147251d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
147351dcc849SKerry Stevens     }
147451d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
147551d315f7SKerry Stevens     job_true.iNumJobThreads--;
147651d315f7SKerry Stevens     job_true.iNumReadyThreads--;
147751d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
147851d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
147951d315f7SKerry Stevens     if(job_true.pdata==NULL) {
148051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
148151dcc849SKerry Stevens     }
148251dcc849SKerry Stevens     else {
148351d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
148451dcc849SKerry Stevens     }
14850ca81413SKerry Stevens     if(iterr!=0) {
14860ca81413SKerry Stevens       ithreaderr = 1;
14870ca81413SKerry Stevens     }
14886c9b208dSKerry Stevens     //printf("Thread %d Finished Job\n",ThreadId);
1489ba61063dSBarry Smith     /* the barrier is necessary BECAUSE: look at job_true.iNumReadyThreads
1490ba61063dSBarry Smith       what happens if a thread finishes before they all start? BAD!
1491ba61063dSBarry Smith      what happens if a thread finishes before any else start? BAD! */
1492ba61063dSBarry Smith     pthread_barrier_wait(job_true.pbarr); /* ensures all threads are finished */
1493ba61063dSBarry Smith     /* reset job */
149451dcc849SKerry Stevens     if(PetscThreadGo) {
149551d315f7SKerry Stevens       pthread_mutex_lock(&job_true.mutex);
149651d315f7SKerry Stevens       job_true.iNumReadyThreads++;
149751d315f7SKerry Stevens       if(job_true.iNumReadyThreads==PetscMaxThreads) {
1498ba61063dSBarry Smith 	/* signal the 'main' thread that the job is done! (only done once) */
149951dcc849SKerry Stevens 	ierr = pthread_cond_signal(&main_cond);
150051dcc849SKerry Stevens       }
150151dcc849SKerry Stevens     }
150251dcc849SKerry Stevens   }
150351dcc849SKerry Stevens   return NULL;
150451dcc849SKerry Stevens }
150551dcc849SKerry Stevens 
1506f09cb4aaSKerry Stevens #undef __FUNCT__
150751d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_True"
150851d315f7SKerry Stevens void* PetscThreadInitialize_True(PetscInt N) {
150951dcc849SKerry Stevens   PetscInt i;
151051dcc849SKerry Stevens   int status;
15110ca81413SKerry Stevens 
1512f09cb4aaSKerry Stevens   pVal = (int*)malloc(N*sizeof(int));
1513ba61063dSBarry Smith   /* allocate memory in the heap for the thread structure */
151451dcc849SKerry Stevens   PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1515ba61063dSBarry Smith   BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
151651d315f7SKerry Stevens   job_true.pdata = (void**)malloc(N*sizeof(void*));
151751dcc849SKerry Stevens   for(i=0; i<N; i++) {
1518f09cb4aaSKerry Stevens     pVal[i] = i;
1519f09cb4aaSKerry Stevens     status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1520ba61063dSBarry Smith     /* error check to ensure proper thread creation */
152151dcc849SKerry Stevens     status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1522ba61063dSBarry Smith     /* should check error */
152351dcc849SKerry Stevens   }
15246c9b208dSKerry Stevens   //printf("Finished True Thread Pool Initialization\n");
152551dcc849SKerry Stevens   return NULL;
152651dcc849SKerry Stevens }
152751dcc849SKerry Stevens 
1528f09cb4aaSKerry Stevens 
1529f09cb4aaSKerry Stevens #undef __FUNCT__
153051d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
153151d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
153251dcc849SKerry Stevens   int i,ierr;
153351dcc849SKerry Stevens   void* jstatus;
153451dcc849SKerry Stevens 
153551dcc849SKerry Stevens   PetscFunctionBegin;
1536cfcfc605SKerry Stevens 
1537ba61063dSBarry Smith   MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1538ba61063dSBarry Smith   /* join the threads */
153951dcc849SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
154051dcc849SKerry Stevens     ierr = pthread_join(PetscThreadPoint[i],&jstatus);
154151dcc849SKerry Stevens   }
154251dcc849SKerry Stevens   free(BarrPoint);
154351dcc849SKerry Stevens   free(PetscThreadPoint);
1544cfcfc605SKerry Stevens 
154551dcc849SKerry Stevens   PetscFunctionReturn(0);
154651dcc849SKerry Stevens }
154751dcc849SKerry Stevens 
1548f09cb4aaSKerry Stevens #undef __FUNCT__
154951d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
155051d315f7SKerry Stevens void MainWait_True() {
155151dcc849SKerry Stevens   int ierr;
15526c9b208dSKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
155351d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
155451d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
155551dcc849SKerry Stevens   }
155651d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
155751dcc849SKerry Stevens }
155851dcc849SKerry Stevens 
1559f09cb4aaSKerry Stevens #undef __FUNCT__
156051d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
156151d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
156251dcc849SKerry Stevens   int ierr;
15630ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
156436d20dc5SKerry Stevens 
15650ca81413SKerry Stevens   MainWait();
156651d315f7SKerry Stevens   job_true.pfunc = pFunc;
156751d315f7SKerry Stevens   job_true.pdata = data;
156851d315f7SKerry Stevens   job_true.pbarr = &BarrPoint[n];
156951d315f7SKerry Stevens   job_true.iNumJobThreads = n;
157051d315f7SKerry Stevens   job_true.startJob = PETSC_TRUE;
157151d315f7SKerry Stevens   ierr = pthread_cond_broadcast(&job_true.cond);
15720ca81413SKerry Stevens   if(pFunc!=FuncFinish) {
1573ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done */
15740ca81413SKerry Stevens   }
157536d20dc5SKerry Stevens 
15760ca81413SKerry Stevens   if(ithreaderr) {
15770ca81413SKerry Stevens     ijoberr = ithreaderr;
15780ca81413SKerry Stevens   }
15790ca81413SKerry Stevens   return ijoberr;
158051dcc849SKerry Stevens }
1581683509dcSKerry Stevens /**** NO THREAD POOL FUNCTION ****/
1582683509dcSKerry Stevens #undef __FUNCT__
1583683509dcSKerry Stevens #define __FUNCT__ "MainJob_Spawn"
1584683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void** data,PetscInt n) {
1585683509dcSKerry Stevens   PetscErrorCode ijoberr = 0;
1586683509dcSKerry Stevens 
1587683509dcSKerry Stevens   pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
1588cfcfc605SKerry Stevens   PetscThreadPoint = apThread; /* point to same place */
1589683509dcSKerry Stevens   PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1590683509dcSKerry Stevens   PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
1591683509dcSKerry Stevens   free(apThread);
1592683509dcSKerry Stevens 
1593683509dcSKerry Stevens   return ijoberr;
1594683509dcSKerry Stevens }
159551d315f7SKerry Stevens /****  ****/
1596ba61063dSBarry Smith #endif
159751dcc849SKerry Stevens 
159851dcc849SKerry Stevens void* FuncFinish(void* arg) {
159951dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16000ca81413SKerry Stevens   return(0);
160151dcc849SKerry Stevens }
1602ba61063dSBarry Smith 
1603ba61063dSBarry Smith #endif
1604