xref: /petsc/src/sys/objects/init.c (revision 4d4e02f7946d71247125f2631c72ecc068aab84b)
1e5c89e4eSSatish Balay /*
2e5c89e4eSSatish Balay 
3e5c89e4eSSatish Balay    This file defines part of the initialization of PETSc
4e5c89e4eSSatish Balay 
5e5c89e4eSSatish Balay   This file uses regular malloc and free because it cannot know
6e5c89e4eSSatish Balay   what malloc is being used until it has already processed the input.
7e5c89e4eSSatish Balay */
8*4d4e02f7SKerry Stevens //#if defined(PETSC_HAVE_SCHED_H) && defined(PETSC_USE_PTHREAD)
98f54b378SBarry Smith #ifndef _GNU_SOURCE
108f54b378SBarry Smith #define _GNU_SOURCE
118f54b378SBarry Smith #endif
128f54b378SBarry Smith #include <sched.h>
13*4d4e02f7SKerry Stevens //#endif
14*4d4e02f7SKerry Stevens #include <petscsys.h>        /*I  "petscsys.h"   I*/
15ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD)
1651dcc849SKerry Stevens #include <pthread.h>
17ba61063dSBarry Smith #endif
18ba61063dSBarry Smith #if defined(PETSC_HAVE_SYS_SYSINFO_H)
1951d315f7SKerry Stevens #include <sys/sysinfo.h>
20ba61063dSBarry Smith #endif
2151d315f7SKerry Stevens #include <unistd.h>
22e5c89e4eSSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
23e5c89e4eSSatish Balay #include <stdlib.h>
24e5c89e4eSSatish Balay #endif
25e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
26e5c89e4eSSatish Balay #include <malloc.h>
27e5c89e4eSSatish Balay #endif
28555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
29555d055bSBarry Smith #include <valgrind/valgrind.h>
30555d055bSBarry Smith #endif
31555d055bSBarry Smith 
32e5c89e4eSSatish Balay /* ------------------------Nasty global variables -------------------------------*/
33e5c89e4eSSatish Balay /*
34e5c89e4eSSatish Balay      Indicates if PETSc started up MPI, or it was
35e5c89e4eSSatish Balay    already started before PETSc was initialized.
36e5c89e4eSSatish Balay */
377087cfbeSBarry Smith PetscBool    PetscBeganMPI         = PETSC_FALSE;
387087cfbeSBarry Smith PetscBool    PetscInitializeCalled = PETSC_FALSE;
397087cfbeSBarry Smith PetscBool    PetscFinalizeCalled   = PETSC_FALSE;
4051dcc849SKerry Stevens PetscBool    PetscUseThreadPool    = PETSC_FALSE;
4151dcc849SKerry Stevens PetscBool    PetscThreadGo         = PETSC_TRUE;
427087cfbeSBarry Smith PetscMPIInt  PetscGlobalRank = -1;
437087cfbeSBarry Smith PetscMPIInt  PetscGlobalSize = -1;
44ba61063dSBarry Smith 
45ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
4651dcc849SKerry Stevens PetscMPIInt  PetscMaxThreads = 2;
4751dcc849SKerry Stevens pthread_t*   PetscThreadPoint;
48af359df3SBarry Smith #define PETSC_HAVE_PTHREAD_BARRIER
49ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
50ba61063dSBarry Smith pthread_barrier_t* BarrPoint;   /* used by 'true' thread pool */
51ba61063dSBarry Smith #endif
5251d315f7SKerry Stevens PetscErrorCode ithreaderr = 0;
53f09cb4aaSKerry Stevens int*         pVal;
5451dcc849SKerry Stevens 
55ba61063dSBarry Smith #define CACHE_LINE_SIZE 64  /* used by 'chain', 'main','tree' thread pools */
5651d315f7SKerry Stevens int* ThreadCoreAffinity;
5751d315f7SKerry Stevens 
58ba61063dSBarry Smith typedef enum {JobInitiated,ThreadsWorking,JobCompleted} estat;  /* used by 'chain','tree' thread pool */
5951d315f7SKerry Stevens 
6051d315f7SKerry Stevens typedef struct {
6151d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
6251d315f7SKerry Stevens   pthread_cond_t**  cond1array;
6351d315f7SKerry Stevens   pthread_cond_t** cond2array;
6451d315f7SKerry Stevens   void* (*pfunc)(void*);
6551d315f7SKerry Stevens   void** pdata;
6651d315f7SKerry Stevens   PetscBool startJob;
6751d315f7SKerry Stevens   estat eJobStat;
6851d315f7SKerry Stevens   PetscBool** arrThreadStarted;
6951d315f7SKerry Stevens   PetscBool** arrThreadReady;
7051d315f7SKerry Stevens } sjob_tree;
7151d315f7SKerry Stevens sjob_tree job_tree;
7251d315f7SKerry Stevens typedef struct {
7351d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
7451d315f7SKerry Stevens   pthread_cond_t**  cond1array;
7551d315f7SKerry Stevens   pthread_cond_t** cond2array;
7651d315f7SKerry Stevens   void* (*pfunc)(void*);
7751d315f7SKerry Stevens   void** pdata;
7851d315f7SKerry Stevens   PetscBool** arrThreadReady;
7951d315f7SKerry Stevens } sjob_main;
8051d315f7SKerry Stevens sjob_main job_main;
8151d315f7SKerry Stevens typedef struct {
8251d315f7SKerry Stevens   pthread_mutex_t** mutexarray;
8351d315f7SKerry Stevens   pthread_cond_t**  cond1array;
8451d315f7SKerry Stevens   pthread_cond_t** cond2array;
8551d315f7SKerry Stevens   void* (*pfunc)(void*);
8651d315f7SKerry Stevens   void** pdata;
8751d315f7SKerry Stevens   PetscBool startJob;
8851d315f7SKerry Stevens   estat eJobStat;
8951d315f7SKerry Stevens   PetscBool** arrThreadStarted;
9051d315f7SKerry Stevens   PetscBool** arrThreadReady;
9151d315f7SKerry Stevens } sjob_chain;
9251d315f7SKerry Stevens sjob_chain job_chain;
93ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
9451dcc849SKerry Stevens typedef struct {
9551dcc849SKerry Stevens   pthread_mutex_t mutex;
9651dcc849SKerry Stevens   pthread_cond_t cond;
9751dcc849SKerry Stevens   void* (*pfunc)(void*);
9851dcc849SKerry Stevens   void** pdata;
9951dcc849SKerry Stevens   pthread_barrier_t* pbarr;
10051dcc849SKerry Stevens   int iNumJobThreads;
10151dcc849SKerry Stevens   int iNumReadyThreads;
10251dcc849SKerry Stevens   PetscBool startJob;
10351d315f7SKerry Stevens } sjob_true;
10451d315f7SKerry Stevens sjob_true job_true = {PTHREAD_MUTEX_INITIALIZER,PTHREAD_COND_INITIALIZER,NULL,NULL,NULL,0,0,PETSC_FALSE};
105ba61063dSBarry Smith #endif
10651dcc849SKerry Stevens 
107ba61063dSBarry Smith pthread_cond_t  main_cond  = PTHREAD_COND_INITIALIZER;  /* used by 'true', 'chain','tree' thread pools */
108ba61063dSBarry Smith char* arrmutex; /* used by 'chain','main','tree' thread pools */
109ba61063dSBarry Smith char* arrcond1; /* used by 'chain','main','tree' thread pools */
110ba61063dSBarry Smith char* arrcond2; /* used by 'chain','main','tree' thread pools */
111ba61063dSBarry Smith char* arrstart; /* used by 'chain','main','tree' thread pools */
112ba61063dSBarry Smith char* arrready; /* used by 'chain','main','tree' thread pools */
11351dcc849SKerry Stevens 
11451d315f7SKerry Stevens /* Function Pointers */
11551d315f7SKerry Stevens void*          (*PetscThreadFunc)(void*) = NULL;
11651d315f7SKerry Stevens void*          (*PetscThreadInitialize)(PetscInt) = NULL;
11751d315f7SKerry Stevens PetscErrorCode (*PetscThreadFinalize)(void) = NULL;
11851d315f7SKerry Stevens void           (*MainWait)(void) = NULL;
11951d315f7SKerry Stevens PetscErrorCode (*MainJob)(void* (*pFunc)(void*),void**,PetscInt) = NULL;
12036d20dc5SKerry Stevens /**** Tree Thread Pool Functions ****/
12151d315f7SKerry Stevens void*          PetscThreadFunc_Tree(void*);
12251d315f7SKerry Stevens void*          PetscThreadInitialize_Tree(PetscInt);
12351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree(void);
12451d315f7SKerry Stevens void           MainWait_Tree(void);
12551d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void**,PetscInt);
12636d20dc5SKerry Stevens /**** Main Thread Pool Functions ****/
12751d315f7SKerry Stevens void*          PetscThreadFunc_Main(void*);
12851d315f7SKerry Stevens void*          PetscThreadInitialize_Main(PetscInt);
12951d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main(void);
13051d315f7SKerry Stevens void           MainWait_Main(void);
13151d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void**,PetscInt);
13236d20dc5SKerry Stevens /**** Chain Thread Pool Functions ****/
13351d315f7SKerry Stevens void*          PetscThreadFunc_Chain(void*);
13451d315f7SKerry Stevens void*          PetscThreadInitialize_Chain(PetscInt);
13551d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain(void);
13651d315f7SKerry Stevens void           MainWait_Chain(void);
13751d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void**,PetscInt);
13836d20dc5SKerry Stevens /**** True Thread Pool Functions ****/
13951d315f7SKerry Stevens void*          PetscThreadFunc_True(void*);
14051d315f7SKerry Stevens void*          PetscThreadInitialize_True(PetscInt);
14151d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True(void);
14251d315f7SKerry Stevens void           MainWait_True(void);
14351d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void**,PetscInt);
14436d20dc5SKerry Stevens /**** NO Thread Pool Function  ****/
145683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void**,PetscInt);
14636d20dc5SKerry Stevens /****  ****/
14751dcc849SKerry Stevens void* FuncFinish(void*);
1480ca81413SKerry Stevens void* PetscThreadRun(MPI_Comm Comm,void* (*pFunc)(void*),int,pthread_t*,void**);
1490ca81413SKerry Stevens void* PetscThreadStop(MPI_Comm Comm,int,pthread_t*);
150ba61063dSBarry Smith #endif
151e5c89e4eSSatish Balay 
152e5c89e4eSSatish Balay #if defined(PETSC_USE_COMPLEX)
153e5c89e4eSSatish Balay #if defined(PETSC_COMPLEX_INSTANTIATE)
154e5c89e4eSSatish Balay template <> class std::complex<double>; /* instantiate complex template class */
155e5c89e4eSSatish Balay #endif
1562c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1577087cfbeSBarry Smith MPI_Datatype   MPI_C_DOUBLE_COMPLEX;
1587087cfbeSBarry Smith MPI_Datatype   MPI_C_COMPLEX;
1592c876bd9SBarry Smith #endif
1607087cfbeSBarry Smith PetscScalar    PETSC_i;
161e5c89e4eSSatish Balay #else
1627087cfbeSBarry Smith PetscScalar    PETSC_i = 0.0;
163e5c89e4eSSatish Balay #endif
164ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
165c90a1750SBarry Smith MPI_Datatype   MPIU___FLOAT128 = 0;
166c90a1750SBarry Smith #endif
1677087cfbeSBarry Smith MPI_Datatype   MPIU_2SCALAR = 0;
1687087cfbeSBarry Smith MPI_Datatype   MPIU_2INT = 0;
16975567043SBarry Smith 
170e5c89e4eSSatish Balay /*
171e5c89e4eSSatish Balay      These are needed by petscbt.h
172e5c89e4eSSatish Balay */
173c6db04a5SJed Brown #include <petscbt.h>
1747087cfbeSBarry Smith char      _BT_mask = ' ';
1757087cfbeSBarry Smith char      _BT_c = ' ';
1767087cfbeSBarry Smith PetscInt  _BT_idx  = 0;
177e5c89e4eSSatish Balay 
178e5c89e4eSSatish Balay /*
179e5c89e4eSSatish Balay        Function that is called to display all error messages
180e5c89e4eSSatish Balay */
1817087cfbeSBarry Smith PetscErrorCode  (*PetscErrorPrintf)(const char [],...)          = PetscErrorPrintfDefault;
1827087cfbeSBarry Smith PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char [],...)  = PetscHelpPrintfDefault;
183238ccf28SShri Abhyankar #if defined(PETSC_HAVE_MATLAB_ENGINE)
1847087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintf_Matlab;
185238ccf28SShri Abhyankar #else
1867087cfbeSBarry Smith PetscErrorCode  (*PetscVFPrintf)(FILE*,const char[],va_list)    = PetscVFPrintfDefault;
187238ccf28SShri Abhyankar #endif
188bab1f7e6SVictor Minden /*
1898154be41SBarry Smith   This is needed to turn on/off cusp synchronization */
1908154be41SBarry Smith PetscBool   synchronizeCUSP = PETSC_FALSE;
191bab1f7e6SVictor Minden 
192e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
193e5c89e4eSSatish Balay /*
194e5c89e4eSSatish Balay    Optional file where all PETSc output from various prints is saved
195e5c89e4eSSatish Balay */
196e5c89e4eSSatish Balay FILE *petsc_history = PETSC_NULL;
197e5c89e4eSSatish Balay 
198e5c89e4eSSatish Balay #undef __FUNCT__
199f3dea69dSBarry Smith #define __FUNCT__ "PetscOpenHistoryFile"
2007087cfbeSBarry Smith PetscErrorCode  PetscOpenHistoryFile(const char filename[],FILE **fd)
201e5c89e4eSSatish Balay {
202e5c89e4eSSatish Balay   PetscErrorCode ierr;
203e5c89e4eSSatish Balay   PetscMPIInt    rank,size;
204e5c89e4eSSatish Balay   char           pfile[PETSC_MAX_PATH_LEN],pname[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN],date[64];
205e5c89e4eSSatish Balay   char           version[256];
206e5c89e4eSSatish Balay 
207e5c89e4eSSatish Balay   PetscFunctionBegin;
208e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
209e5c89e4eSSatish Balay   if (!rank) {
210e5c89e4eSSatish Balay     char        arch[10];
211f56c2debSBarry Smith     int         err;
21288c29154SBarry Smith     PetscViewer viewer;
213f56c2debSBarry Smith 
214e5c89e4eSSatish Balay     ierr = PetscGetArchType(arch,10);CHKERRQ(ierr);
215e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
216a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
217e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
218e5c89e4eSSatish Balay     if (filename) {
219e5c89e4eSSatish Balay       ierr = PetscFixFilename(filename,fname);CHKERRQ(ierr);
220e5c89e4eSSatish Balay     } else {
221e5c89e4eSSatish Balay       ierr = PetscGetHomeDirectory(pfile,240);CHKERRQ(ierr);
222e5c89e4eSSatish Balay       ierr = PetscStrcat(pfile,"/.petschistory");CHKERRQ(ierr);
223e5c89e4eSSatish Balay       ierr = PetscFixFilename(pfile,fname);CHKERRQ(ierr);
224e5c89e4eSSatish Balay     }
225e5c89e4eSSatish Balay 
226e32f2f54SBarry Smith     *fd = fopen(fname,"a"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file: %s",fname);
227e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
228e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s %s\n",version,date);CHKERRQ(ierr);
229e5c89e4eSSatish Balay     ierr = PetscGetProgramName(pname,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
230e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"%s on a %s, %d proc. with options:\n",pname,arch,size);CHKERRQ(ierr);
23188c29154SBarry Smith     ierr = PetscViewerASCIIOpenWithFILE(PETSC_COMM_WORLD,*fd,&viewer);CHKERRQ(ierr);
23288c29154SBarry Smith     ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
2336bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
234e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
235f56c2debSBarry Smith     err = fflush(*fd);
236e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
237e5c89e4eSSatish Balay   }
238e5c89e4eSSatish Balay   PetscFunctionReturn(0);
239e5c89e4eSSatish Balay }
240e5c89e4eSSatish Balay 
241e5c89e4eSSatish Balay #undef __FUNCT__
242f3dea69dSBarry Smith #define __FUNCT__ "PetscCloseHistoryFile"
2437087cfbeSBarry Smith PetscErrorCode  PetscCloseHistoryFile(FILE **fd)
244e5c89e4eSSatish Balay {
245e5c89e4eSSatish Balay   PetscErrorCode ierr;
246e5c89e4eSSatish Balay   PetscMPIInt    rank;
247e5c89e4eSSatish Balay   char           date[64];
248f56c2debSBarry Smith   int            err;
249e5c89e4eSSatish Balay 
250e5c89e4eSSatish Balay   PetscFunctionBegin;
251e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
252e5c89e4eSSatish Balay   if (!rank) {
253e5c89e4eSSatish Balay     ierr = PetscGetDate(date,64);CHKERRQ(ierr);
254e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
255e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"Finished at %s\n",date);CHKERRQ(ierr);
256e5c89e4eSSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,*fd,"---------------------------------------------------------\n");CHKERRQ(ierr);
257f56c2debSBarry Smith     err = fflush(*fd);
258e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
259f56c2debSBarry Smith     err = fclose(*fd);
260e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
261e5c89e4eSSatish Balay   }
262e5c89e4eSSatish Balay   PetscFunctionReturn(0);
263e5c89e4eSSatish Balay }
264e5c89e4eSSatish Balay 
265e5c89e4eSSatish Balay /* ------------------------------------------------------------------------------*/
266e5c89e4eSSatish Balay 
267e5c89e4eSSatish Balay /*
268e5c89e4eSSatish Balay    This is ugly and probably belongs somewhere else, but I want to
269e5c89e4eSSatish Balay   be able to put a true MPI abort error handler with command line args.
270e5c89e4eSSatish Balay 
271e5c89e4eSSatish Balay     This is so MPI errors in the debugger will leave all the stack
2723c311c98SBarry Smith   frames. The default MP_Abort() cleans up and exits thus providing no useful information
2733c311c98SBarry Smith   in the debugger hence we call abort() instead of MPI_Abort().
274e5c89e4eSSatish Balay */
275e5c89e4eSSatish Balay 
276e5c89e4eSSatish Balay #undef __FUNCT__
277e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_AbortOnError"
278e5c89e4eSSatish Balay void Petsc_MPI_AbortOnError(MPI_Comm *comm,PetscMPIInt *flag)
279e5c89e4eSSatish Balay {
280e5c89e4eSSatish Balay   PetscFunctionBegin;
2813c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
282e5c89e4eSSatish Balay   abort();
283e5c89e4eSSatish Balay }
284e5c89e4eSSatish Balay 
285e5c89e4eSSatish Balay #undef __FUNCT__
286e5c89e4eSSatish Balay #define __FUNCT__ "Petsc_MPI_DebuggerOnError"
287e5c89e4eSSatish Balay void Petsc_MPI_DebuggerOnError(MPI_Comm *comm,PetscMPIInt *flag)
288e5c89e4eSSatish Balay {
289e5c89e4eSSatish Balay   PetscErrorCode ierr;
290e5c89e4eSSatish Balay 
291e5c89e4eSSatish Balay   PetscFunctionBegin;
2923c311c98SBarry Smith   (*PetscErrorPrintf)("MPI error %d\n",*flag);
293e5c89e4eSSatish Balay   ierr = PetscAttachDebugger();
294e5c89e4eSSatish Balay   if (ierr) { /* hopeless so get out */
2953c311c98SBarry Smith     MPI_Abort(*comm,*flag);
296e5c89e4eSSatish Balay   }
297e5c89e4eSSatish Balay }
298e5c89e4eSSatish Balay 
299e5c89e4eSSatish Balay #undef __FUNCT__
300e5c89e4eSSatish Balay #define __FUNCT__ "PetscEnd"
301e5c89e4eSSatish Balay /*@C
302e5c89e4eSSatish Balay    PetscEnd - Calls PetscFinalize() and then ends the program. This is useful if one
303e5c89e4eSSatish Balay      wishes a clean exit somewhere deep in the program.
304e5c89e4eSSatish Balay 
305e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
306e5c89e4eSSatish Balay 
307e5c89e4eSSatish Balay    Options Database Keys are the same as for PetscFinalize()
308e5c89e4eSSatish Balay 
309e5c89e4eSSatish Balay    Level: advanced
310e5c89e4eSSatish Balay 
311e5c89e4eSSatish Balay    Note:
312e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
313e5c89e4eSSatish Balay 
31488c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscFinalize()
315e5c89e4eSSatish Balay @*/
3167087cfbeSBarry Smith PetscErrorCode  PetscEnd(void)
317e5c89e4eSSatish Balay {
318e5c89e4eSSatish Balay   PetscFunctionBegin;
319e5c89e4eSSatish Balay   PetscFinalize();
320e5c89e4eSSatish Balay   exit(0);
321e5c89e4eSSatish Balay   return 0;
322e5c89e4eSSatish Balay }
323e5c89e4eSSatish Balay 
324ace3abfcSBarry Smith PetscBool    PetscOptionsPublish = PETSC_FALSE;
32509573ac7SBarry Smith extern PetscErrorCode        PetscSetUseTrMalloc_Private(void);
326ace3abfcSBarry Smith extern PetscBool  petscsetmallocvisited;
327e5c89e4eSSatish Balay static char       emacsmachinename[256];
328e5c89e4eSSatish Balay 
329e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalVersionFunction)(MPI_Comm) = 0;
330e5c89e4eSSatish Balay PetscErrorCode (*PetscExternalHelpFunction)(MPI_Comm)    = 0;
331e5c89e4eSSatish Balay 
332e5c89e4eSSatish Balay #undef __FUNCT__
333e5c89e4eSSatish Balay #define __FUNCT__ "PetscSetHelpVersionFunctions"
334e5c89e4eSSatish Balay /*@C
335e5c89e4eSSatish Balay    PetscSetHelpVersionFunctions - Sets functions that print help and version information
336e5c89e4eSSatish Balay    before the PETSc help and version information is printed. Must call BEFORE PetscInitialize().
337e5c89e4eSSatish Balay    This routine enables a "higher-level" package that uses PETSc to print its messages first.
338e5c89e4eSSatish Balay 
339e5c89e4eSSatish Balay    Input Parameter:
340e5c89e4eSSatish Balay +  help - the help function (may be PETSC_NULL)
341da93591fSBarry Smith -  version - the version function (may be PETSC_NULL)
342e5c89e4eSSatish Balay 
343e5c89e4eSSatish Balay    Level: developer
344e5c89e4eSSatish Balay 
345e5c89e4eSSatish Balay    Concepts: package help message
346e5c89e4eSSatish Balay 
347e5c89e4eSSatish Balay @*/
3487087cfbeSBarry Smith PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*help)(MPI_Comm),PetscErrorCode (*version)(MPI_Comm))
349e5c89e4eSSatish Balay {
350e5c89e4eSSatish Balay   PetscFunctionBegin;
351e5c89e4eSSatish Balay   PetscExternalHelpFunction    = help;
352e5c89e4eSSatish Balay   PetscExternalVersionFunction = version;
353e5c89e4eSSatish Balay   PetscFunctionReturn(0);
354e5c89e4eSSatish Balay }
355e5c89e4eSSatish Balay 
356e5c89e4eSSatish Balay #undef __FUNCT__
357e5c89e4eSSatish Balay #define __FUNCT__ "PetscOptionsCheckInitial_Private"
3587087cfbeSBarry Smith PetscErrorCode  PetscOptionsCheckInitial_Private(void)
359e5c89e4eSSatish Balay {
360e5c89e4eSSatish Balay   char           string[64],mname[PETSC_MAX_PATH_LEN],*f;
361e5c89e4eSSatish Balay   MPI_Comm       comm = PETSC_COMM_WORLD;
362ace3abfcSBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flag,flgz,flgzout;
363e5c89e4eSSatish Balay   PetscErrorCode ierr;
364a6d0e24fSJed Brown   PetscReal      si;
365e5c89e4eSSatish Balay   int            i;
366e5c89e4eSSatish Balay   PetscMPIInt    rank;
367e5c89e4eSSatish Balay   char           version[256];
368e5c89e4eSSatish Balay 
369e5c89e4eSSatish Balay   PetscFunctionBegin;
370e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
371e5c89e4eSSatish Balay 
372e5c89e4eSSatish Balay   /*
373e5c89e4eSSatish Balay       Setup the memory management; support for tracing malloc() usage
374e5c89e4eSSatish Balay   */
3758bb29257SSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-malloc_log",&flg3);CHKERRQ(ierr);
37681b192fdSBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
377acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg1,&flg2);CHKERRQ(ierr);
378e5c89e4eSSatish Balay   if ((!flg2 || flg1) && !petscsetmallocvisited) {
379555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
380555d055bSBarry Smith     if (flg2 || !(RUNNING_ON_VALGRIND)) {
381555d055bSBarry Smith       /* turn off default -malloc if valgrind is being used */
382555d055bSBarry Smith #endif
383e5c89e4eSSatish Balay       ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
384555d055bSBarry Smith #if defined(PETSC_HAVE_VALGRIND)
385555d055bSBarry Smith     }
386555d055bSBarry Smith #endif
387e5c89e4eSSatish Balay   }
388e5c89e4eSSatish Balay #else
389acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_dump",&flg1,PETSC_NULL);CHKERRQ(ierr);
390acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc",&flg2,PETSC_NULL);CHKERRQ(ierr);
391e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3) {ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);}
392e5c89e4eSSatish Balay #endif
393e5c89e4eSSatish Balay   if (flg3) {
394e5c89e4eSSatish Balay     ierr = PetscMallocSetDumpLog();CHKERRQ(ierr);
395e5c89e4eSSatish Balay   }
39690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
397acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_debug",&flg1,PETSC_NULL);CHKERRQ(ierr);
398e5c89e4eSSatish Balay   if (flg1) {
399e5c89e4eSSatish Balay     ierr = PetscSetUseTrMalloc_Private();CHKERRQ(ierr);
400e5c89e4eSSatish Balay     ierr = PetscMallocDebug(PETSC_TRUE);CHKERRQ(ierr);
401e5c89e4eSSatish Balay   }
402e5c89e4eSSatish Balay 
40390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
404acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-malloc_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4057783f70dSSatish Balay   if (!flg1) {
40690d69ab7SBarry Smith     flg1 = PETSC_FALSE;
407acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-memory_info",&flg1,PETSC_NULL);CHKERRQ(ierr);
4087783f70dSSatish Balay   }
409e5c89e4eSSatish Balay   if (flg1) {
410e5c89e4eSSatish Balay     ierr = PetscMemorySetGetMaximumUsage();CHKERRQ(ierr);
411e5c89e4eSSatish Balay   }
412e5c89e4eSSatish Balay 
413e5c89e4eSSatish Balay   /*
414e5c89e4eSSatish Balay       Set the display variable for graphics
415e5c89e4eSSatish Balay   */
416e5c89e4eSSatish Balay   ierr = PetscSetDisplay();CHKERRQ(ierr);
417e5c89e4eSSatish Balay 
41851dcc849SKerry Stevens 
41951dcc849SKerry Stevens   /*
420e5c89e4eSSatish Balay       Print the PETSc version information
421e5c89e4eSSatish Balay   */
422e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-v",&flg1);CHKERRQ(ierr);
423e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-version",&flg2);CHKERRQ(ierr);
424e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg3);CHKERRQ(ierr);
425e5c89e4eSSatish Balay   if (flg1 || flg2 || flg3){
426e5c89e4eSSatish Balay 
427e5c89e4eSSatish Balay     /*
428e5c89e4eSSatish Balay        Print "higher-level" package version message
429e5c89e4eSSatish Balay     */
430e5c89e4eSSatish Balay     if (PetscExternalVersionFunction) {
431e5c89e4eSSatish Balay       ierr = (*PetscExternalVersionFunction)(comm);CHKERRQ(ierr);
432e5c89e4eSSatish Balay     }
433e5c89e4eSSatish Balay 
434a523d312SBarry Smith     ierr = PetscGetVersion(version,256);CHKERRQ(ierr);
435e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
436e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
437e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s\n",version);CHKERRQ(ierr);
438e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"%s",PETSC_AUTHOR_INFO);CHKERRQ(ierr);
439e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/changes/index.html for recent updates.\n");CHKERRQ(ierr);
44084e42920SBarry Smith     ierr = (*PetscHelpPrintf)(comm,"See docs/faq.html for problems.\n");CHKERRQ(ierr);
441e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"See docs/manualpages/index.html for help. \n");CHKERRQ(ierr);
442e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Libraries linked from %s\n",PETSC_LIB_DIR);CHKERRQ(ierr);
443e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"--------------------------------------------\
444e5c89e4eSSatish Balay ------------------------------\n");CHKERRQ(ierr);
445e5c89e4eSSatish Balay   }
446e5c89e4eSSatish Balay 
447e5c89e4eSSatish Balay   /*
448e5c89e4eSSatish Balay        Print "higher-level" package help message
449e5c89e4eSSatish Balay   */
450e5c89e4eSSatish Balay   if (flg3){
451e5c89e4eSSatish Balay     if (PetscExternalHelpFunction) {
452e5c89e4eSSatish Balay       ierr = (*PetscExternalHelpFunction)(comm);CHKERRQ(ierr);
453e5c89e4eSSatish Balay     }
454e5c89e4eSSatish Balay   }
455e5c89e4eSSatish Balay 
456e5c89e4eSSatish Balay   /*
457e5c89e4eSSatish Balay       Setup the error handling
458e5c89e4eSSatish Balay   */
45990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
460acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_abort",&flg1,PETSC_NULL);CHKERRQ(ierr);
461cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);}
46290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
463acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-on_error_mpiabort",&flg1,PETSC_NULL);CHKERRQ(ierr);
464cb9801acSJed Brown   if (flg1) { ierr = PetscPushErrorHandler(PetscMPIAbortErrorHandler,0);CHKERRQ(ierr);}
46590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
466acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mpi_return_on_error",&flg1,PETSC_NULL);CHKERRQ(ierr);
467e5c89e4eSSatish Balay   if (flg1) {
468e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,MPI_ERRORS_RETURN);CHKERRQ(ierr);
469e5c89e4eSSatish Balay   }
47090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
471acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-no_signal_handler",&flg1,PETSC_NULL);CHKERRQ(ierr);
472cb9801acSJed Brown   if (!flg1) {ierr = PetscPushSignalHandler(PetscDefaultSignalHandler,(void*)0);CHKERRQ(ierr);}
47396cc47afSJed Brown   flg1 = PETSC_FALSE;
474acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-fp_trap",&flg1,PETSC_NULL);CHKERRQ(ierr);
47596cc47afSJed Brown   if (flg1) {ierr = PetscSetFPTrap(PETSC_FP_TRAP_ON);CHKERRQ(ierr);}
476e5c89e4eSSatish Balay 
477e5c89e4eSSatish Balay   /*
478e5c89e4eSSatish Balay       Setup debugger information
479e5c89e4eSSatish Balay   */
480e5c89e4eSSatish Balay   ierr = PetscSetDefaultDebugger();CHKERRQ(ierr);
481e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_attach_debugger",string,64,&flg1);CHKERRQ(ierr);
482e5c89e4eSSatish Balay   if (flg1) {
483e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
484e5c89e4eSSatish Balay 
485e5c89e4eSSatish Balay     ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
486e5c89e4eSSatish Balay     ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_DebuggerOnError,&err_handler);CHKERRQ(ierr);
487e5c89e4eSSatish Balay     ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
488e5c89e4eSSatish Balay     ierr = PetscPushErrorHandler(PetscAttachDebuggerErrorHandler,0);CHKERRQ(ierr);
489e5c89e4eSSatish Balay   }
4905e96ac45SJed Brown   ierr = PetscOptionsGetString(PETSC_NULL,"-debug_terminal",string,64,&flg1);CHKERRQ(ierr);
4915e96ac45SJed Brown   if (flg1) { ierr = PetscSetDebugTerminal(string);CHKERRQ(ierr); }
492e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-start_in_debugger",string,64,&flg1);CHKERRQ(ierr);
493e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-stop_for_debugger",string,64,&flg2);CHKERRQ(ierr);
494e5c89e4eSSatish Balay   if (flg1 || flg2) {
495e5c89e4eSSatish Balay     PetscMPIInt    size;
496e5c89e4eSSatish Balay     PetscInt       lsize,*nodes;
497e5c89e4eSSatish Balay     MPI_Errhandler err_handler;
498e5c89e4eSSatish Balay     /*
499e5c89e4eSSatish Balay        we have to make sure that all processors have opened
500e5c89e4eSSatish Balay        connections to all other processors, otherwise once the
501e5c89e4eSSatish Balay        debugger has stated it is likely to receive a SIGUSR1
502e5c89e4eSSatish Balay        and kill the program.
503e5c89e4eSSatish Balay     */
504e5c89e4eSSatish Balay     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
505e5c89e4eSSatish Balay     if (size > 2) {
506533163c2SBarry Smith       PetscMPIInt dummy = 0;
507e5c89e4eSSatish Balay       MPI_Status  status;
508e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
509e5c89e4eSSatish Balay         if (rank != i) {
510e5c89e4eSSatish Balay           ierr = MPI_Send(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD);CHKERRQ(ierr);
511e5c89e4eSSatish Balay         }
512e5c89e4eSSatish Balay       }
513e5c89e4eSSatish Balay       for (i=0; i<size; i++) {
514e5c89e4eSSatish Balay         if (rank != i) {
515e5c89e4eSSatish Balay           ierr = MPI_Recv(&dummy,1,MPI_INT,i,109,PETSC_COMM_WORLD,&status);CHKERRQ(ierr);
516e5c89e4eSSatish Balay         }
517e5c89e4eSSatish Balay       }
518e5c89e4eSSatish Balay     }
519e5c89e4eSSatish Balay     /* check if this processor node should be in debugger */
520e5c89e4eSSatish Balay     ierr  = PetscMalloc(size*sizeof(PetscInt),&nodes);CHKERRQ(ierr);
521e5c89e4eSSatish Balay     lsize = size;
522e5c89e4eSSatish Balay     ierr  = PetscOptionsGetIntArray(PETSC_NULL,"-debugger_nodes",nodes,&lsize,&flag);CHKERRQ(ierr);
523e5c89e4eSSatish Balay     if (flag) {
524e5c89e4eSSatish Balay       for (i=0; i<lsize; i++) {
525e5c89e4eSSatish Balay         if (nodes[i] == rank) { flag = PETSC_FALSE; break; }
526e5c89e4eSSatish Balay       }
527e5c89e4eSSatish Balay     }
528e5c89e4eSSatish Balay     if (!flag) {
529e5c89e4eSSatish Balay       ierr = PetscSetDebuggerFromString(string);CHKERRQ(ierr);
530e5c89e4eSSatish Balay       ierr = PetscPushErrorHandler(PetscAbortErrorHandler,0);CHKERRQ(ierr);
531e5c89e4eSSatish Balay       if (flg1) {
532e5c89e4eSSatish Balay         ierr = PetscAttachDebugger();CHKERRQ(ierr);
533e5c89e4eSSatish Balay       } else {
534e5c89e4eSSatish Balay         ierr = PetscStopForDebugger();CHKERRQ(ierr);
535e5c89e4eSSatish Balay       }
536e5c89e4eSSatish Balay       ierr = MPI_Errhandler_create((MPI_Handler_function*)Petsc_MPI_AbortOnError,&err_handler);CHKERRQ(ierr);
537e5c89e4eSSatish Balay       ierr = MPI_Errhandler_set(comm,err_handler);CHKERRQ(ierr);
538e5c89e4eSSatish Balay     }
539e5c89e4eSSatish Balay     ierr = PetscFree(nodes);CHKERRQ(ierr);
540e5c89e4eSSatish Balay   }
541e5c89e4eSSatish Balay 
542e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-on_error_emacs",emacsmachinename,128,&flg1);CHKERRQ(ierr);
543cb9801acSJed Brown   if (flg1 && !rank) {ierr = PetscPushErrorHandler(PetscEmacsClientErrorHandler,emacsmachinename);CHKERRQ(ierr);}
544e5c89e4eSSatish Balay 
54593ba235fSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
54622b84c2fSbcordonn   /*
54722b84c2fSbcordonn     Activates new sockets for zope if needed
54822b84c2fSbcordonn   */
54984ab5442Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-zope", &flgz);CHKERRQ(ierr);
550d8c6e182Sbcordonn   ierr = PetscOptionsHasName(PETSC_NULL,"-nostdout", &flgzout);CHKERRQ(ierr);
5516dc8fec2Sbcordonn   if (flgz){
55222b84c2fSbcordonn     int  sockfd;
553f1384234SBarry Smith     char hostname[256];
55422b84c2fSbcordonn     char username[256];
5556dc8fec2Sbcordonn     int  remoteport = 9999;
5569c4c166aSBarry Smith 
55784ab5442Sbcordonn     ierr = PetscOptionsGetString(PETSC_NULL, "-zope", hostname, 256, &flgz);CHKERRQ(ierr);
55884ab5442Sbcordonn     if (!hostname[0]){
5599c4c166aSBarry Smith       ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
5609c4c166aSBarry Smith     }
56122b84c2fSbcordonn     ierr = PetscOpenSocket(hostname, remoteport, &sockfd);CHKERRQ(ierr);
5629c4c166aSBarry Smith     ierr = PetscGetUserName(username, 256);CHKERRQ(ierr);
56322b84c2fSbcordonn     PETSC_ZOPEFD = fdopen(sockfd, "w");
56422b84c2fSbcordonn     if (flgzout){
56522b84c2fSbcordonn       PETSC_STDOUT = PETSC_ZOPEFD;
566606f100bSbcordonn       fprintf(PETSC_STDOUT, "<<<user>>> %s\n",username);
5676dc8fec2Sbcordonn       fprintf(PETSC_STDOUT, "<<<start>>>");
5689c4c166aSBarry Smith     } else {
569d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<user>>> %s\n",username);
570d8c6e182Sbcordonn       fprintf(PETSC_ZOPEFD, "<<<start>>>");
5719c4c166aSBarry Smith     }
5729c4c166aSBarry Smith   }
57393ba235fSBarry Smith #endif
574ffc871a5SBarry Smith #if defined(PETSC_USE_SERVER)
575ffc871a5SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-server", &flgz);CHKERRQ(ierr);
576ffc871a5SBarry Smith   if (flgz){
577ffc871a5SBarry Smith     PetscInt port = PETSC_DECIDE;
578ffc871a5SBarry Smith     ierr = PetscOptionsGetInt(PETSC_NULL,"-server",&port,PETSC_NULL);CHKERRQ(ierr);
579ffc871a5SBarry Smith     ierr = PetscWebServe(PETSC_COMM_WORLD,(int)port);CHKERRQ(ierr);
580ffc871a5SBarry Smith   }
581ffc871a5SBarry Smith #endif
5826dc8fec2Sbcordonn 
583e5c89e4eSSatish Balay   /*
584e5c89e4eSSatish Balay         Setup profiling and logging
585e5c89e4eSSatish Balay   */
5866cf91177SBarry Smith #if defined (PETSC_USE_INFO)
5878bb29257SSatish Balay   {
588e5c89e4eSSatish Balay     char logname[PETSC_MAX_PATH_LEN]; logname[0] = 0;
5896cf91177SBarry Smith     ierr = PetscOptionsGetString(PETSC_NULL,"-info",logname,250,&flg1);CHKERRQ(ierr);
5908bb29257SSatish Balay     if (flg1 && logname[0]) {
591fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,logname);CHKERRQ(ierr);
5928bb29257SSatish Balay     } else if (flg1) {
593fcc2139eSBarry Smith       ierr = PetscInfoAllow(PETSC_TRUE,PETSC_NULL);CHKERRQ(ierr);
594e5c89e4eSSatish Balay     }
595e5c89e4eSSatish Balay   }
596865f6aa8SSatish Balay #endif
597865f6aa8SSatish Balay #if defined(PETSC_USE_LOG)
598865f6aa8SSatish Balay   mname[0] = 0;
599f3dea69dSBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-history",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
600865f6aa8SSatish Balay   if (flg1) {
601865f6aa8SSatish Balay     if (mname[0]) {
602f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(mname,&petsc_history);CHKERRQ(ierr);
603865f6aa8SSatish Balay     } else {
604f3dea69dSBarry Smith       ierr = PetscOpenHistoryFile(0,&petsc_history);CHKERRQ(ierr);
605865f6aa8SSatish Balay     }
606865f6aa8SSatish Balay   }
607e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
60890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
609fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_mpe",&flg1);CHKERRQ(ierr);
610e5c89e4eSSatish Balay   if (flg1) PetscLogMPEBegin();
611e5c89e4eSSatish Balay #endif
61290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
61390d69ab7SBarry Smith   flg2 = PETSC_FALSE;
61490d69ab7SBarry Smith   flg3 = PETSC_FALSE;
615acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log_all",&flg1,PETSC_NULL);CHKERRQ(ierr);
616acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-log",&flg2,PETSC_NULL);CHKERRQ(ierr);
617d44e083bSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
6189f7b6320SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary_python",&flg4);CHKERRQ(ierr);
619e5c89e4eSSatish Balay   if (flg1)                      {  ierr = PetscLogAllBegin();CHKERRQ(ierr); }
6209f7b6320SBarry Smith   else if (flg2 || flg3 || flg4) {  ierr = PetscLogBegin();CHKERRQ(ierr);}
621e5c89e4eSSatish Balay 
622e5c89e4eSSatish Balay   ierr = PetscOptionsGetString(PETSC_NULL,"-log_trace",mname,250,&flg1);CHKERRQ(ierr);
623e5c89e4eSSatish Balay   if (flg1) {
624e5c89e4eSSatish Balay     char name[PETSC_MAX_PATH_LEN],fname[PETSC_MAX_PATH_LEN];
625e5c89e4eSSatish Balay     FILE *file;
626e5c89e4eSSatish Balay     if (mname[0]) {
627e5c89e4eSSatish Balay       sprintf(name,"%s.%d",mname,rank);
628e5c89e4eSSatish Balay       ierr = PetscFixFilename(name,fname);CHKERRQ(ierr);
629e5c89e4eSSatish Balay       file = fopen(fname,"w");
630f3dea69dSBarry Smith       if (!file) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Unable to open trace file: %s",fname);
631e5c89e4eSSatish Balay     } else {
632da9f1d6bSBarry Smith       file = PETSC_STDOUT;
633e5c89e4eSSatish Balay     }
634e5c89e4eSSatish Balay     ierr = PetscLogTraceBegin(file);CHKERRQ(ierr);
635e5c89e4eSSatish Balay   }
636e5c89e4eSSatish Balay #endif
637e5c89e4eSSatish Balay 
638e5c89e4eSSatish Balay   /*
639e5c89e4eSSatish Balay       Setup building of stack frames for all function calls
640e5c89e4eSSatish Balay   */
64163d6bff0SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_USE_PTHREAD)
642e5c89e4eSSatish Balay   ierr = PetscStackCreate();CHKERRQ(ierr);
643e5c89e4eSSatish Balay #endif
644e5c89e4eSSatish Balay 
645acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-options_gui",&PetscOptionsPublish,PETSC_NULL);CHKERRQ(ierr);
646e5c89e4eSSatish Balay 
647af359df3SBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
648af359df3SBarry Smith   /*
649af359df3SBarry Smith       Determine whether user specified maximum number of threads
650af359df3SBarry Smith    */
651af359df3SBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-thread_max",&PetscMaxThreads,PETSC_NULL);CHKERRQ(ierr);
652af359df3SBarry Smith 
653b154f58aSKerry Stevens   ierr = PetscOptionsHasName(PETSC_NULL,"-main",&flg1);CHKERRQ(ierr);
654b154f58aSKerry Stevens   if(flg1) {
655b154f58aSKerry Stevens     cpu_set_t mset;
656b154f58aSKerry Stevens     int icorr,ncorr = get_nprocs();
657b154f58aSKerry Stevens     ierr = PetscOptionsGetInt(PETSC_NULL,"-main",&icorr,PETSC_NULL);CHKERRQ(ierr);
658b154f58aSKerry Stevens     CPU_ZERO(&mset);
659b154f58aSKerry Stevens     CPU_SET(icorr%ncorr,&mset);
660b154f58aSKerry Stevens     sched_setaffinity(0,sizeof(cpu_set_t),&mset);
661b154f58aSKerry Stevens   }
662b154f58aSKerry Stevens 
663af359df3SBarry Smith   /*
664af359df3SBarry Smith       Determine whether to use thread pool
665af359df3SBarry Smith    */
666af359df3SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-use_thread_pool",&flg1);CHKERRQ(ierr);
667af359df3SBarry Smith   if (flg1) {
668af359df3SBarry Smith     PetscUseThreadPool = PETSC_TRUE;
669af359df3SBarry Smith     PetscInt N_CORES = get_nprocs();
670af359df3SBarry Smith     ThreadCoreAffinity = (int*)malloc(N_CORES*sizeof(int));
671af359df3SBarry Smith     char tstr[9];
672af359df3SBarry Smith     char tbuf[2];
673af359df3SBarry Smith     strcpy(tstr,"-thread");
674af359df3SBarry Smith     for(i=0;i<PetscMaxThreads;i++) {
675af359df3SBarry Smith       ThreadCoreAffinity[i] = i;
676af359df3SBarry Smith       sprintf(tbuf,"%d",i);
677af359df3SBarry Smith       strcat(tstr,tbuf);
678af359df3SBarry Smith       ierr = PetscOptionsHasName(PETSC_NULL,tstr,&flg1);CHKERRQ(ierr);
679af359df3SBarry Smith       if(flg1) {
680af359df3SBarry Smith         ierr = PetscOptionsGetInt(PETSC_NULL,tstr,&ThreadCoreAffinity[i],PETSC_NULL);CHKERRQ(ierr);
681af359df3SBarry Smith         ThreadCoreAffinity[i] = ThreadCoreAffinity[i]%N_CORES; /* check on the user */
682af359df3SBarry Smith       }
683af359df3SBarry Smith       tstr[7] = '\0';
684af359df3SBarry Smith     }
685af359df3SBarry Smith     /* get the thread pool type */
686af359df3SBarry Smith     PetscInt ipool = 0;
687af359df3SBarry Smith     const char *choices[4] = {"true","tree","main","chain"};
688af359df3SBarry Smith 
689af359df3SBarry Smith     ierr = PetscOptionsGetEList(PETSC_NULL,"-use_thread_pool",choices,4,&ipool,PETSC_NULL);CHKERRQ(ierr);
690af359df3SBarry Smith     switch(ipool) {
691af359df3SBarry Smith     case 1:
692af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Tree;
693af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Tree;
694af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Tree;
695af359df3SBarry Smith       MainWait              = &MainWait_Tree;
696af359df3SBarry Smith       MainJob               = &MainJob_Tree;
697af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using tree thread pool\n");
698af359df3SBarry Smith       break;
699af359df3SBarry Smith     case 2:
700af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Main;
701af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Main;
702af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Main;
703af359df3SBarry Smith       MainWait              = &MainWait_Main;
704af359df3SBarry Smith       MainJob               = &MainJob_Main;
705af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using main thread pool\n");
706af359df3SBarry Smith       break;
707af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
708af359df3SBarry Smith     case 3:
709af359df3SBarry Smith #else
710af359df3SBarry Smith     default:
711af359df3SBarry Smith #endif
712af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_Chain;
713af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_Chain;
714af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_Chain;
715af359df3SBarry Smith       MainWait              = &MainWait_Chain;
716af359df3SBarry Smith       MainJob               = &MainJob_Chain;
717af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using chain thread pool\n");
718af359df3SBarry Smith       break;
719af359df3SBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
720af359df3SBarry Smith     default:
721af359df3SBarry Smith       PetscThreadFunc       = &PetscThreadFunc_True;
722af359df3SBarry Smith       PetscThreadInitialize = &PetscThreadInitialize_True;
723af359df3SBarry Smith       PetscThreadFinalize   = &PetscThreadFinalize_True;
724af359df3SBarry Smith       MainWait              = &MainWait_True;
725af359df3SBarry Smith       MainJob               = &MainJob_True;
726af359df3SBarry Smith       PetscInfo(PETSC_NULL,"Using true thread pool\n");
727af359df3SBarry Smith       break;
728af359df3SBarry Smith #endif
729af359df3SBarry Smith     }
730af359df3SBarry Smith     PetscThreadInitialize(PetscMaxThreads);
731683509dcSKerry Stevens   } else {
732683509dcSKerry Stevens     //need to define these in the case on 'no threads' or 'thread create/destroy'
733683509dcSKerry Stevens     //could take any of the above versions
734683509dcSKerry Stevens     MainJob               = &MainJob_Spawn;
735af359df3SBarry Smith   }
736af359df3SBarry Smith #endif
737e5c89e4eSSatish Balay   /*
738e5c89e4eSSatish Balay        Print basic help message
739e5c89e4eSSatish Balay   */
740e5c89e4eSSatish Balay   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg1);CHKERRQ(ierr);
741e5c89e4eSSatish Balay   if (flg1) {
742e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"Options for all PETSc programs:\n");CHKERRQ(ierr);
743301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -help: prints help method for each option\n");CHKERRQ(ierr);
744301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm," -on_error_abort: cause an abort when an error is detected. Useful \n ");CHKERRQ(ierr);
745301d30feSBarry Smith     ierr = (*PetscHelpPrintf)(comm,"       only when run in the debugger\n");CHKERRQ(ierr);
746e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_attach_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
747e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start the debugger in new xterm\n");CHKERRQ(ierr);
748e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       unless noxterm is given\n");CHKERRQ(ierr);
749e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -start_in_debugger [gdb,dbx,xxgdb,ups,noxterm]\n");CHKERRQ(ierr);
750e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"       start all processes in the debugger\n");CHKERRQ(ierr);
751e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -on_error_emacs <machinename>\n");CHKERRQ(ierr);
752e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"    emacs jumps to error file\n");CHKERRQ(ierr);
753e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_nodes [n1,n2,..] Nodes to start in debugger\n");CHKERRQ(ierr);
754e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -debugger_pause [m] : delay (in seconds) to attach debugger\n");CHKERRQ(ierr);
755e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -stop_for_debugger : prints message on how to attach debugger manually\n");CHKERRQ(ierr);
756e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"                      waits the delay for you to attach\n");CHKERRQ(ierr);
757e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -display display: Location where graphics and debuggers are displayed\n");CHKERRQ(ierr);
758e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -no_signal_handler: do not trap error signals\n");CHKERRQ(ierr);
759e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -mpi_return_on_error: MPI returns error code, rather than abort on internal error\n");CHKERRQ(ierr);
760e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -fp_trap: stop on floating point exceptions\n");CHKERRQ(ierr);
761e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"           note on IBM RS6000 this slows run greatly\n");CHKERRQ(ierr);
762e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_dump <optional filename>: dump list of unfreed memory at conclusion\n");CHKERRQ(ierr);
763e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc: use our error checking malloc\n");CHKERRQ(ierr);
764e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc no: don't use error checking malloc\n");CHKERRQ(ierr);
7654161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_info: prints total memory usage\n");CHKERRQ(ierr);
7664161f2a3SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -malloc_log: keeps log of all memory allocations\n");CHKERRQ(ierr);
767e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -malloc_debug: enables extended checking for memory corruption\n");CHKERRQ(ierr);
768e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_table: dump list of options inputted\n");CHKERRQ(ierr);
769e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left: dump list of unused options\n");CHKERRQ(ierr);
770e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_left no: don't dump list of unused options\n");CHKERRQ(ierr);
771e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -tmp tmpdir: alternative /tmp directory\n");CHKERRQ(ierr);
772e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -shared_tmp: tmp directory is shared by all processors\n");CHKERRQ(ierr);
773a8c7a070SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -not_shared_tmp: each processor has separate tmp directory\n");CHKERRQ(ierr);
774e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -memory_info: print memory usage at end of run\n");CHKERRQ(ierr);
77540ab9619SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -server <port>: Run PETSc webserver (default port is 8080) see PetscWebServe()\n");CHKERRQ(ierr);
776e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
777e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -get_total_flops: total flops over all processors\n");CHKERRQ(ierr);
7787ef452c0SMatthew G Knepley     ierr = (*PetscHelpPrintf)(comm," -log[_all _summary _summary_python]: logging objects and events\n");CHKERRQ(ierr);
779e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_trace [filename]: prints trace of all PETSc calls\n");CHKERRQ(ierr);
780e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
781e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -log_mpe: Also create logfile viewable through upshot\n");CHKERRQ(ierr);
782e5c89e4eSSatish Balay #endif
7836cf91177SBarry Smith     ierr = (*PetscHelpPrintf)(comm," -info <optional filename>: print informative messages about the calculations\n");CHKERRQ(ierr);
784e5c89e4eSSatish Balay #endif
785e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -v: prints PETSc version number and release date\n");CHKERRQ(ierr);
786e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -options_file <file>: reads options from file\n");CHKERRQ(ierr);
787e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm," -petsc_sleep n: sleeps n seconds before running program\n");CHKERRQ(ierr);
788e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
789e5c89e4eSSatish Balay   }
790e5c89e4eSSatish Balay 
791a6d0e24fSJed Brown   ierr = PetscOptionsGetReal(PETSC_NULL,"-petsc_sleep",&si,&flg1);CHKERRQ(ierr);
792e5c89e4eSSatish Balay   if (flg1) {
793e5c89e4eSSatish Balay     ierr = PetscSleep(si);CHKERRQ(ierr);
794e5c89e4eSSatish Balay   }
795e5c89e4eSSatish Balay 
7966cf91177SBarry Smith   ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
797e5c89e4eSSatish Balay   ierr = PetscStrstr(mname,"null",&f);CHKERRQ(ierr);
798e5c89e4eSSatish Balay   if (f) {
7996cf91177SBarry Smith     ierr = PetscInfoDeactivateClass(PETSC_NULL);CHKERRQ(ierr);
800e5c89e4eSSatish Balay   }
801827f890bSBarry Smith 
8028154be41SBarry Smith #if defined(PETSC_HAVE_CUSP)
803c97f9302SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-log_summary",&flg3);CHKERRQ(ierr);
80473113deaSBarry Smith   if (flg3) flg1 = PETSC_TRUE;
80573113deaSBarry Smith   else flg1 = PETSC_FALSE;
8068154be41SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-cusp_synchronize",&flg1,PETSC_NULL);CHKERRQ(ierr);
8078154be41SBarry Smith   if (flg1) synchronizeCUSP = PETSC_TRUE;
808bab1f7e6SVictor Minden #endif
809192daf7cSBarry Smith 
810e5c89e4eSSatish Balay   PetscFunctionReturn(0);
811e5c89e4eSSatish Balay }
812df413903SBarry Smith 
813ba61063dSBarry Smith #if defined(PETSC_USE_PTHREAD_CLASSES)
814ba61063dSBarry Smith 
81551d315f7SKerry Stevens /**** 'Tree' Thread Pool Functions ****/
81651d315f7SKerry Stevens void* PetscThreadFunc_Tree(void* arg) {
81751d315f7SKerry Stevens   PetscErrorCode iterr;
81851d315f7SKerry Stevens   int icorr,ierr;
81951d315f7SKerry Stevens   int* pId = (int*)arg;
82051d315f7SKerry Stevens   int ThreadId = *pId,Mary = 2,i,SubWorker;
82151d315f7SKerry Stevens   PetscBool PeeOn;
82251d315f7SKerry Stevens   cpu_set_t mset;
8239e800a48SKerry Stevens   //printf("Thread %d In Tree Thread Function\n",ThreadId);
82451d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
82551d315f7SKerry Stevens   CPU_ZERO(&mset);
82651d315f7SKerry Stevens   CPU_SET(icorr,&mset);
82751d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
82851d315f7SKerry Stevens 
82951d315f7SKerry Stevens   if((Mary*ThreadId+1)>(PetscMaxThreads-1)) {
83051d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
83151d315f7SKerry Stevens   }
83251d315f7SKerry Stevens   else {
83351d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
83451d315f7SKerry Stevens   }
83551d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
836ba61063dSBarry Smith     /* check your subordinates, wait for them to be ready */
83751d315f7SKerry Stevens     for(i=1;i<=Mary;i++) {
83851d315f7SKerry Stevens       SubWorker = Mary*ThreadId+i;
83951d315f7SKerry Stevens       if(SubWorker<PetscMaxThreads) {
84051d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
84151d315f7SKerry Stevens         while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE) {
842ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
843ba61063dSBarry Smith            upon return, has the lock */
84451d315f7SKerry Stevens           ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
84551d315f7SKerry Stevens         }
84651d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
84751d315f7SKerry Stevens       }
84851d315f7SKerry Stevens     }
849ba61063dSBarry Smith     /* your subordinates are now ready */
85051d315f7SKerry Stevens   }
85151d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
852ba61063dSBarry Smith   /* update your ready status */
85351d315f7SKerry Stevens   *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
85451d315f7SKerry Stevens   if(ThreadId==0) {
85551d315f7SKerry Stevens     job_tree.eJobStat = JobCompleted;
856ba61063dSBarry Smith     /* ignal main */
85751d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
85851d315f7SKerry Stevens   }
85951d315f7SKerry Stevens   else {
860ba61063dSBarry Smith     /* tell your boss that you're ready to work */
86151d315f7SKerry Stevens     ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
86251d315f7SKerry Stevens   }
863ba61063dSBarry Smith   /* the while loop needs to have an exit
864ba61063dSBarry Smith   the 'main' thread can terminate all the threads by performing a broadcast
865ba61063dSBarry Smith    and calling FuncFinish */
86651d315f7SKerry Stevens   while(PetscThreadGo) {
867ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
868ba61063dSBarry Smith       waiting when you don't have to causes problems
869ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
87051d315f7SKerry Stevens     while(*(job_tree.arrThreadReady[ThreadId])==PETSC_TRUE) {
871ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
872ba61063dSBarry Smith        upon return, has the lock */
87351d315f7SKerry Stevens         ierr = pthread_cond_wait(job_tree.cond2array[ThreadId],job_tree.mutexarray[ThreadId]);
87451d315f7SKerry Stevens 	*(job_tree.arrThreadStarted[ThreadId]) = PETSC_TRUE;
87551d315f7SKerry Stevens 	*(job_tree.arrThreadReady[ThreadId])   = PETSC_FALSE;
87651d315f7SKerry Stevens     }
87751d315f7SKerry Stevens     if(ThreadId==0) {
87851d315f7SKerry Stevens       job_tree.startJob = PETSC_FALSE;
87951d315f7SKerry Stevens       job_tree.eJobStat = ThreadsWorking;
88051d315f7SKerry Stevens     }
88151d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_tree.mutexarray[ThreadId]);
88251d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
883ba61063dSBarry Smith       /* tell your subordinates it's time to get to work */
88451d315f7SKerry Stevens       for(i=1; i<=Mary; i++) {
88551d315f7SKerry Stevens 	SubWorker = Mary*ThreadId+i;
88651d315f7SKerry Stevens         if(SubWorker<PetscMaxThreads) {
88751d315f7SKerry Stevens           ierr = pthread_cond_signal(job_tree.cond2array[SubWorker]);
88851d315f7SKerry Stevens         }
88951d315f7SKerry Stevens       }
89051d315f7SKerry Stevens     }
891ba61063dSBarry Smith     /* do your job */
89251d315f7SKerry Stevens     if(job_tree.pdata==NULL) {
89351d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata);
89451d315f7SKerry Stevens     }
89551d315f7SKerry Stevens     else {
89651d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_tree.pfunc(job_tree.pdata[ThreadId]);
89751d315f7SKerry Stevens     }
89851d315f7SKerry Stevens     if(iterr!=0) {
89951d315f7SKerry Stevens       ithreaderr = 1;
90051d315f7SKerry Stevens     }
90151d315f7SKerry Stevens     if(PetscThreadGo) {
902ba61063dSBarry Smith       /* reset job, get ready for more */
90351d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
904ba61063dSBarry Smith         /* check your subordinates, waiting for them to be ready
905ba61063dSBarry Smith          how do you know for a fact that a given subordinate has actually started? */
90651d315f7SKerry Stevens 	for(i=1;i<=Mary;i++) {
90751d315f7SKerry Stevens 	  SubWorker = Mary*ThreadId+i;
90851d315f7SKerry Stevens           if(SubWorker<PetscMaxThreads) {
90951d315f7SKerry Stevens             ierr = pthread_mutex_lock(job_tree.mutexarray[SubWorker]);
91051d315f7SKerry Stevens             while(*(job_tree.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_tree.arrThreadStarted[SubWorker])==PETSC_FALSE) {
911ba61063dSBarry Smith               /* upon entry, automically releases the lock and blocks
912ba61063dSBarry Smith                upon return, has the lock */
91351d315f7SKerry Stevens               ierr = pthread_cond_wait(job_tree.cond1array[SubWorker],job_tree.mutexarray[SubWorker]);
91451d315f7SKerry Stevens             }
91551d315f7SKerry Stevens             ierr = pthread_mutex_unlock(job_tree.mutexarray[SubWorker]);
91651d315f7SKerry Stevens           }
91751d315f7SKerry Stevens 	}
918ba61063dSBarry Smith         /* your subordinates are now ready */
91951d315f7SKerry Stevens       }
92051d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_tree.mutexarray[ThreadId]);
92151d315f7SKerry Stevens       *(job_tree.arrThreadReady[ThreadId]) = PETSC_TRUE;
92251d315f7SKerry Stevens       if(ThreadId==0) {
923ba61063dSBarry Smith 	job_tree.eJobStat = JobCompleted; /* oot thread: last thread to complete, guaranteed! */
924ba61063dSBarry Smith         /* root thread signals 'main' */
92551d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
92651d315f7SKerry Stevens       }
92751d315f7SKerry Stevens       else {
928ba61063dSBarry Smith         /* signal your boss before you go to sleep */
92951d315f7SKerry Stevens         ierr = pthread_cond_signal(job_tree.cond1array[ThreadId]);
93051d315f7SKerry Stevens       }
93151d315f7SKerry Stevens     }
93251d315f7SKerry Stevens   }
93351d315f7SKerry Stevens   return NULL;
93451d315f7SKerry Stevens }
93551d315f7SKerry Stevens 
93651d315f7SKerry Stevens #undef __FUNCT__
93751d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Tree"
93851d315f7SKerry Stevens void* PetscThreadInitialize_Tree(PetscInt N) {
93951d315f7SKerry Stevens   PetscInt i,ierr;
94051d315f7SKerry Stevens   int status;
94151d315f7SKerry Stevens 
94251d315f7SKerry Stevens   if(PetscUseThreadPool) {
94351d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
94451d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
94551d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
94651d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
94751d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
94851d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
94951d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
95051d315f7SKerry Stevens     job_tree.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
95151d315f7SKerry Stevens     job_tree.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95251d315f7SKerry Stevens     job_tree.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
95351d315f7SKerry Stevens     job_tree.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
95451d315f7SKerry Stevens     job_tree.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
955ba61063dSBarry Smith     /* initialize job structure */
95651d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
95751d315f7SKerry Stevens       job_tree.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
95851d315f7SKerry Stevens       job_tree.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
95951d315f7SKerry Stevens       job_tree.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
96051d315f7SKerry Stevens       job_tree.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
96151d315f7SKerry Stevens       job_tree.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
96251d315f7SKerry Stevens     }
96351d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
96451d315f7SKerry Stevens       ierr = pthread_mutex_init(job_tree.mutexarray[i],NULL);
96551d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond1array[i],NULL);
96651d315f7SKerry Stevens       ierr = pthread_cond_init(job_tree.cond2array[i],NULL);
96751d315f7SKerry Stevens       *(job_tree.arrThreadStarted[i])  = PETSC_FALSE;
96851d315f7SKerry Stevens       *(job_tree.arrThreadReady[i])    = PETSC_FALSE;
96951d315f7SKerry Stevens     }
97051d315f7SKerry Stevens     job_tree.pfunc = NULL;
97151d315f7SKerry Stevens     job_tree.pdata = (void**)malloc(N*sizeof(void*));
97251d315f7SKerry Stevens     job_tree.startJob = PETSC_FALSE;
97351d315f7SKerry Stevens     job_tree.eJobStat = JobInitiated;
97451d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
975ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
97651d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
977ba61063dSBarry Smith     /* create threads */
97851d315f7SKerry Stevens     for(i=0; i<N; i++) {
97951d315f7SKerry Stevens       pVal[i] = i;
98051d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
981ba61063dSBarry Smith       /* should check status */
98251d315f7SKerry Stevens     }
98351d315f7SKerry Stevens   }
98451d315f7SKerry Stevens   return NULL;
98551d315f7SKerry Stevens }
98651d315f7SKerry Stevens 
98751d315f7SKerry Stevens #undef __FUNCT__
98851d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Tree"
98951d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Tree() {
99051d315f7SKerry Stevens   int i,ierr;
99151d315f7SKerry Stevens   void* jstatus;
99251d315f7SKerry Stevens 
99351d315f7SKerry Stevens   PetscFunctionBegin;
99451d315f7SKerry Stevens 
99551d315f7SKerry Stevens   if(PetscUseThreadPool) {
996ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
997ba61063dSBarry Smith     /* join the threads */
99851d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
99951d315f7SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1000ba61063dSBarry Smith       /* do error checking*/
100151d315f7SKerry Stevens     }
100251d315f7SKerry Stevens     free(PetscThreadPoint);
100351d315f7SKerry Stevens     free(arrmutex);
100451d315f7SKerry Stevens     free(arrcond1);
100551d315f7SKerry Stevens     free(arrcond2);
100651d315f7SKerry Stevens     free(arrstart);
100751d315f7SKerry Stevens     free(arrready);
100851d315f7SKerry Stevens     free(job_tree.pdata);
100951d315f7SKerry Stevens     free(pVal);
101051d315f7SKerry Stevens   }
101151d315f7SKerry Stevens   else {
101251d315f7SKerry Stevens   }
101351d315f7SKerry Stevens   PetscFunctionReturn(0);
101451d315f7SKerry Stevens }
101551d315f7SKerry Stevens 
101651d315f7SKerry Stevens #undef __FUNCT__
101751d315f7SKerry Stevens #define __FUNCT__ "MainWait_Tree"
101851d315f7SKerry Stevens void MainWait_Tree() {
101951d315f7SKerry Stevens   int ierr;
102051d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_tree.mutexarray[0]);
102151d315f7SKerry Stevens   while(job_tree.eJobStat<JobCompleted||job_tree.startJob==PETSC_TRUE) {
102251d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_tree.mutexarray[0]);
102351d315f7SKerry Stevens   }
102451d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_tree.mutexarray[0]);
102551d315f7SKerry Stevens }
102651d315f7SKerry Stevens 
102751d315f7SKerry Stevens #undef __FUNCT__
102851d315f7SKerry Stevens #define __FUNCT__ "MainJob_Tree"
102951d315f7SKerry Stevens PetscErrorCode MainJob_Tree(void* (*pFunc)(void*),void** data,PetscInt n) {
103051d315f7SKerry Stevens   int i,ierr;
103151d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
103236d20dc5SKerry Stevens 
103351d315f7SKerry Stevens   MainWait();
103451d315f7SKerry Stevens   job_tree.pfunc = pFunc;
103551d315f7SKerry Stevens   job_tree.pdata = data;
103651d315f7SKerry Stevens   job_tree.startJob = PETSC_TRUE;
103751d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
103851d315f7SKerry Stevens     *(job_tree.arrThreadStarted[i]) = PETSC_FALSE;
103951d315f7SKerry Stevens   }
104051d315f7SKerry Stevens   job_tree.eJobStat = JobInitiated;
104151d315f7SKerry Stevens   ierr = pthread_cond_signal(job_tree.cond2array[0]);
104251d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1043ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
104451d315f7SKerry Stevens   }
104536d20dc5SKerry Stevens 
104651d315f7SKerry Stevens   if(ithreaderr) {
104751d315f7SKerry Stevens     ijoberr = ithreaderr;
104851d315f7SKerry Stevens   }
104951d315f7SKerry Stevens   return ijoberr;
105051d315f7SKerry Stevens }
105151d315f7SKerry Stevens /****  ****/
105251d315f7SKerry Stevens 
105351d315f7SKerry Stevens /**** 'Main' Thread Pool Functions ****/
105451d315f7SKerry Stevens void* PetscThreadFunc_Main(void* arg) {
105551d315f7SKerry Stevens   PetscErrorCode iterr;
105651d315f7SKerry Stevens   int icorr,ierr;
105751d315f7SKerry Stevens   int* pId = (int*)arg;
105851d315f7SKerry Stevens   int ThreadId = *pId;
105951d315f7SKerry Stevens   cpu_set_t mset;
10609e800a48SKerry Stevens   //printf("Thread %d In Main Thread Function\n",ThreadId);
106151d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
106251d315f7SKerry Stevens   CPU_ZERO(&mset);
106351d315f7SKerry Stevens   CPU_SET(icorr,&mset);
106451d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
106551d315f7SKerry Stevens 
106651d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
1067ba61063dSBarry Smith   /* update your ready status */
106851d315f7SKerry Stevens   *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1069ba61063dSBarry Smith   /* tell the BOSS that you're ready to work before you go to sleep */
107051d315f7SKerry Stevens   ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
107151d315f7SKerry Stevens 
1072ba61063dSBarry Smith   /* the while loop needs to have an exit
1073ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1074ba61063dSBarry Smith      and calling FuncFinish */
107551d315f7SKerry Stevens   while(PetscThreadGo) {
1076ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1077ba61063dSBarry Smith        waiting when you don't have to causes problems
1078ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
107951d315f7SKerry Stevens     while(*(job_main.arrThreadReady[ThreadId])==PETSC_TRUE) {
1080ba61063dSBarry Smith       /* upon entry, atomically releases the lock and blocks
1081ba61063dSBarry Smith        upon return, has the lock */
108251d315f7SKerry Stevens         ierr = pthread_cond_wait(job_main.cond2array[ThreadId],job_main.mutexarray[ThreadId]);
1083ba61063dSBarry Smith 	/* (job_main.arrThreadReady[ThreadId])   = PETSC_FALSE; */
108451d315f7SKerry Stevens     }
108551d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[ThreadId]);
108651d315f7SKerry Stevens     if(job_main.pdata==NULL) {
108751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata);
108851d315f7SKerry Stevens     }
108951d315f7SKerry Stevens     else {
109051d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_main.pfunc(job_main.pdata[ThreadId]);
109151d315f7SKerry Stevens     }
109251d315f7SKerry Stevens     if(iterr!=0) {
109351d315f7SKerry Stevens       ithreaderr = 1;
109451d315f7SKerry Stevens     }
109551d315f7SKerry Stevens     if(PetscThreadGo) {
1096ba61063dSBarry Smith       /* reset job, get ready for more */
109751d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_main.mutexarray[ThreadId]);
109851d315f7SKerry Stevens       *(job_main.arrThreadReady[ThreadId]) = PETSC_TRUE;
1099ba61063dSBarry Smith       /* tell the BOSS that you're ready to work before you go to sleep */
110051d315f7SKerry Stevens       ierr = pthread_cond_signal(job_main.cond1array[ThreadId]);
110151d315f7SKerry Stevens     }
110251d315f7SKerry Stevens   }
110351d315f7SKerry Stevens   return NULL;
110451d315f7SKerry Stevens }
110551d315f7SKerry Stevens 
110651d315f7SKerry Stevens #undef __FUNCT__
110751d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Main"
110851d315f7SKerry Stevens void* PetscThreadInitialize_Main(PetscInt N) {
110951d315f7SKerry Stevens   PetscInt i,ierr;
111051d315f7SKerry Stevens   int status;
111151d315f7SKerry Stevens 
111251d315f7SKerry Stevens   if(PetscUseThreadPool) {
111351d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
111451d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
111551d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
111651d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
111751d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
111851d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
111951d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
112051d315f7SKerry Stevens     job_main.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
112151d315f7SKerry Stevens     job_main.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112251d315f7SKerry Stevens     job_main.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
112351d315f7SKerry Stevens     job_main.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1124ba61063dSBarry Smith     /* initialize job structure */
112551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
112651d315f7SKerry Stevens       job_main.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
112751d315f7SKerry Stevens       job_main.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
112851d315f7SKerry Stevens       job_main.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
112951d315f7SKerry Stevens       job_main.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
113051d315f7SKerry Stevens     }
113151d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
113251d315f7SKerry Stevens       ierr = pthread_mutex_init(job_main.mutexarray[i],NULL);
113351d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond1array[i],NULL);
113451d315f7SKerry Stevens       ierr = pthread_cond_init(job_main.cond2array[i],NULL);
113551d315f7SKerry Stevens       *(job_main.arrThreadReady[i])    = PETSC_FALSE;
113651d315f7SKerry Stevens     }
113751d315f7SKerry Stevens     job_main.pfunc = NULL;
113851d315f7SKerry Stevens     job_main.pdata = (void**)malloc(N*sizeof(void*));
113951d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1140ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
114151d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1142ba61063dSBarry Smith     /* create threads */
114351d315f7SKerry Stevens     for(i=0; i<N; i++) {
114451d315f7SKerry Stevens       pVal[i] = i;
114551d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1146ba61063dSBarry Smith       /* error check */
114751d315f7SKerry Stevens     }
114851d315f7SKerry Stevens   }
114951d315f7SKerry Stevens   else {
115051d315f7SKerry Stevens   }
115151d315f7SKerry Stevens   return NULL;
115251d315f7SKerry Stevens }
115351d315f7SKerry Stevens 
115451d315f7SKerry Stevens #undef __FUNCT__
115551d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Main"
115651d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Main() {
115751d315f7SKerry Stevens   int i,ierr;
115851d315f7SKerry Stevens   void* jstatus;
115951d315f7SKerry Stevens 
116051d315f7SKerry Stevens   PetscFunctionBegin;
116151d315f7SKerry Stevens 
116251d315f7SKerry Stevens   if(PetscUseThreadPool) {
1163ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1164ba61063dSBarry Smith     /* join the threads */
116551d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
1166ba61063dSBarry Smith       ierr = pthread_join(PetscThreadPoint[i],&jstatus);CHKERRQ(ierr);
116751d315f7SKerry Stevens     }
116851d315f7SKerry Stevens     free(PetscThreadPoint);
116951d315f7SKerry Stevens     free(arrmutex);
117051d315f7SKerry Stevens     free(arrcond1);
117151d315f7SKerry Stevens     free(arrcond2);
117251d315f7SKerry Stevens     free(arrstart);
117351d315f7SKerry Stevens     free(arrready);
117451d315f7SKerry Stevens     free(job_main.pdata);
117551d315f7SKerry Stevens     free(pVal);
117651d315f7SKerry Stevens   }
117751d315f7SKerry Stevens   PetscFunctionReturn(0);
117851d315f7SKerry Stevens }
117951d315f7SKerry Stevens 
118051d315f7SKerry Stevens #undef __FUNCT__
118151d315f7SKerry Stevens #define __FUNCT__ "MainWait_Main"
118251d315f7SKerry Stevens void MainWait_Main() {
118351d315f7SKerry Stevens   int i,ierr;
118451d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
118551d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_main.mutexarray[i]);
118651d315f7SKerry Stevens     while(*(job_main.arrThreadReady[i])==PETSC_FALSE) {
118751d315f7SKerry Stevens       ierr = pthread_cond_wait(job_main.cond1array[i],job_main.mutexarray[i]);
118851d315f7SKerry Stevens     }
118951d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_main.mutexarray[i]);
119051d315f7SKerry Stevens   }
119151d315f7SKerry Stevens }
119251d315f7SKerry Stevens 
119351d315f7SKerry Stevens #undef __FUNCT__
119451d315f7SKerry Stevens #define __FUNCT__ "MainJob_Main"
119551d315f7SKerry Stevens PetscErrorCode MainJob_Main(void* (*pFunc)(void*),void** data,PetscInt n) {
119651d315f7SKerry Stevens   int i,ierr;
119751d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
119836d20dc5SKerry Stevens 
1199ba61063dSBarry Smith   MainWait(); /* you know everyone is waiting to be signalled! */
120051d315f7SKerry Stevens   job_main.pfunc = pFunc;
120151d315f7SKerry Stevens   job_main.pdata = data;
120251d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
1203ba61063dSBarry Smith     *(job_main.arrThreadReady[i]) = PETSC_FALSE; /* why do this?  suppose you get into MainWait first */
120451d315f7SKerry Stevens   }
1205ba61063dSBarry Smith   /* tell the threads to go to work */
120651d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
120751d315f7SKerry Stevens     ierr = pthread_cond_signal(job_main.cond2array[i]);
120851d315f7SKerry Stevens   }
120951d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1210ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
121151d315f7SKerry Stevens   }
121236d20dc5SKerry Stevens 
121351d315f7SKerry Stevens   if(ithreaderr) {
121451d315f7SKerry Stevens     ijoberr = ithreaderr;
121551d315f7SKerry Stevens   }
121651d315f7SKerry Stevens   return ijoberr;
121751d315f7SKerry Stevens }
121851d315f7SKerry Stevens /****  ****/
121951d315f7SKerry Stevens 
122051d315f7SKerry Stevens /**** Chain Thread Functions ****/
122151d315f7SKerry Stevens void* PetscThreadFunc_Chain(void* arg) {
122251d315f7SKerry Stevens   PetscErrorCode iterr;
122351d315f7SKerry Stevens   int icorr,ierr;
122451d315f7SKerry Stevens   int* pId = (int*)arg;
122551d315f7SKerry Stevens   int ThreadId = *pId;
122651d315f7SKerry Stevens   int SubWorker = ThreadId + 1;
122751d315f7SKerry Stevens   PetscBool PeeOn;
122851d315f7SKerry Stevens   cpu_set_t mset;
12299e800a48SKerry Stevens   //printf("Thread %d In Chain Thread Function\n",ThreadId);
123051d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
123151d315f7SKerry Stevens   CPU_ZERO(&mset);
123251d315f7SKerry Stevens   CPU_SET(icorr,&mset);
123351d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
123451d315f7SKerry Stevens 
123551d315f7SKerry Stevens   if(ThreadId==(PetscMaxThreads-1)) {
123651d315f7SKerry Stevens     PeeOn = PETSC_TRUE;
123751d315f7SKerry Stevens   }
123851d315f7SKerry Stevens   else {
123951d315f7SKerry Stevens     PeeOn = PETSC_FALSE;
124051d315f7SKerry Stevens   }
124151d315f7SKerry Stevens   if(PeeOn==PETSC_FALSE) {
1242ba61063dSBarry Smith     /* check your subordinate, wait for him to be ready */
124351d315f7SKerry Stevens     ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
124451d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE) {
1245ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1246ba61063dSBarry Smith        upon return, has the lock */
124751d315f7SKerry Stevens       ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
124851d315f7SKerry Stevens     }
124951d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1250ba61063dSBarry Smith     /* your subordinate is now ready*/
125151d315f7SKerry Stevens   }
125251d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
1253ba61063dSBarry Smith   /* update your ready status */
125451d315f7SKerry Stevens   *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
125551d315f7SKerry Stevens   if(ThreadId==0) {
125651d315f7SKerry Stevens     job_chain.eJobStat = JobCompleted;
1257ba61063dSBarry Smith     /* signal main */
125851d315f7SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
125951d315f7SKerry Stevens   }
126051d315f7SKerry Stevens   else {
1261ba61063dSBarry Smith     /* tell your boss that you're ready to work */
126251d315f7SKerry Stevens     ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
126351d315f7SKerry Stevens   }
1264ba61063dSBarry Smith   /*  the while loop needs to have an exit
1265ba61063dSBarry Smith      the 'main' thread can terminate all the threads by performing a broadcast
1266ba61063dSBarry Smith    and calling FuncFinish */
126751d315f7SKerry Stevens   while(PetscThreadGo) {
1268ba61063dSBarry Smith     /* need to check the condition to ensure we don't have to wait
1269ba61063dSBarry Smith        waiting when you don't have to causes problems
1270ba61063dSBarry Smith      also need to check the condition to ensure proper handling of spurious wakeups */
127151d315f7SKerry Stevens     while(*(job_chain.arrThreadReady[ThreadId])==PETSC_TRUE) {
1272ba61063dSBarry Smith       /*upon entry, automically releases the lock and blocks
1273ba61063dSBarry Smith        upon return, has the lock */
127451d315f7SKerry Stevens         ierr = pthread_cond_wait(job_chain.cond2array[ThreadId],job_chain.mutexarray[ThreadId]);
127551d315f7SKerry Stevens 	*(job_chain.arrThreadStarted[ThreadId]) = PETSC_TRUE;
127651d315f7SKerry Stevens 	*(job_chain.arrThreadReady[ThreadId])   = PETSC_FALSE;
127751d315f7SKerry Stevens     }
127851d315f7SKerry Stevens     if(ThreadId==0) {
127951d315f7SKerry Stevens       job_chain.startJob = PETSC_FALSE;
128051d315f7SKerry Stevens       job_chain.eJobStat = ThreadsWorking;
128151d315f7SKerry Stevens     }
128251d315f7SKerry Stevens     ierr = pthread_mutex_unlock(job_chain.mutexarray[ThreadId]);
128351d315f7SKerry Stevens     if(PeeOn==PETSC_FALSE) {
1284ba61063dSBarry Smith       /* tell your subworker it's time to get to work */
128551d315f7SKerry Stevens       ierr = pthread_cond_signal(job_chain.cond2array[SubWorker]);
128651d315f7SKerry Stevens     }
1287ba61063dSBarry Smith     /* do your job */
128851d315f7SKerry Stevens     if(job_chain.pdata==NULL) {
128951d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata);
129051d315f7SKerry Stevens     }
129151d315f7SKerry Stevens     else {
129251d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_chain.pfunc(job_chain.pdata[ThreadId]);
129351d315f7SKerry Stevens     }
129451d315f7SKerry Stevens     if(iterr!=0) {
129551d315f7SKerry Stevens       ithreaderr = 1;
129651d315f7SKerry Stevens     }
129751d315f7SKerry Stevens     if(PetscThreadGo) {
1298ba61063dSBarry Smith       /* reset job, get ready for more */
129951d315f7SKerry Stevens       if(PeeOn==PETSC_FALSE) {
1300ba61063dSBarry Smith         /* check your subordinate, wait for him to be ready
1301ba61063dSBarry Smith          how do you know for a fact that your subordinate has actually started? */
130251d315f7SKerry Stevens         ierr = pthread_mutex_lock(job_chain.mutexarray[SubWorker]);
130351d315f7SKerry Stevens         while(*(job_chain.arrThreadReady[SubWorker])==PETSC_FALSE||*(job_chain.arrThreadStarted[SubWorker])==PETSC_FALSE) {
1304ba61063dSBarry Smith           /* upon entry, automically releases the lock and blocks
1305ba61063dSBarry Smith            upon return, has the lock */
130651d315f7SKerry Stevens           ierr = pthread_cond_wait(job_chain.cond1array[SubWorker],job_chain.mutexarray[SubWorker]);
130751d315f7SKerry Stevens         }
130851d315f7SKerry Stevens         ierr = pthread_mutex_unlock(job_chain.mutexarray[SubWorker]);
1309ba61063dSBarry Smith         /* your subordinate is now ready */
131051d315f7SKerry Stevens       }
131151d315f7SKerry Stevens       ierr = pthread_mutex_lock(job_chain.mutexarray[ThreadId]);
131251d315f7SKerry Stevens       *(job_chain.arrThreadReady[ThreadId]) = PETSC_TRUE;
131351d315f7SKerry Stevens       if(ThreadId==0) {
1314ba61063dSBarry Smith 	job_chain.eJobStat = JobCompleted; /* foreman: last thread to complete, guaranteed! */
1315ba61063dSBarry Smith         /* root thread (foreman) signals 'main' */
131651d315f7SKerry Stevens         ierr = pthread_cond_signal(&main_cond);
131751d315f7SKerry Stevens       }
131851d315f7SKerry Stevens       else {
1319ba61063dSBarry Smith         /* signal your boss before you go to sleep */
132051d315f7SKerry Stevens         ierr = pthread_cond_signal(job_chain.cond1array[ThreadId]);
132151d315f7SKerry Stevens       }
132251d315f7SKerry Stevens     }
132351d315f7SKerry Stevens   }
132451d315f7SKerry Stevens   return NULL;
132551d315f7SKerry Stevens }
132651d315f7SKerry Stevens 
132751d315f7SKerry Stevens #undef __FUNCT__
132851d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_Chain"
132951d315f7SKerry Stevens void* PetscThreadInitialize_Chain(PetscInt N) {
133051d315f7SKerry Stevens   PetscInt i,ierr;
133151d315f7SKerry Stevens   int status;
133251d315f7SKerry Stevens 
133351d315f7SKerry Stevens   if(PetscUseThreadPool) {
133451d315f7SKerry Stevens     size_t Val1 = (size_t)CACHE_LINE_SIZE;
133551d315f7SKerry Stevens     size_t Val2 = (size_t)PetscMaxThreads*CACHE_LINE_SIZE;
133651d315f7SKerry Stevens     arrmutex = (char*)memalign(Val1,Val2);
133751d315f7SKerry Stevens     arrcond1 = (char*)memalign(Val1,Val2);
133851d315f7SKerry Stevens     arrcond2 = (char*)memalign(Val1,Val2);
133951d315f7SKerry Stevens     arrstart = (char*)memalign(Val1,Val2);
134051d315f7SKerry Stevens     arrready = (char*)memalign(Val1,Val2);
134151d315f7SKerry Stevens     job_chain.mutexarray       = (pthread_mutex_t**)malloc(PetscMaxThreads*sizeof(pthread_mutex_t*));
134251d315f7SKerry Stevens     job_chain.cond1array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134351d315f7SKerry Stevens     job_chain.cond2array       = (pthread_cond_t**)malloc(PetscMaxThreads*sizeof(pthread_cond_t*));
134451d315f7SKerry Stevens     job_chain.arrThreadStarted = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
134551d315f7SKerry Stevens     job_chain.arrThreadReady   = (PetscBool**)malloc(PetscMaxThreads*sizeof(PetscBool*));
1346ba61063dSBarry Smith     /* initialize job structure */
134751d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
134851d315f7SKerry Stevens       job_chain.mutexarray[i]        = (pthread_mutex_t*)(arrmutex+CACHE_LINE_SIZE*i);
134951d315f7SKerry Stevens       job_chain.cond1array[i]        = (pthread_cond_t*)(arrcond1+CACHE_LINE_SIZE*i);
135051d315f7SKerry Stevens       job_chain.cond2array[i]        = (pthread_cond_t*)(arrcond2+CACHE_LINE_SIZE*i);
135151d315f7SKerry Stevens       job_chain.arrThreadStarted[i]  = (PetscBool*)(arrstart+CACHE_LINE_SIZE*i);
135251d315f7SKerry Stevens       job_chain.arrThreadReady[i]    = (PetscBool*)(arrready+CACHE_LINE_SIZE*i);
135351d315f7SKerry Stevens     }
135451d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
135551d315f7SKerry Stevens       ierr = pthread_mutex_init(job_chain.mutexarray[i],NULL);
135651d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond1array[i],NULL);
135751d315f7SKerry Stevens       ierr = pthread_cond_init(job_chain.cond2array[i],NULL);
135851d315f7SKerry Stevens       *(job_chain.arrThreadStarted[i])  = PETSC_FALSE;
135951d315f7SKerry Stevens       *(job_chain.arrThreadReady[i])    = PETSC_FALSE;
136051d315f7SKerry Stevens     }
136151d315f7SKerry Stevens     job_chain.pfunc = NULL;
136251d315f7SKerry Stevens     job_chain.pdata = (void**)malloc(N*sizeof(void*));
136351d315f7SKerry Stevens     job_chain.startJob = PETSC_FALSE;
136451d315f7SKerry Stevens     job_chain.eJobStat = JobInitiated;
136551d315f7SKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1366ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
136751d315f7SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1368ba61063dSBarry Smith     /* create threads */
136951d315f7SKerry Stevens     for(i=0; i<N; i++) {
137051d315f7SKerry Stevens       pVal[i] = i;
137151d315f7SKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1372ba61063dSBarry Smith       /* should check error */
137351d315f7SKerry Stevens     }
137451d315f7SKerry Stevens   }
137551d315f7SKerry Stevens   else {
137651d315f7SKerry Stevens   }
137751d315f7SKerry Stevens   return NULL;
137851d315f7SKerry Stevens }
137951d315f7SKerry Stevens 
138051d315f7SKerry Stevens 
138151d315f7SKerry Stevens #undef __FUNCT__
138251d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_Chain"
138351d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_Chain() {
138451d315f7SKerry Stevens   int i,ierr;
138551d315f7SKerry Stevens   void* jstatus;
138651d315f7SKerry Stevens 
138751d315f7SKerry Stevens   PetscFunctionBegin;
138851d315f7SKerry Stevens 
138951d315f7SKerry Stevens   if(PetscUseThreadPool) {
1390ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1391ba61063dSBarry Smith     /* join the threads */
139251d315f7SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
139351d315f7SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1394ba61063dSBarry Smith       /* should check error */
139551d315f7SKerry Stevens     }
139651d315f7SKerry Stevens     free(PetscThreadPoint);
139751d315f7SKerry Stevens     free(arrmutex);
139851d315f7SKerry Stevens     free(arrcond1);
139951d315f7SKerry Stevens     free(arrcond2);
140051d315f7SKerry Stevens     free(arrstart);
140151d315f7SKerry Stevens     free(arrready);
140251d315f7SKerry Stevens     free(job_chain.pdata);
140351d315f7SKerry Stevens     free(pVal);
140451d315f7SKerry Stevens   }
140551d315f7SKerry Stevens   else {
140651d315f7SKerry Stevens   }
140751d315f7SKerry Stevens   PetscFunctionReturn(0);
140851d315f7SKerry Stevens }
140951d315f7SKerry Stevens 
141051d315f7SKerry Stevens #undef __FUNCT__
141151d315f7SKerry Stevens #define __FUNCT__ "MainWait_Chain"
141251d315f7SKerry Stevens void MainWait_Chain() {
141351d315f7SKerry Stevens   int ierr;
141451d315f7SKerry Stevens   ierr = pthread_mutex_lock(job_chain.mutexarray[0]);
141551d315f7SKerry Stevens   while(job_chain.eJobStat<JobCompleted||job_chain.startJob==PETSC_TRUE) {
141651d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,job_chain.mutexarray[0]);
141751d315f7SKerry Stevens   }
141851d315f7SKerry Stevens   ierr = pthread_mutex_unlock(job_chain.mutexarray[0]);
141951d315f7SKerry Stevens }
142051d315f7SKerry Stevens 
142151d315f7SKerry Stevens #undef __FUNCT__
142251d315f7SKerry Stevens #define __FUNCT__ "MainJob_Chain"
142351d315f7SKerry Stevens PetscErrorCode MainJob_Chain(void* (*pFunc)(void*),void** data,PetscInt n) {
142451d315f7SKerry Stevens   int i,ierr;
142551d315f7SKerry Stevens   PetscErrorCode ijoberr = 0;
142636d20dc5SKerry Stevens 
142751d315f7SKerry Stevens   MainWait();
142851d315f7SKerry Stevens   job_chain.pfunc = pFunc;
142951d315f7SKerry Stevens   job_chain.pdata = data;
143051d315f7SKerry Stevens   job_chain.startJob = PETSC_TRUE;
143151d315f7SKerry Stevens   for(i=0; i<PetscMaxThreads; i++) {
143251d315f7SKerry Stevens     *(job_chain.arrThreadStarted[i]) = PETSC_FALSE;
143351d315f7SKerry Stevens   }
143451d315f7SKerry Stevens   job_chain.eJobStat = JobInitiated;
143551d315f7SKerry Stevens   ierr = pthread_cond_signal(job_chain.cond2array[0]);
143651d315f7SKerry Stevens   if(pFunc!=FuncFinish) {
1437ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done before proceeding with result collection (if any) */
143851d315f7SKerry Stevens   }
143936d20dc5SKerry Stevens 
144051d315f7SKerry Stevens   if(ithreaderr) {
144151d315f7SKerry Stevens     ijoberr = ithreaderr;
144251d315f7SKerry Stevens   }
144351d315f7SKerry Stevens   return ijoberr;
144451d315f7SKerry Stevens }
144551d315f7SKerry Stevens /****  ****/
144651d315f7SKerry Stevens 
1447ba61063dSBarry Smith #if defined(PETSC_HAVE_PTHREAD_BARRIER)
144851d315f7SKerry Stevens /**** True Thread Functions ****/
144951d315f7SKerry Stevens void* PetscThreadFunc_True(void* arg) {
145051d315f7SKerry Stevens   int icorr,ierr,iVal;
145151dcc849SKerry Stevens   int* pId = (int*)arg;
145251dcc849SKerry Stevens   int ThreadId = *pId;
14530ca81413SKerry Stevens   PetscErrorCode iterr;
145451d315f7SKerry Stevens   cpu_set_t mset;
14559e800a48SKerry Stevens   //printf("Thread %d In True Pool Thread Function\n",ThreadId);
145651d315f7SKerry Stevens   icorr = ThreadCoreAffinity[ThreadId];
145751d315f7SKerry Stevens   CPU_ZERO(&mset);
145851d315f7SKerry Stevens   CPU_SET(icorr,&mset);
145951d315f7SKerry Stevens   sched_setaffinity(0,sizeof(cpu_set_t),&mset);
146051d315f7SKerry Stevens 
146151d315f7SKerry Stevens   ierr = pthread_mutex_lock(&job_true.mutex);
146251d315f7SKerry Stevens   job_true.iNumReadyThreads++;
146351d315f7SKerry Stevens   if(job_true.iNumReadyThreads==PetscMaxThreads) {
146451dcc849SKerry Stevens     ierr = pthread_cond_signal(&main_cond);
146551dcc849SKerry Stevens   }
1466ba61063dSBarry Smith   /*the while loop needs to have an exit
1467ba61063dSBarry Smith     the 'main' thread can terminate all the threads by performing a broadcast
1468ba61063dSBarry Smith    and calling FuncFinish */
146951dcc849SKerry Stevens   while(PetscThreadGo) {
1470ba61063dSBarry Smith     /*need to check the condition to ensure we don't have to wait
1471ba61063dSBarry Smith       waiting when you don't have to causes problems
1472ba61063dSBarry Smith      also need to wait if another thread sneaks in and messes with the predicate */
147351d315f7SKerry Stevens     while(job_true.startJob==PETSC_FALSE&&job_true.iNumJobThreads==0) {
1474ba61063dSBarry Smith       /* upon entry, automically releases the lock and blocks
1475ba61063dSBarry Smith        upon return, has the lock */
147651d315f7SKerry Stevens       ierr = pthread_cond_wait(&job_true.cond,&job_true.mutex);
147751dcc849SKerry Stevens     }
147851d315f7SKerry Stevens     job_true.startJob = PETSC_FALSE;
147951d315f7SKerry Stevens     job_true.iNumJobThreads--;
148051d315f7SKerry Stevens     job_true.iNumReadyThreads--;
148151d315f7SKerry Stevens     iVal = PetscMaxThreads-job_true.iNumReadyThreads-1;
148251d315f7SKerry Stevens     pthread_mutex_unlock(&job_true.mutex);
148351d315f7SKerry Stevens     if(job_true.pdata==NULL) {
148451d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata);
148551dcc849SKerry Stevens     }
148651dcc849SKerry Stevens     else {
148751d315f7SKerry Stevens       iterr = (PetscErrorCode)(long int)job_true.pfunc(job_true.pdata[iVal]);
148851dcc849SKerry Stevens     }
14890ca81413SKerry Stevens     if(iterr!=0) {
14900ca81413SKerry Stevens       ithreaderr = 1;
14910ca81413SKerry Stevens     }
1492ba61063dSBarry Smith     /* the barrier is necessary BECAUSE: look at job_true.iNumReadyThreads
1493ba61063dSBarry Smith       what happens if a thread finishes before they all start? BAD!
1494ba61063dSBarry Smith      what happens if a thread finishes before any else start? BAD! */
1495ba61063dSBarry Smith     pthread_barrier_wait(job_true.pbarr); /* ensures all threads are finished */
1496ba61063dSBarry Smith     /* reset job */
149751dcc849SKerry Stevens     if(PetscThreadGo) {
149851d315f7SKerry Stevens       pthread_mutex_lock(&job_true.mutex);
149951d315f7SKerry Stevens       job_true.iNumReadyThreads++;
150051d315f7SKerry Stevens       if(job_true.iNumReadyThreads==PetscMaxThreads) {
1501ba61063dSBarry Smith 	/* signal the 'main' thread that the job is done! (only done once) */
150251dcc849SKerry Stevens 	ierr = pthread_cond_signal(&main_cond);
150351dcc849SKerry Stevens       }
150451dcc849SKerry Stevens     }
150551dcc849SKerry Stevens   }
150651dcc849SKerry Stevens   return NULL;
150751dcc849SKerry Stevens }
150851dcc849SKerry Stevens 
1509f09cb4aaSKerry Stevens #undef __FUNCT__
151051d315f7SKerry Stevens #define __FUNCT__ "PetscThreadInitialize_True"
151151d315f7SKerry Stevens void* PetscThreadInitialize_True(PetscInt N) {
151251dcc849SKerry Stevens   PetscInt i;
151351dcc849SKerry Stevens   int status;
15140ca81413SKerry Stevens 
15150ca81413SKerry Stevens   if(PetscUseThreadPool) {
1516f09cb4aaSKerry Stevens     pVal = (int*)malloc(N*sizeof(int));
1517ba61063dSBarry Smith     /* allocate memory in the heap for the thread structure */
151851dcc849SKerry Stevens     PetscThreadPoint = (pthread_t*)malloc(N*sizeof(pthread_t));
1519ba61063dSBarry Smith     BarrPoint = (pthread_barrier_t*)malloc((N+1)*sizeof(pthread_barrier_t)); /* BarrPoint[0] makes no sense, don't use it! */
152051d315f7SKerry Stevens     job_true.pdata = (void**)malloc(N*sizeof(void*));
152151dcc849SKerry Stevens     for(i=0; i<N; i++) {
1522f09cb4aaSKerry Stevens       pVal[i] = i;
1523f09cb4aaSKerry Stevens       status = pthread_create(&PetscThreadPoint[i],NULL,PetscThreadFunc,&pVal[i]);
1524ba61063dSBarry Smith       /* error check to ensure proper thread creation */
152551dcc849SKerry Stevens       status = pthread_barrier_init(&BarrPoint[i+1],NULL,i+1);
1526ba61063dSBarry Smith       /* should check error */
152751dcc849SKerry Stevens     }
15280ca81413SKerry Stevens   }
15290ca81413SKerry Stevens   else {
15300ca81413SKerry Stevens   }
153151dcc849SKerry Stevens   return NULL;
153251dcc849SKerry Stevens }
153351dcc849SKerry Stevens 
1534f09cb4aaSKerry Stevens 
1535f09cb4aaSKerry Stevens #undef __FUNCT__
153651d315f7SKerry Stevens #define __FUNCT__ "PetscThreadFinalize_True"
153751d315f7SKerry Stevens PetscErrorCode PetscThreadFinalize_True() {
153851dcc849SKerry Stevens   int i,ierr;
153951dcc849SKerry Stevens   void* jstatus;
154051dcc849SKerry Stevens 
154151dcc849SKerry Stevens   PetscFunctionBegin;
15420ca81413SKerry Stevens 
15430ca81413SKerry Stevens   if(PetscUseThreadPool) {
1544ba61063dSBarry Smith     MainJob(FuncFinish,NULL,PetscMaxThreads);  /* set up job and broadcast work */
1545ba61063dSBarry Smith     /* join the threads */
154651dcc849SKerry Stevens     for(i=0; i<PetscMaxThreads; i++) {
154751dcc849SKerry Stevens       ierr = pthread_join(PetscThreadPoint[i],&jstatus);
1548ba61063dSBarry Smith       /* should check error */
154951dcc849SKerry Stevens     }
155051dcc849SKerry Stevens     free(BarrPoint);
155151dcc849SKerry Stevens     free(PetscThreadPoint);
15520ca81413SKerry Stevens   }
15530ca81413SKerry Stevens   else {
15540ca81413SKerry Stevens   }
155551dcc849SKerry Stevens   PetscFunctionReturn(0);
155651dcc849SKerry Stevens }
155751dcc849SKerry Stevens 
1558f09cb4aaSKerry Stevens #undef __FUNCT__
155951d315f7SKerry Stevens #define __FUNCT__ "MainWait_True"
156051d315f7SKerry Stevens void MainWait_True() {
156151dcc849SKerry Stevens   int ierr;
156251d315f7SKerry Stevens   while(job_true.iNumReadyThreads<PetscMaxThreads||job_true.startJob==PETSC_TRUE) {
156351d315f7SKerry Stevens     ierr = pthread_cond_wait(&main_cond,&job_true.mutex);
156451dcc849SKerry Stevens   }
156551d315f7SKerry Stevens   ierr = pthread_mutex_unlock(&job_true.mutex);
156651dcc849SKerry Stevens }
156751dcc849SKerry Stevens 
1568f09cb4aaSKerry Stevens #undef __FUNCT__
156951d315f7SKerry Stevens #define __FUNCT__ "MainJob_True"
157051d315f7SKerry Stevens PetscErrorCode MainJob_True(void* (*pFunc)(void*),void** data,PetscInt n) {
157151dcc849SKerry Stevens   int ierr;
15720ca81413SKerry Stevens   PetscErrorCode ijoberr = 0;
157336d20dc5SKerry Stevens 
15740ca81413SKerry Stevens   MainWait();
157551d315f7SKerry Stevens   job_true.pfunc = pFunc;
157651d315f7SKerry Stevens   job_true.pdata = data;
157751d315f7SKerry Stevens   job_true.pbarr = &BarrPoint[n];
157851d315f7SKerry Stevens   job_true.iNumJobThreads = n;
157951d315f7SKerry Stevens   job_true.startJob = PETSC_TRUE;
158051d315f7SKerry Stevens   ierr = pthread_cond_broadcast(&job_true.cond);
15810ca81413SKerry Stevens   if(pFunc!=FuncFinish) {
1582ba61063dSBarry Smith     MainWait(); /* why wait after? guarantees that job gets done */
15830ca81413SKerry Stevens   }
158436d20dc5SKerry Stevens 
15850ca81413SKerry Stevens   if(ithreaderr) {
15860ca81413SKerry Stevens     ijoberr = ithreaderr;
15870ca81413SKerry Stevens   }
15880ca81413SKerry Stevens   return ijoberr;
158951dcc849SKerry Stevens }
1590683509dcSKerry Stevens /**** NO THREAD POOL FUNCTION ****/
1591683509dcSKerry Stevens #undef __FUNCT__
1592683509dcSKerry Stevens #define __FUNCT__ "MainJob_Spawn"
1593683509dcSKerry Stevens PetscErrorCode MainJob_Spawn(void* (*pFunc)(void*),void** data,PetscInt n) {
1594683509dcSKerry Stevens   PetscErrorCode ijoberr = 0;
1595683509dcSKerry Stevens 
1596683509dcSKerry Stevens   pthread_t* apThread = (pthread_t*)malloc(n*sizeof(pthread_t));
1597683509dcSKerry Stevens   PetscThreadRun(MPI_COMM_WORLD,pFunc,n,apThread,data);
1598683509dcSKerry Stevens   PetscThreadStop(MPI_COMM_WORLD,n,apThread); /* ensures that all threads are finished with the job */
1599683509dcSKerry Stevens   free(apThread);
1600683509dcSKerry Stevens 
1601683509dcSKerry Stevens   return ijoberr;
1602683509dcSKerry Stevens }
160351d315f7SKerry Stevens /****  ****/
1604ba61063dSBarry Smith #endif
160551dcc849SKerry Stevens 
160651dcc849SKerry Stevens void* FuncFinish(void* arg) {
160751dcc849SKerry Stevens   PetscThreadGo = PETSC_FALSE;
16080ca81413SKerry Stevens   return(0);
160951dcc849SKerry Stevens }
1610ba61063dSBarry Smith 
1611ba61063dSBarry Smith #endif
1612