xref: /petsc/src/sys/objects/init.c (revision 36d20dc55d5999f7e3262f63e358d6f8a5a9893c)
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 */
8e5c89e4eSSatish Balay 
951d315f7SKerry Stevens #define _GNU_SOURCE
1051d315f7SKerry Stevens #include <sched.h>
11c6db04a5SJed Brown #include <petscsys.h>        /*I  "petscsys.h"   I*/
12ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD)
1351dcc849SKerry Stevens #include <pthread.h>
14ba61063dSBarry Smith #endif
15ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
1651d315f7SKerry Stevens #include <sys/sysinfo.h>
17ba61063dSBarry Smith #endif
1851d315f7SKerry Stevens #include <unistd.h>
19e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
20e5c89e4eSSatish Balay #include <stdlib.h>
21e5c89e4eSSatish Balay #endif
22e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
23e5c89e4eSSatish Balay #include <malloc.h>
24e5c89e4eSSatish Balay #endif
25555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
26555d055bSBarry Smith #include <valgrind/valgrind.h>
27555d055bSBarry Smith #endif
28555d055bSBarry Smith 
29e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
30e5c89e4eSSatish Balay /*
31e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
32e5c89e4eSSatish Balay    already started before PETSc was initialized.
33e5c89e4eSSatish Balay */
347087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
357087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
367087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
3751dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
3851dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
397087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
407087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
41ba61063dSBarry Smith 
42ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
4351dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
4451dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
45af359df3SBarry Smith #define PETSC_HAVE_PTHREAD_BARRIER
46ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
47ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
48ba61063dSBarry Smith #endif
4951d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
50f09cb4aaSKerry Stevens int*         pVal;
5151dcc849SKerry Stevens 
52ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
5351d315f7SKerry Stevens int* ThreadCoreAffinity;
5451d315f7SKerry Stevens 
55ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
5651d315f7SKerry Stevens 
5751d315f7SKerry Stevens typedef struct {
5851d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
5951d315f7SKerry Stevens   pthread_cond_t**  cond1array;
6051d315f7SKerry Stevens   pthread_cond_t** cond2array;
6151d315f7SKerry Stevens   void* (*pfunc)(void*);
6251d315f7SKerry Stevens   void** pdata;
6351d315f7SKerry Stevens   PetscBool startJob;
6451d315f7SKerry Stevens   estat eJobStat;
6551d315f7SKerry Stevens   PetscBool** arrThreadStarted;
6651d315f7SKerry Stevens   PetscBool** arrThreadReady;
6751d315f7SKerry Stevens } sjob_tree;
6851d315f7SKerry Stevens sjob_tree job_tree;
6951d315f7SKerry Stevens typedef struct {
7051d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
7151d315f7SKerry Stevens   pthread_cond_t**  cond1array;
7251d315f7SKerry Stevens   pthread_cond_t** cond2array;
7351d315f7SKerry Stevens   void* (*pfunc)(void*);
7451d315f7SKerry Stevens   void** pdata;
7551d315f7SKerry Stevens   PetscBool** arrThreadReady;
7651d315f7SKerry Stevens } sjob_main;
7751d315f7SKerry Stevens sjob_main job_main;
7851d315f7SKerry Stevens typedef struct {
7951d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8051d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8151d315f7SKerry Stevens   pthread_cond_t** cond2array;
8251d315f7SKerry Stevens   void* (*pfunc)(void*);
8351d315f7SKerry Stevens   void** pdata;
8451d315f7SKerry Stevens   PetscBool startJob;
8551d315f7SKerry Stevens   estat eJobStat;
8651d315f7SKerry Stevens   PetscBool** arrThreadStarted;
8751d315f7SKerry Stevens   PetscBool** arrThreadReady;
8851d315f7SKerry Stevens } sjob_chain;
8951d315f7SKerry Stevens sjob_chain job_chain;
90ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
9151dcc849SKerry Stevens typedef struct {
9251dcc849SKerry Stevens   pthread_mutex_t mutex;
9351dcc849SKerry Stevens   pthread_cond_t cond;
9451dcc849SKerry Stevens   void* (*pfunc)(void*);
9551dcc849SKerry Stevens   void** pdata;
9651dcc849SKerry Stevens   pthread_barrier_t* pbarr;
9751dcc849SKerry Stevens   int iNumJobThreads;
9851dcc849SKerry Stevens   int iNumReadyThreads;
9951dcc849SKerry Stevens   PetscBool startJob;
10051d315f7SKerry Stevens } sjob_true;
10151d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
102ba61063dSBarry Smith #endif
10351dcc849SKerry Stevens 
104ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
105ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
106ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
107ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
108ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
109ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
11051dcc849SKerry Stevens 
11151d315f7SKerry Stevens /* Function Pointers */
11251d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
11351d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
11451d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
11551d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
11651d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
117*36d20dc5SKerry Stevens /**** Tree Thread Pool Functions ****/
11851d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
11951d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
12051d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
12151d315f7SKerry Stevens void           MainWait_Tree(void);
12251d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
123*36d20dc5SKerry Stevens /**** Main Thread Pool Functions ****/
12451d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
12551d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
12651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
12751d315f7SKerry Stevens void           MainWait_Main(void);
12851d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
129*36d20dc5SKerry Stevens /**** Chain Thread Pool Functions ****/
13051d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
13151d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
13251d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
13351d315f7SKerry Stevens void           MainWait_Chain(void);
13451d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
135*36d20dc5SKerry Stevens /**** True Thread Pool Functions ****/
13651d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
13751d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
13851d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
13951d315f7SKerry Stevens void           MainWait_True(void);
14051d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
141*36d20dc5SKerry Stevens /**** NO Thread Pool Function  ****/
142683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void**,PetscInt);
143*36d20dc5SKerry Stevens /****  ****/
14451dcc849SKerry Stevens void* FuncFinish(void*);
1450ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1460ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
147ba61063dSBarry Smith #endif
148e5c89e4eSSatish Balay 
149e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
150e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
151e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
152e5c89e4eSSatish Balay #endif
1532c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1547087cfbeSBarry Smith MPI_Datatype   MPI_C_DOUBLE_COMPLEX;
1557087cfbeSBarry Smith MPI_Datatype   MPI_C_COMPLEX;
1562c876bd9SBarry Smith #endif
1577087cfbeSBarry Smith PetscScalar    PETSC_i;
158e5c89e4eSSatish Balay #else
1597087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
160e5c89e4eSSatish Balay #endif
161ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
162c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
163c90a1750SBarry Smith #endif
1647087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1657087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
16675567043SBarry Smith 
167e5c89e4eSSatish Balay /*
168e5c89e4eSSatish Balay      These are needed by petscbt.h
169e5c89e4eSSatish Balay */
170c6db04a5SJed Brown #include <petscbt.h>
1717087cfbeSBarry Smith char      _BT_mask = ' ';
1727087cfbeSBarry Smith char      _BT_c = ' ';
1737087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
174e5c89e4eSSatish Balay 
175e5c89e4eSSatish Balay /*
176e5c89e4eSSatish Balay        Function that is called to display all error messages
177e5c89e4eSSatish Balay */
1787087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1797087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
180238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1817087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
182238ccf28SShri Abhyankar #else
1837087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
184238ccf28SShri Abhyankar #endif
185bab1f7e6SVictor Minden /*
1868154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1878154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
188bab1f7e6SVictor Minden 
189e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
190e5c89e4eSSatish Balay /*
191e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
192e5c89e4eSSatish Balay */
193e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
194e5c89e4eSSatish Balay 
195e5c89e4eSSatish Balay #undef __FUNCT__
196f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
1977087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
198e5c89e4eSSatish Balay {
199e5c89e4eSSatish Balay   PetscErrorCode ierr;
200e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
201e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
202e5c89e4eSSatish Balay   char           version[256];
203e5c89e4eSSatish Balay 
204e5c89e4eSSatish Balay   PetscFunctionBegin;
205e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
206e5c89e4eSSatish Balay   if (!rank) {
207e5c89e4eSSatish Balay     char        arch[10];
208f56c2debSBarry Smith     int         err;
20988c29154SBarry Smith     PetscViewer viewer;
210f56c2debSBarry Smith 
211e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
212e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
213a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
214e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
215e5c89e4eSSatish Balay     if (filename) {
216e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
217e5c89e4eSSatish Balay     } else {
218e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
219e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
220e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
221e5c89e4eSSatish Balay     }
222e5c89e4eSSatish Balay 
223e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
224e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
225e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
226e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
227e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
22888c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
22988c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2306bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
231e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
232f56c2debSBarry Smith     err = fflush(*fd);
233e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
234e5c89e4eSSatish Balay   }
235e5c89e4eSSatish Balay   PetscFunctionReturn(0);
236e5c89e4eSSatish Balay }
237e5c89e4eSSatish Balay 
238e5c89e4eSSatish Balay #undef __FUNCT__
239f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2407087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
241e5c89e4eSSatish Balay {
242e5c89e4eSSatish Balay   PetscErrorCode ierr;
243e5c89e4eSSatish Balay   PetscMPIInt    rank;
244e5c89e4eSSatish Balay   char           date[64];
245f56c2debSBarry Smith   int            err;
246e5c89e4eSSatish Balay 
247e5c89e4eSSatish Balay   PetscFunctionBegin;
248e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
249e5c89e4eSSatish Balay   if (!rank) {
250e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
251e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
252e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
253e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
254f56c2debSBarry Smith     err = fflush(*fd);
255e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
256f56c2debSBarry Smith     err = fclose(*fd);
257e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
258e5c89e4eSSatish Balay   }
259e5c89e4eSSatish Balay   PetscFunctionReturn(0);
260e5c89e4eSSatish Balay }
261e5c89e4eSSatish Balay 
262e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
263e5c89e4eSSatish Balay 
264e5c89e4eSSatish Balay /*
265e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
266e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
267e5c89e4eSSatish Balay 
268e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2693c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2703c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
271e5c89e4eSSatish Balay */
272e5c89e4eSSatish Balay 
273e5c89e4eSSatish Balay #undef __FUNCT__
274e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
275e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
276e5c89e4eSSatish Balay {
277e5c89e4eSSatish Balay   PetscFunctionBegin;
2783c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
279e5c89e4eSSatish Balay   abort();
280e5c89e4eSSatish Balay }
281e5c89e4eSSatish Balay 
282e5c89e4eSSatish Balay #undef __FUNCT__
283e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
284e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
285e5c89e4eSSatish Balay {
286e5c89e4eSSatish Balay   PetscErrorCode ierr;
287e5c89e4eSSatish Balay 
288e5c89e4eSSatish Balay   PetscFunctionBegin;
2893c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
290e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
291e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
2923c311c98SBarry Smith     MPI_Abort(*comm,*flag);
293e5c89e4eSSatish Balay   }
294e5c89e4eSSatish Balay }
295e5c89e4eSSatish Balay 
296e5c89e4eSSatish Balay #undef __FUNCT__
297e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
298e5c89e4eSSatish Balay /*@C
299e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
300e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
301e5c89e4eSSatish Balay 
302e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
303e5c89e4eSSatish Balay 
304e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
305e5c89e4eSSatish Balay 
306e5c89e4eSSatish Balay    Level: advanced
307e5c89e4eSSatish Balay 
308e5c89e4eSSatish Balay    Note:
309e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
310e5c89e4eSSatish Balay 
31188c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
312e5c89e4eSSatish Balay @*/
3137087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
314e5c89e4eSSatish Balay {
315e5c89e4eSSatish Balay   PetscFunctionBegin;
316e5c89e4eSSatish Balay   PetscFinalize();
317e5c89e4eSSatish Balay   exit(0);
318e5c89e4eSSatish Balay   return 0;
319e5c89e4eSSatish Balay }
320e5c89e4eSSatish Balay 
321ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
32209573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
323ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
324e5c89e4eSSatish Balay static char       emacsmachinename[256];
325e5c89e4eSSatish Balay 
326e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
327e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
328e5c89e4eSSatish Balay 
329e5c89e4eSSatish Balay #undef __FUNCT__
330e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
331e5c89e4eSSatish Balay /*@C
332e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
333e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
334e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
335e5c89e4eSSatish Balay 
336e5c89e4eSSatish Balay    Input Parameter:
337e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
338da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
339e5c89e4eSSatish Balay 
340e5c89e4eSSatish Balay    Level: developer
341e5c89e4eSSatish Balay 
342e5c89e4eSSatish Balay    Concepts: package help message
343e5c89e4eSSatish Balay 
344e5c89e4eSSatish Balay @*/
3457087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
346e5c89e4eSSatish Balay {
347e5c89e4eSSatish Balay   PetscFunctionBegin;
348e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
349e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
350e5c89e4eSSatish Balay   PetscFunctionReturn(0);
351e5c89e4eSSatish Balay }
352e5c89e4eSSatish Balay 
353e5c89e4eSSatish Balay #undef __FUNCT__
354e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3557087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
356e5c89e4eSSatish Balay {
357e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
358e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
359ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
360e5c89e4eSSatish Balay   PetscErrorCode ierr;
361a6d0e24fSJed Brown   PetscReal      si;
362e5c89e4eSSatish Balay   int            i;
363e5c89e4eSSatish Balay   PetscMPIInt    rank;
364e5c89e4eSSatish Balay   char           version[256];
365e5c89e4eSSatish Balay 
366e5c89e4eSSatish Balay   PetscFunctionBegin;
367e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
368e5c89e4eSSatish Balay 
369e5c89e4eSSatish Balay   /*
370e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
371e5c89e4eSSatish Balay   */
3728bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
37381b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
374acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
375e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
376555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
377555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
378555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
379555d055bSBarry Smith #endif
380e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
381555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
382555d055bSBarry Smith     }
383555d055bSBarry Smith #endif
384e5c89e4eSSatish Balay   }
385e5c89e4eSSatish Balay #else
386acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
387acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
388e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
389e5c89e4eSSatish Balay #endif
390e5c89e4eSSatish Balay   if (flg3) {
391e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
392e5c89e4eSSatish Balay   }
39390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
394acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
395e5c89e4eSSatish Balay   if (flg1) {
396e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
397e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
398e5c89e4eSSatish Balay   }
399e5c89e4eSSatish Balay 
40090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
401acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4027783f70dSSatish Balay   if (!flg1) {
40390d69ab7SBarry Smith     flg1 = PETSC_FALSE;
404acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4057783f70dSSatish Balay   }
406e5c89e4eSSatish Balay   if (flg1) {
407e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
408e5c89e4eSSatish Balay   }
409e5c89e4eSSatish Balay 
410e5c89e4eSSatish Balay   /*
411e5c89e4eSSatish Balay       Set the display variable for graphics
412e5c89e4eSSatish Balay   */
413e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
414e5c89e4eSSatish Balay 
41551dcc849SKerry Stevens 
41651dcc849SKerry Stevens   /*
417e5c89e4eSSatish Balay       Print the PETSc version information
418e5c89e4eSSatish Balay   */
419e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
420e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
421e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
422e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
423e5c89e4eSSatish Balay 
424e5c89e4eSSatish Balay     /*
425e5c89e4eSSatish Balay        Print "higher-level" package version message
426e5c89e4eSSatish Balay     */
427e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
428e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
429e5c89e4eSSatish Balay     }
430e5c89e4eSSatish Balay 
431a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
432e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
433e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
434e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
435e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
436e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
43784e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
438e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
439e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
440e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
441e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
442e5c89e4eSSatish Balay   }
443e5c89e4eSSatish Balay 
444e5c89e4eSSatish Balay   /*
445e5c89e4eSSatish Balay        Print "higher-level" package help message
446e5c89e4eSSatish Balay   */
447e5c89e4eSSatish Balay   if (flg3){
448e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
449e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
450e5c89e4eSSatish Balay     }
451e5c89e4eSSatish Balay   }
452e5c89e4eSSatish Balay 
453e5c89e4eSSatish Balay   /*
454e5c89e4eSSatish Balay       Setup the error handling
455e5c89e4eSSatish Balay   */
45690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
457acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
458cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
45990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
460acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
461cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
46290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
463acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
464e5c89e4eSSatish Balay   if (flg1) {
465e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
466e5c89e4eSSatish Balay   }
46790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
468acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
469cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
47096cc47afSJed Brown   flg1 = PETSC_FALSE;
471acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
47296cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
473e5c89e4eSSatish Balay 
474e5c89e4eSSatish Balay   /*
475e5c89e4eSSatish Balay       Setup debugger information
476e5c89e4eSSatish Balay   */
477e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
478e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
479e5c89e4eSSatish Balay   if (flg1) {
480e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
481e5c89e4eSSatish Balay 
482e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
483e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
484e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
485e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
486e5c89e4eSSatish Balay   }
4875e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
4885e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
489e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
490e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
491e5c89e4eSSatish Balay   if (flg1 || flg2) {
492e5c89e4eSSatish Balay     PetscMPIInt    size;
493e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
494e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
495e5c89e4eSSatish Balay     /*
496e5c89e4eSSatish Balay        we have to make sure that all processors have opened
497e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
498e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
499e5c89e4eSSatish Balay        and kill the program.
500e5c89e4eSSatish Balay     */
501e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
502e5c89e4eSSatish Balay     if (size > 2) {
503533163c2SBarry Smith       PetscMPIInt dummy = 0;
504e5c89e4eSSatish Balay       MPI_Status  status;
505e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
506e5c89e4eSSatish Balay         if (rank != i) {
507e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
508e5c89e4eSSatish Balay         }
509e5c89e4eSSatish Balay       }
510e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
511e5c89e4eSSatish Balay         if (rank != i) {
512e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
513e5c89e4eSSatish Balay         }
514e5c89e4eSSatish Balay       }
515e5c89e4eSSatish Balay     }
516e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
517e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
518e5c89e4eSSatish Balay     lsize = size;
519e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
520e5c89e4eSSatish Balay     if (flag) {
521e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
522e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
523e5c89e4eSSatish Balay       }
524e5c89e4eSSatish Balay     }
525e5c89e4eSSatish Balay     if (!flag) {
526e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
527e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
528e5c89e4eSSatish Balay       if (flg1) {
529e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
530e5c89e4eSSatish Balay       } else {
531e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
532e5c89e4eSSatish Balay       }
533e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
534e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
535e5c89e4eSSatish Balay     }
536e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
537e5c89e4eSSatish Balay   }
538e5c89e4eSSatish Balay 
539e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
540cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
541e5c89e4eSSatish Balay 
54293ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
54322b84c2fSbcordonn   /*
54422b84c2fSbcordonn     Activates new sockets for zope if needed
54522b84c2fSbcordonn   */
54684ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
547d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5486dc8fec2Sbcordonn   if (flgz){
54922b84c2fSbcordonn     int  sockfd;
550f1384234SBarry Smith     char hostname[256];
55122b84c2fSbcordonn     char username[256];
5526dc8fec2Sbcordonn     int  remoteport = 9999;
5539c4c166aSBarry Smith 
55484ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
55584ab5442Sbcordonn     if (!hostname[0]){
5569c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5579c4c166aSBarry Smith     }
55822b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5599c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
56022b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
56122b84c2fSbcordonn     if (flgzout){
56222b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
563606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5646dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5659c4c166aSBarry Smith     } else {
566d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
567d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5689c4c166aSBarry Smith     }
5699c4c166aSBarry Smith   }
57093ba235fSBarry Smith #endif
571ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
572ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
573ffc871a5SBarry Smith   if (flgz){
574ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
575ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
576ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
577ffc871a5SBarry Smith   }
578ffc871a5SBarry Smith #endif
5796dc8fec2Sbcordonn 
580e5c89e4eSSatish Balay   /*
581e5c89e4eSSatish Balay         Setup profiling and logging
582e5c89e4eSSatish Balay   */
5836cf91177SBarry Smith #if defined (PETSC_USE_INFO)
5848bb29257SSatish Balay   {
585e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
5866cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
5878bb29257SSatish Balay     if (flg1 && logname[0]) {
588fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
5898bb29257SSatish Balay     } else if (flg1) {
590fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
591e5c89e4eSSatish Balay     }
592e5c89e4eSSatish Balay   }
593865f6aa8SSatish Balay #endif
594865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
595865f6aa8SSatish Balay   mname[0] = 0;
596f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
597865f6aa8SSatish Balay   if (flg1) {
598865f6aa8SSatish Balay     if (mname[0]) {
599f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
600865f6aa8SSatish Balay     } else {
601f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
602865f6aa8SSatish Balay     }
603865f6aa8SSatish Balay   }
604e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
60590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
606fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
607e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
608e5c89e4eSSatish Balay #endif
60990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
61090d69ab7SBarry Smith   flg2 = PETSC_FALSE;
61190d69ab7SBarry Smith   flg3 = PETSC_FALSE;
612acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
613acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
614d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6159f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
616e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6179f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
618e5c89e4eSSatish Balay 
619e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
620e5c89e4eSSatish Balay   if (flg1) {
621e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
622e5c89e4eSSatish Balay     FILE *file;
623e5c89e4eSSatish Balay     if (mname[0]) {
624e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
625e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
626e5c89e4eSSatish Balay       file = fopen(fname,"w");
627f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
628e5c89e4eSSatish Balay     } else {
629da9f1d6bSBarry Smith       file = PETSC_STDOUT;
630e5c89e4eSSatish Balay     }
631e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
632e5c89e4eSSatish Balay   }
633e5c89e4eSSatish Balay #endif
634e5c89e4eSSatish Balay 
635e5c89e4eSSatish Balay   /*
636e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
637e5c89e4eSSatish Balay   */
63863d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
639e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
640e5c89e4eSSatish Balay #endif
641e5c89e4eSSatish Balay 
642acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
643e5c89e4eSSatish Balay 
644af359df3SBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
645af359df3SBarry Smith   /*
646af359df3SBarry Smith       Determine whether user specified maximum number of threads
647af359df3SBarry Smith    */
648af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
649af359df3SBarry Smith 
650b154f58aSKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-main",&flg1);CHKERRQ(ierr);
651b154f58aSKerry Stevens   if(flg1) {
652b154f58aSKerry Stevens     cpu_set_t mset;
653b154f58aSKerry Stevens     int icorr,ncorr = get_nprocs();
654b154f58aSKerry Stevens     ierr = PetscOptionsGetInt(PETSC_NULL,"-main",&icorr,PETSC_NULL);CHKERRQ(ierr);
655b154f58aSKerry Stevens     CPU_ZERO(&mset);
656b154f58aSKerry Stevens     CPU_SET(icorr%ncorr,&mset);
657b154f58aSKerry Stevens     sched_setaffinity(0,sizeof(cpu_set_t),&mset);
658b154f58aSKerry Stevens   }
659b154f58aSKerry Stevens 
660af359df3SBarry Smith   /*
661af359df3SBarry Smith       Determine whether to use thread pool
662af359df3SBarry Smith    */
663af359df3SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
664af359df3SBarry Smith   if (flg1) {
665af359df3SBarry Smith     PetscUseThreadPool = PETSC_TRUE;
666af359df3SBarry Smith     PetscInt N_CORES = get_nprocs();
667af359df3SBarry Smith     ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
668af359df3SBarry Smith     char tstr[9];
669af359df3SBarry Smith     char tbuf[2];
670af359df3SBarry Smith     strcpy(tstr,"-thread");
671af359df3SBarry Smith     for(i=0;i<PetscMaxThreads;i++) {
672af359df3SBarry Smith       ThreadCoreAffinity[i] = i;
673af359df3SBarry Smith       sprintf(tbuf,"%d",i);
674af359df3SBarry Smith       strcat(tstr,tbuf);
675af359df3SBarry Smith       ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
676af359df3SBarry Smith       if(flg1) {
677af359df3SBarry Smith         ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
678af359df3SBarry Smith         ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
679af359df3SBarry Smith       }
680af359df3SBarry Smith       tstr[7] = '\0';
681af359df3SBarry Smith     }
682af359df3SBarry Smith     /* get the thread pool type */
683af359df3SBarry Smith     PetscInt ipool = 0;
684af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
685af359df3SBarry Smith 
686af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
687af359df3SBarry Smith     switch(ipool) {
688af359df3SBarry Smith     case 1:
689af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
690af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
691af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
692af359df3SBarry Smith       MainWait              = &MainWait_Tree;
693af359df3SBarry Smith       MainJob               = &MainJob_Tree;
694af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
695af359df3SBarry Smith       break;
696af359df3SBarry Smith     case 2:
697af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
698af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
699af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
700af359df3SBarry Smith       MainWait              = &MainWait_Main;
701af359df3SBarry Smith       MainJob               = &MainJob_Main;
702af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
703af359df3SBarry Smith       break;
704af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
705af359df3SBarry Smith     case 3:
706af359df3SBarry Smith #else
707af359df3SBarry Smith     default:
708af359df3SBarry Smith #endif
709af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
710af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
711af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
712af359df3SBarry Smith       MainWait              = &MainWait_Chain;
713af359df3SBarry Smith       MainJob               = &MainJob_Chain;
714af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
715af359df3SBarry Smith       break;
716af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
717af359df3SBarry Smith     default:
718af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
719af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
720af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
721af359df3SBarry Smith       MainWait              = &MainWait_True;
722af359df3SBarry Smith       MainJob               = &MainJob_True;
723af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
724af359df3SBarry Smith       break;
725af359df3SBarry Smith #endif
726af359df3SBarry Smith     }
727af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
728683509dcSKerry Stevens   } else {
729683509dcSKerry Stevens     //need to define these in the case on 'no threads' or 'thread create/destroy'
730683509dcSKerry Stevens     //could take any of the above versions
731683509dcSKerry Stevens     MainJob               = &MainJob_Spawn;
732af359df3SBarry Smith   }
733af359df3SBarry Smith #endif
734e5c89e4eSSatish Balay   /*
735e5c89e4eSSatish Balay        Print basic help message
736e5c89e4eSSatish Balay   */
737e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
738e5c89e4eSSatish Balay   if (flg1) {
739e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
740301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
741301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
742301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
743e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
744e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
745e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
746e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
747e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
748e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
749e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
750e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
751e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
752e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
753e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
754e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
755e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
756e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
757e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
758e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
759e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
760e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
761e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7624161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7634161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
764e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
765e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
766e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
767e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
768e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
769e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
770a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
771e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
77240ab9619SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -server <port>: Run PETSc webserver (default port is 8080) see PetscWebServe()\n");CHKERRQ(ierr);
773e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
774e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
7757ef452c0SMatthew G Knepley     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary _summary_python]: logging objects and events\n");CHKERRQ(ierr);
776e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
777e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
778e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
779e5c89e4eSSatish Balay #endif
7806cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
781e5c89e4eSSatish Balay #endif
782e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
783e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
784e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
785e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
786e5c89e4eSSatish Balay   }
787e5c89e4eSSatish Balay 
788a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
789e5c89e4eSSatish Balay   if (flg1) {
790e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
791e5c89e4eSSatish Balay   }
792e5c89e4eSSatish Balay 
7936cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
794e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
795e5c89e4eSSatish Balay   if (f) {
7966cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
797e5c89e4eSSatish Balay   }
798827f890bSBarry Smith 
7998154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
800c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
80173113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
80273113deaSBarry Smith   else flg1 = PETSC_FALSE;
8038154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
8048154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
805bab1f7e6SVictor Minden #endif
806192daf7cSBarry Smith 
807e5c89e4eSSatish Balay   PetscFunctionReturn(0);
808e5c89e4eSSatish Balay }
809df413903SBarry Smith 
810ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
811ba61063dSBarry Smith 
81251d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
81351d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
81451d315f7SKerry Stevens   PetscErrorCode iterr;
81551d315f7SKerry Stevens   int icorr,ierr;
81651d315f7SKerry Stevens   int* pId = (int*)arg;
81751d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
81851d315f7SKerry Stevens   PetscBool PeeOn;
81951d315f7SKerry Stevens   cpu_set_t mset;
8209e800a48SKerry Stevens   //printf("Thread %d In Tree Thread Function\n",ThreadId);
82151d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
82251d315f7SKerry Stevens   CPU_ZERO(&mset);
82351d315f7SKerry Stevens   CPU_SET(icorr,&mset);
82451d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
82551d315f7SKerry Stevens 
82651d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
82751d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
82851d315f7SKerry Stevens   }
82951d315f7SKerry Stevens   else {
83051d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
83151d315f7SKerry Stevens   }
83251d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
833ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
83451d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
83551d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
83651d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
83751d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
83851d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
839ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
840ba61063dSBarry Smith            upon return, has the lock */
84151d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
84251d315f7SKerry Stevens         }
84351d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
84451d315f7SKerry Stevens       }
84551d315f7SKerry Stevens     }
846ba61063dSBarry Smith     /* your subordinates are now ready */
84751d315f7SKerry Stevens   }
84851d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
849ba61063dSBarry Smith   /* update your ready status */
85051d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
85151d315f7SKerry Stevens   if(ThreadId==0) {
85251d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
853ba61063dSBarry Smith     /* ignal main */
85451d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
85551d315f7SKerry Stevens   }
85651d315f7SKerry Stevens   else {
857ba61063dSBarry Smith     /* tell your boss that you're ready to work */
85851d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
85951d315f7SKerry Stevens   }
860ba61063dSBarry Smith   /* the while loop needs to have an exit
861ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
862ba61063dSBarry Smith    and calling FuncFinish */
86351d315f7SKerry Stevens   while(PetscThreadGo) {
864ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
865ba61063dSBarry Smith       waiting when you don't have to causes problems
866ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
86751d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
868ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
869ba61063dSBarry Smith        upon return, has the lock */
87051d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
87151d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
87251d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
87351d315f7SKerry Stevens     }
87451d315f7SKerry Stevens     if(ThreadId==0) {
87551d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
87651d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
87751d315f7SKerry Stevens     }
87851d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
87951d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
880ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
88151d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
88251d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
88351d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
88451d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
88551d315f7SKerry Stevens         }
88651d315f7SKerry Stevens       }
88751d315f7SKerry Stevens     }
888ba61063dSBarry Smith     /* do your job */
88951d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
89051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
89151d315f7SKerry Stevens     }
89251d315f7SKerry Stevens     else {
89351d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
89451d315f7SKerry Stevens     }
89551d315f7SKerry Stevens     if(iterr!=0) {
89651d315f7SKerry Stevens       ithreaderr = 1;
89751d315f7SKerry Stevens     }
89851d315f7SKerry Stevens     if(PetscThreadGo) {
899ba61063dSBarry Smith       /* reset job, get ready for more */
90051d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
901ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
902ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
90351d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
90451d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
90551d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
90651d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
90751d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
908ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
909ba61063dSBarry Smith                upon return, has the lock */
91051d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
91151d315f7SKerry Stevens             }
91251d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
91351d315f7SKerry Stevens           }
91451d315f7SKerry Stevens 	}
915ba61063dSBarry Smith         /* your subordinates are now ready */
91651d315f7SKerry Stevens       }
91751d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
91851d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
91951d315f7SKerry Stevens       if(ThreadId==0) {
920ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
921ba61063dSBarry Smith         /* root thread signals 'main' */
92251d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
92351d315f7SKerry Stevens       }
92451d315f7SKerry Stevens       else {
925ba61063dSBarry Smith         /* signal your boss before you go to sleep */
92651d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
92751d315f7SKerry Stevens       }
92851d315f7SKerry Stevens     }
92951d315f7SKerry Stevens   }
93051d315f7SKerry Stevens   return NULL;
93151d315f7SKerry Stevens }
93251d315f7SKerry Stevens 
93351d315f7SKerry Stevens #undef __FUNCT__
93451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
93551d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
93651d315f7SKerry Stevens   PetscInt i,ierr;
93751d315f7SKerry Stevens   int status;
93851d315f7SKerry Stevens 
93951d315f7SKerry Stevens   if(PetscUseThreadPool) {
94051d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
94151d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
94251d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
94351d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
94451d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
94551d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
94651d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
94751d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
94851d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
94951d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95051d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
95151d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
952ba61063dSBarry Smith     /* initialize job structure */
95351d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
95451d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
95551d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
95651d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
95751d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
95851d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
95951d315f7SKerry Stevens     }
96051d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
96151d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
96251d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
96351d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
96451d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
96551d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
96651d315f7SKerry Stevens     }
96751d315f7SKerry Stevens     job_tree.pfunc = NULL;
96851d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
96951d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
97051d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
97151d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
972ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
97351d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
974ba61063dSBarry Smith     /* create threads */
97551d315f7SKerry Stevens     for(i=0; i<N; i++) {
97651d315f7SKerry Stevens       pVal[i] = i;
97751d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
978ba61063dSBarry Smith       /* should check status */
97951d315f7SKerry Stevens     }
98051d315f7SKerry Stevens   }
98151d315f7SKerry Stevens   return NULL;
98251d315f7SKerry Stevens }
98351d315f7SKerry Stevens 
98451d315f7SKerry Stevens #undef __FUNCT__
98551d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
98651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
98751d315f7SKerry Stevens   int i,ierr;
98851d315f7SKerry Stevens   void* jstatus;
98951d315f7SKerry Stevens 
99051d315f7SKerry Stevens   PetscFunctionBegin;
99151d315f7SKerry Stevens 
99251d315f7SKerry Stevens   if(PetscUseThreadPool) {
993ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
994ba61063dSBarry Smith     /* join the threads */
99551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
99651d315f7SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
997ba61063dSBarry Smith       /* do error checking*/
99851d315f7SKerry Stevens     }
99951d315f7SKerry Stevens     free(PetscThreadPoint);
100051d315f7SKerry Stevens     free(arrmutex);
100151d315f7SKerry Stevens     free(arrcond1);
100251d315f7SKerry Stevens     free(arrcond2);
100351d315f7SKerry Stevens     free(arrstart);
100451d315f7SKerry Stevens     free(arrready);
100551d315f7SKerry Stevens     free(job_tree.pdata);
100651d315f7SKerry Stevens     free(pVal);
100751d315f7SKerry Stevens   }
100851d315f7SKerry Stevens   else {
100951d315f7SKerry Stevens   }
101051d315f7SKerry Stevens   PetscFunctionReturn(0);
101151d315f7SKerry Stevens }
101251d315f7SKerry Stevens 
101351d315f7SKerry Stevens #undef __FUNCT__
101451d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
101551d315f7SKerry Stevens void MainWait_Tree() {
101651d315f7SKerry Stevens   int ierr;
101751d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
101851d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
101951d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
102051d315f7SKerry Stevens   }
102151d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
102251d315f7SKerry Stevens }
102351d315f7SKerry Stevens 
102451d315f7SKerry Stevens #undef __FUNCT__
102551d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
102651d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
102751d315f7SKerry Stevens   int i,ierr;
102851d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
1029*36d20dc5SKerry Stevens 
103051d315f7SKerry Stevens   MainWait();
103151d315f7SKerry Stevens   job_tree.pfunc = pFunc;
103251d315f7SKerry Stevens   job_tree.pdata = data;
103351d315f7SKerry Stevens   job_tree.startJob = PETSC_TRUE;
103451d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
103551d315f7SKerry Stevens     *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
103651d315f7SKerry Stevens   }
103751d315f7SKerry Stevens   job_tree.eJobStat = JobInitiated;
103851d315f7SKerry Stevens   ierr = pthread_cond_signal(job_tree.cond2array[0]);
103951d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1040ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
104151d315f7SKerry Stevens   }
1042*36d20dc5SKerry Stevens 
104351d315f7SKerry Stevens   if(ithreaderr) {
104451d315f7SKerry Stevens     ijoberr = ithreaderr;
104551d315f7SKerry Stevens   }
104651d315f7SKerry Stevens   return ijoberr;
104751d315f7SKerry Stevens }
104851d315f7SKerry Stevens /****  ****/
104951d315f7SKerry Stevens 
105051d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
105151d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
105251d315f7SKerry Stevens   PetscErrorCode iterr;
105351d315f7SKerry Stevens   int icorr,ierr;
105451d315f7SKerry Stevens   int* pId = (int*)arg;
105551d315f7SKerry Stevens   int ThreadId = *pId;
105651d315f7SKerry Stevens   cpu_set_t mset;
10579e800a48SKerry Stevens   //printf("Thread %d In Main Thread Function\n",ThreadId);
105851d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
105951d315f7SKerry Stevens   CPU_ZERO(&mset);
106051d315f7SKerry Stevens   CPU_SET(icorr,&mset);
106151d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
106251d315f7SKerry Stevens 
106351d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1064ba61063dSBarry Smith   /* update your ready status */
106551d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1066ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
106751d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
106851d315f7SKerry Stevens 
1069ba61063dSBarry Smith   /* the while loop needs to have an exit
1070ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1071ba61063dSBarry Smith      and calling FuncFinish */
107251d315f7SKerry Stevens   while(PetscThreadGo) {
1073ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1074ba61063dSBarry Smith        waiting when you don't have to causes problems
1075ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
107651d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1077ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1078ba61063dSBarry Smith        upon return, has the lock */
107951d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1080ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
108151d315f7SKerry Stevens     }
108251d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
108351d315f7SKerry Stevens     if(job_main.pdata==NULL) {
108451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
108551d315f7SKerry Stevens     }
108651d315f7SKerry Stevens     else {
108751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
108851d315f7SKerry Stevens     }
108951d315f7SKerry Stevens     if(iterr!=0) {
109051d315f7SKerry Stevens       ithreaderr = 1;
109151d315f7SKerry Stevens     }
109251d315f7SKerry Stevens     if(PetscThreadGo) {
1093ba61063dSBarry Smith       /* reset job, get ready for more */
109451d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
109551d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1096ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
109751d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
109851d315f7SKerry Stevens     }
109951d315f7SKerry Stevens   }
110051d315f7SKerry Stevens   return NULL;
110151d315f7SKerry Stevens }
110251d315f7SKerry Stevens 
110351d315f7SKerry Stevens #undef __FUNCT__
110451d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
110551d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
110651d315f7SKerry Stevens   PetscInt i,ierr;
110751d315f7SKerry Stevens   int status;
110851d315f7SKerry Stevens 
110951d315f7SKerry Stevens   if(PetscUseThreadPool) {
111051d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
111151d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
111251d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
111351d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
111451d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
111551d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
111651d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
111751d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
111851d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
111951d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112051d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1121ba61063dSBarry Smith     /* initialize job structure */
112251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
112351d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
112451d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
112551d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
112651d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
112751d315f7SKerry Stevens     }
112851d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
112951d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
113051d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
113151d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
113251d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
113351d315f7SKerry Stevens     }
113451d315f7SKerry Stevens     job_main.pfunc = NULL;
113551d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
113651d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1137ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
113851d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1139ba61063dSBarry Smith     /* create threads */
114051d315f7SKerry Stevens     for(i=0; i<N; i++) {
114151d315f7SKerry Stevens       pVal[i] = i;
114251d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1143ba61063dSBarry Smith       /* error check */
114451d315f7SKerry Stevens     }
114551d315f7SKerry Stevens   }
114651d315f7SKerry Stevens   else {
114751d315f7SKerry Stevens   }
114851d315f7SKerry Stevens   return NULL;
114951d315f7SKerry Stevens }
115051d315f7SKerry Stevens 
115151d315f7SKerry Stevens #undef __FUNCT__
115251d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
115351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
115451d315f7SKerry Stevens   int i,ierr;
115551d315f7SKerry Stevens   void* jstatus;
115651d315f7SKerry Stevens 
115751d315f7SKerry Stevens   PetscFunctionBegin;
115851d315f7SKerry Stevens 
115951d315f7SKerry Stevens   if(PetscUseThreadPool) {
1160ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1161ba61063dSBarry Smith     /* join the threads */
116251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
1163ba61063dSBarry Smith       ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
116451d315f7SKerry Stevens     }
116551d315f7SKerry Stevens     free(PetscThreadPoint);
116651d315f7SKerry Stevens     free(arrmutex);
116751d315f7SKerry Stevens     free(arrcond1);
116851d315f7SKerry Stevens     free(arrcond2);
116951d315f7SKerry Stevens     free(arrstart);
117051d315f7SKerry Stevens     free(arrready);
117151d315f7SKerry Stevens     free(job_main.pdata);
117251d315f7SKerry Stevens     free(pVal);
117351d315f7SKerry Stevens   }
117451d315f7SKerry Stevens   PetscFunctionReturn(0);
117551d315f7SKerry Stevens }
117651d315f7SKerry Stevens 
117751d315f7SKerry Stevens #undef __FUNCT__
117851d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
117951d315f7SKerry Stevens void MainWait_Main() {
118051d315f7SKerry Stevens   int i,ierr;
118151d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
118251d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
118351d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
118451d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
118551d315f7SKerry Stevens     }
118651d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
118751d315f7SKerry Stevens   }
118851d315f7SKerry Stevens }
118951d315f7SKerry Stevens 
119051d315f7SKerry Stevens #undef __FUNCT__
119151d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
119251d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
119351d315f7SKerry Stevens   int i,ierr;
119451d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
1195*36d20dc5SKerry Stevens 
1196ba61063dSBarry Smith   MainWait(); /* you know everyone is waiting to be signalled! */
119751d315f7SKerry Stevens   job_main.pfunc = pFunc;
119851d315f7SKerry Stevens   job_main.pdata = data;
119951d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1200ba61063dSBarry Smith     *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
120151d315f7SKerry Stevens   }
1202ba61063dSBarry Smith   /* tell the threads to go to work */
120351d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
120451d315f7SKerry Stevens     ierr = pthread_cond_signal(job_main.cond2array[i]);
120551d315f7SKerry Stevens   }
120651d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1207ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
120851d315f7SKerry Stevens   }
1209*36d20dc5SKerry Stevens 
121051d315f7SKerry Stevens   if(ithreaderr) {
121151d315f7SKerry Stevens     ijoberr = ithreaderr;
121251d315f7SKerry Stevens   }
121351d315f7SKerry Stevens   return ijoberr;
121451d315f7SKerry Stevens }
121551d315f7SKerry Stevens /****  ****/
121651d315f7SKerry Stevens 
121751d315f7SKerry Stevens /**** Chain Thread Functions ****/
121851d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
121951d315f7SKerry Stevens   PetscErrorCode iterr;
122051d315f7SKerry Stevens   int icorr,ierr;
122151d315f7SKerry Stevens   int* pId = (int*)arg;
122251d315f7SKerry Stevens   int ThreadId = *pId;
122351d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
122451d315f7SKerry Stevens   PetscBool PeeOn;
122551d315f7SKerry Stevens   cpu_set_t mset;
12269e800a48SKerry Stevens   //printf("Thread %d In Chain Thread Function\n",ThreadId);
122751d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
122851d315f7SKerry Stevens   CPU_ZERO(&mset);
122951d315f7SKerry Stevens   CPU_SET(icorr,&mset);
123051d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
123151d315f7SKerry Stevens 
123251d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
123351d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
123451d315f7SKerry Stevens   }
123551d315f7SKerry Stevens   else {
123651d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
123751d315f7SKerry Stevens   }
123851d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1239ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
124051d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
124151d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1242ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1243ba61063dSBarry Smith        upon return, has the lock */
124451d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
124551d315f7SKerry Stevens     }
124651d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1247ba61063dSBarry Smith     /* your subordinate is now ready*/
124851d315f7SKerry Stevens   }
124951d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1250ba61063dSBarry Smith   /* update your ready status */
125151d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
125251d315f7SKerry Stevens   if(ThreadId==0) {
125351d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1254ba61063dSBarry Smith     /* signal main */
125551d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
125651d315f7SKerry Stevens   }
125751d315f7SKerry Stevens   else {
1258ba61063dSBarry Smith     /* tell your boss that you're ready to work */
125951d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
126051d315f7SKerry Stevens   }
1261ba61063dSBarry Smith   /*  the while loop needs to have an exit
1262ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1263ba61063dSBarry Smith    and calling FuncFinish */
126451d315f7SKerry Stevens   while(PetscThreadGo) {
1265ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1266ba61063dSBarry Smith        waiting when you don't have to causes problems
1267ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
126851d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1269ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1270ba61063dSBarry Smith        upon return, has the lock */
127151d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
127251d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
127351d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
127451d315f7SKerry Stevens     }
127551d315f7SKerry Stevens     if(ThreadId==0) {
127651d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
127751d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
127851d315f7SKerry Stevens     }
127951d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
128051d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1281ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
128251d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
128351d315f7SKerry Stevens     }
1284ba61063dSBarry Smith     /* do your job */
128551d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
128651d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
128751d315f7SKerry Stevens     }
128851d315f7SKerry Stevens     else {
128951d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
129051d315f7SKerry Stevens     }
129151d315f7SKerry Stevens     if(iterr!=0) {
129251d315f7SKerry Stevens       ithreaderr = 1;
129351d315f7SKerry Stevens     }
129451d315f7SKerry Stevens     if(PetscThreadGo) {
1295ba61063dSBarry Smith       /* reset job, get ready for more */
129651d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1297ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1298ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
129951d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
130051d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1301ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1302ba61063dSBarry Smith            upon return, has the lock */
130351d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
130451d315f7SKerry Stevens         }
130551d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1306ba61063dSBarry Smith         /* your subordinate is now ready */
130751d315f7SKerry Stevens       }
130851d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
130951d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
131051d315f7SKerry Stevens       if(ThreadId==0) {
1311ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1312ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
131351d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
131451d315f7SKerry Stevens       }
131551d315f7SKerry Stevens       else {
1316ba61063dSBarry Smith         /* signal your boss before you go to sleep */
131751d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
131851d315f7SKerry Stevens       }
131951d315f7SKerry Stevens     }
132051d315f7SKerry Stevens   }
132151d315f7SKerry Stevens   return NULL;
132251d315f7SKerry Stevens }
132351d315f7SKerry Stevens 
132451d315f7SKerry Stevens #undef __FUNCT__
132551d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
132651d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
132751d315f7SKerry Stevens   PetscInt i,ierr;
132851d315f7SKerry Stevens   int status;
132951d315f7SKerry Stevens 
133051d315f7SKerry Stevens   if(PetscUseThreadPool) {
133151d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
133251d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
133351d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
133451d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
133551d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
133651d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
133751d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
133851d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
133951d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134051d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134151d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
134251d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1343ba61063dSBarry Smith     /* initialize job structure */
134451d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
134551d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
134651d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
134751d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
134851d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
134951d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
135051d315f7SKerry Stevens     }
135151d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
135251d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
135351d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
135451d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
135551d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
135651d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
135751d315f7SKerry Stevens     }
135851d315f7SKerry Stevens     job_chain.pfunc = NULL;
135951d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
136051d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
136151d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
136251d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1363ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
136451d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1365ba61063dSBarry Smith     /* create threads */
136651d315f7SKerry Stevens     for(i=0; i<N; i++) {
136751d315f7SKerry Stevens       pVal[i] = i;
136851d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1369ba61063dSBarry Smith       /* should check error */
137051d315f7SKerry Stevens     }
137151d315f7SKerry Stevens   }
137251d315f7SKerry Stevens   else {
137351d315f7SKerry Stevens   }
137451d315f7SKerry Stevens   return NULL;
137551d315f7SKerry Stevens }
137651d315f7SKerry Stevens 
137751d315f7SKerry Stevens 
137851d315f7SKerry Stevens #undef __FUNCT__
137951d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
138051d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
138151d315f7SKerry Stevens   int i,ierr;
138251d315f7SKerry Stevens   void* jstatus;
138351d315f7SKerry Stevens 
138451d315f7SKerry Stevens   PetscFunctionBegin;
138551d315f7SKerry Stevens 
138651d315f7SKerry Stevens   if(PetscUseThreadPool) {
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);
140151d315f7SKerry Stevens   }
140251d315f7SKerry Stevens   else {
140351d315f7SKerry Stevens   }
140451d315f7SKerry Stevens   PetscFunctionReturn(0);
140551d315f7SKerry Stevens }
140651d315f7SKerry Stevens 
140751d315f7SKerry Stevens #undef __FUNCT__
140851d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
140951d315f7SKerry Stevens void MainWait_Chain() {
141051d315f7SKerry Stevens   int ierr;
141151d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
141251d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
141351d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
141451d315f7SKerry Stevens   }
141551d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
141651d315f7SKerry Stevens }
141751d315f7SKerry Stevens 
141851d315f7SKerry Stevens #undef __FUNCT__
141951d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
142051d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
142151d315f7SKerry Stevens   int i,ierr;
142251d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
1423*36d20dc5SKerry Stevens 
142451d315f7SKerry Stevens   MainWait();
142551d315f7SKerry Stevens   job_chain.pfunc = pFunc;
142651d315f7SKerry Stevens   job_chain.pdata = data;
142751d315f7SKerry Stevens   job_chain.startJob = PETSC_TRUE;
142851d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
142951d315f7SKerry Stevens     *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
143051d315f7SKerry Stevens   }
143151d315f7SKerry Stevens   job_chain.eJobStat = JobInitiated;
143251d315f7SKerry Stevens   ierr = pthread_cond_signal(job_chain.cond2array[0]);
143351d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1434ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
143551d315f7SKerry Stevens   }
1436*36d20dc5SKerry Stevens 
143751d315f7SKerry Stevens   if(ithreaderr) {
143851d315f7SKerry Stevens     ijoberr = ithreaderr;
143951d315f7SKerry Stevens   }
144051d315f7SKerry Stevens   return ijoberr;
144151d315f7SKerry Stevens }
144251d315f7SKerry Stevens /****  ****/
144351d315f7SKerry Stevens 
1444ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
144551d315f7SKerry Stevens /**** True Thread Functions ****/
144651d315f7SKerry Stevens void* PetscThreadFunc_True(void* arg) {
144751d315f7SKerry Stevens   int icorr,ierr,iVal;
144851dcc849SKerry Stevens   int* pId = (int*)arg;
144951dcc849SKerry Stevens   int ThreadId = *pId;
14500ca81413SKerry Stevens   PetscErrorCode iterr;
145151d315f7SKerry Stevens   cpu_set_t mset;
14529e800a48SKerry Stevens   //printf("Thread %d In True Pool Thread Function\n",ThreadId);
145351d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
145451d315f7SKerry Stevens   CPU_ZERO(&mset);
145551d315f7SKerry Stevens   CPU_SET(icorr,&mset);
145651d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
145751d315f7SKerry Stevens 
145851d315f7SKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
145951d315f7SKerry Stevens   job_true.iNumReadyThreads++;
146051d315f7SKerry Stevens   if(job_true.iNumReadyThreads==PetscMaxThreads) {
146151dcc849SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
146251dcc849SKerry Stevens   }
1463ba61063dSBarry Smith   /*the while loop needs to have an exit
1464ba61063dSBarry Smith     the 'main' thread can terminate all the threads by performing a broadcast
1465ba61063dSBarry Smith    and calling FuncFinish */
146651dcc849SKerry Stevens   while(PetscThreadGo) {
1467ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
1468ba61063dSBarry Smith       waiting when you don't have to causes problems
1469ba61063dSBarry Smith      also need to wait if another thread sneaks in and messes with the predicate */
147051d315f7SKerry Stevens     while(job_true.startJob==PETSC_FALSE&&job_true.iNumJobThreads==0) {
1471ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1472ba61063dSBarry Smith        upon return, has the lock */
147351d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
147451dcc849SKerry Stevens     }
147551d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
147651d315f7SKerry Stevens     job_true.iNumJobThreads--;
147751d315f7SKerry Stevens     job_true.iNumReadyThreads--;
147851d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
147951d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
148051d315f7SKerry Stevens     if(job_true.pdata==NULL) {
148151d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
148251dcc849SKerry Stevens     }
148351dcc849SKerry Stevens     else {
148451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
148551dcc849SKerry Stevens     }
14860ca81413SKerry Stevens     if(iterr!=0) {
14870ca81413SKerry Stevens       ithreaderr = 1;
14880ca81413SKerry Stevens     }
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 
15120ca81413SKerry Stevens   if(PetscUseThreadPool) {
1513f09cb4aaSKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1514ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
151551dcc849SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1516ba61063dSBarry Smith     BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
151751d315f7SKerry Stevens     job_true.pdata = (void**)malloc(N*sizeof(void*));
151851dcc849SKerry Stevens     for(i=0; i<N; i++) {
1519f09cb4aaSKerry Stevens       pVal[i] = i;
1520f09cb4aaSKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1521ba61063dSBarry Smith       /* error check to ensure proper thread creation */
152251dcc849SKerry Stevens       status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1523ba61063dSBarry Smith       /* should check error */
152451dcc849SKerry Stevens     }
15250ca81413SKerry Stevens   }
15260ca81413SKerry Stevens   else {
15270ca81413SKerry Stevens   }
152851dcc849SKerry Stevens   return NULL;
152951dcc849SKerry Stevens }
153051dcc849SKerry Stevens 
1531f09cb4aaSKerry Stevens 
1532f09cb4aaSKerry Stevens #undef __FUNCT__
153351d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
153451d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
153551dcc849SKerry Stevens   int i,ierr;
153651dcc849SKerry Stevens   void* jstatus;
153751dcc849SKerry Stevens 
153851dcc849SKerry Stevens   PetscFunctionBegin;
15390ca81413SKerry Stevens 
15400ca81413SKerry Stevens   if(PetscUseThreadPool) {
1541ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1542ba61063dSBarry Smith     /* join the threads */
154351dcc849SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
154451dcc849SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1545ba61063dSBarry Smith       /* should check error */
154651dcc849SKerry Stevens     }
154751dcc849SKerry Stevens     free(BarrPoint);
154851dcc849SKerry Stevens     free(PetscThreadPoint);
15490ca81413SKerry Stevens   }
15500ca81413SKerry Stevens   else {
15510ca81413SKerry Stevens   }
155251dcc849SKerry Stevens   PetscFunctionReturn(0);
155351dcc849SKerry Stevens }
155451dcc849SKerry Stevens 
1555f09cb4aaSKerry Stevens #undef __FUNCT__
155651d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
155751d315f7SKerry Stevens void MainWait_True() {
155851dcc849SKerry Stevens   int ierr;
155951d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
156051d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
156151dcc849SKerry Stevens   }
156251d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
156351dcc849SKerry Stevens }
156451dcc849SKerry Stevens 
1565f09cb4aaSKerry Stevens #undef __FUNCT__
156651d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
156751d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
156851dcc849SKerry Stevens   int ierr;
15690ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
1570*36d20dc5SKerry Stevens 
15710ca81413SKerry Stevens   MainWait();
157251d315f7SKerry Stevens   job_true.pfunc = pFunc;
157351d315f7SKerry Stevens   job_true.pdata = data;
157451d315f7SKerry Stevens   job_true.pbarr = &BarrPoint[n];
157551d315f7SKerry Stevens   job_true.iNumJobThreads = n;
157651d315f7SKerry Stevens   job_true.startJob = PETSC_TRUE;
157751d315f7SKerry Stevens   ierr = pthread_cond_broadcast(&job_true.cond);
15780ca81413SKerry Stevens   if(pFunc!=FuncFinish) {
1579ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done */
15800ca81413SKerry Stevens   }
1581*36d20dc5SKerry Stevens 
15820ca81413SKerry Stevens   if(ithreaderr) {
15830ca81413SKerry Stevens     ijoberr = ithreaderr;
15840ca81413SKerry Stevens   }
15850ca81413SKerry Stevens   return ijoberr;
158651dcc849SKerry Stevens }
1587683509dcSKerry Stevens /**** NO THREAD POOL FUNCTION ****/
1588683509dcSKerry Stevens #undef __FUNCT__
1589683509dcSKerry Stevens #define __FUNCT__ "MainJob_Spawn"
1590683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void** data,PetscInt n) {
1591683509dcSKerry Stevens   PetscErrorCode ijoberr = 0;
1592683509dcSKerry Stevens 
1593683509dcSKerry Stevens   pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
1594683509dcSKerry Stevens   PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1595683509dcSKerry Stevens   PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
1596683509dcSKerry Stevens   free(apThread);
1597683509dcSKerry Stevens 
1598683509dcSKerry Stevens   return ijoberr;
1599683509dcSKerry Stevens }
160051d315f7SKerry Stevens /****  ****/
1601ba61063dSBarry Smith #endif
160251dcc849SKerry Stevens 
160351dcc849SKerry Stevens void* FuncFinish(void* arg) {
160451dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16050ca81413SKerry Stevens   return(0);
160651dcc849SKerry Stevens }
1607ba61063dSBarry Smith 
1608ba61063dSBarry Smith #endif
1609