xref: /petsc/src/sys/objects/pinit.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>        /*I  "petscsys.h"   I*/
6665c2dedSJed Brown #include <petscviewer.h>
7a0e72f99SJunchao Zhang 
86edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
96edef35eSSatish Balay #include <petsc/private/valgrind/valgrind.h>
106edef35eSSatish Balay #endif
116edef35eSSatish Balay 
1285649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
1385649d77SJunchao Zhang #include <petsc/private/fortranimpl.h>
1485649d77SJunchao Zhang #endif
1585649d77SJunchao Zhang 
1656883f60SBarry Smith #if defined(PETSC_USE_GCOV)
1756883f60SBarry Smith EXTERN_C_BEGIN
1856883f60SBarry Smith void  __gcov_flush(void);
1956883f60SBarry Smith EXTERN_C_END
2056883f60SBarry Smith #endif
218101f56cSMatthew Knepley 
222d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
2395c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
242d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
252d53ad75SBarry Smith #endif
262d53ad75SBarry Smith 
2727104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
2816ad0300SBarry Smith #include <petscviewersaws.h>
29a6790183SBarry Smith #endif
30eb02082bSJunchao Zhang 
31e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
32e5c89e4eSSatish Balay 
3395c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
34e5c89e4eSSatish Balay 
3595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
3795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
410069ddf5SShri Abhyankar 
426de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
43e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
456de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
466de5d289SStefano Zampini #else
476de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
486de5d289SStefano Zampini #endif
49e5c89e4eSSatish Balay 
50480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
51480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
52480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
535f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
54480cf27aSJed Brown 
55e5c89e4eSSatish Balay /*
56e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
57e5c89e4eSSatish Balay */
5802c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
5902c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
60e5c89e4eSSatish Balay 
61ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
62ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
630f8e0872SSatish Balay 
64a2f94806SJed Brown PetscInt PetscHotRegionDepth;
65a2f94806SJed Brown 
666edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
676edef35eSSatish Balay 
68b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
69b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
70b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
71b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
72b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
73b22622e2STadeu Manoel #endif
74b22622e2STadeu Manoel 
75e5c89e4eSSatish Balay /*
76945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
7772a42c3cSBarry Smith 
7872a42c3cSBarry Smith    Collective
7972a42c3cSBarry Smith 
8072a42c3cSBarry Smith    Level: advanced
8172a42c3cSBarry Smith 
8295452b02SPatrick Sanan     Notes:
83a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
840f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
85a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
860f11a792SBarry Smith 
87a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
881ea5a559SBarry Smith 
8972a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
900f11a792SBarry Smith */
91945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
9272a42c3cSBarry Smith {
9372a42c3cSBarry Smith   PetscErrorCode ierr;
9472a42c3cSBarry Smith   int            myargc   = argc;
9572a42c3cSBarry Smith   char           **myargs = args;
9672a42c3cSBarry Smith 
9772a42c3cSBarry Smith   PetscFunctionBegin;
9827104ee2SJacob Faibussowitsch   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) PetscFunctionReturn(ierr);
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPopSignalHandler());
100df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
10127104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
10272a42c3cSBarry Smith }
10372a42c3cSBarry Smith 
104f0865b08SBarry Smith /*
105a56f64adSBarry Smith       Used by Julia interface to get communicator
106f0865b08SBarry Smith */
107945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
108f0865b08SBarry Smith {
109f0865b08SBarry Smith   PetscFunctionBegin;
11027104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidPointer(comm,1);
111f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
112f0865b08SBarry Smith   PetscFunctionReturn(0);
113f0865b08SBarry Smith }
114f0865b08SBarry Smith 
115e5c89e4eSSatish Balay /*@C
116e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
117e5c89e4eSSatish Balay         the command line arguments.
118e5c89e4eSSatish Balay 
119e5c89e4eSSatish Balay    Collective
120e5c89e4eSSatish Balay 
121e5c89e4eSSatish Balay    Level: advanced
122e5c89e4eSSatish Balay 
123e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
124e5c89e4eSSatish Balay @*/
1257087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
126e5c89e4eSSatish Balay {
127e5c89e4eSSatish Balay   PetscErrorCode ierr;
128e5c89e4eSSatish Balay   int            argc   = 0;
12902c9f0b5SLisandro Dalcin   char           **args = NULL;
130e5c89e4eSSatish Balay 
131e5c89e4eSSatish Balay   PetscFunctionBegin;
1320298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
133e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
134e5c89e4eSSatish Balay }
135e5c89e4eSSatish Balay 
136e5c89e4eSSatish Balay /*@
137e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
138e5c89e4eSSatish Balay 
13993b6d2d1SJed Brown    Level: beginner
140e5c89e4eSSatish Balay 
141e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
142e5c89e4eSSatish Balay @*/
1437087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized)
144e5c89e4eSSatish Balay {
14527104ee2SJacob Faibussowitsch   PetscFunctionBegin;
14627104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized,1);
147e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
14827104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
149e5c89e4eSSatish Balay }
150e5c89e4eSSatish Balay 
151e5c89e4eSSatish Balay /*@
152e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
153e5c89e4eSSatish Balay 
154e5c89e4eSSatish Balay    Level: developer
155e5c89e4eSSatish Balay 
156e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
157e5c89e4eSSatish Balay @*/
1587087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
159e5c89e4eSSatish Balay {
16027104ee2SJacob Faibussowitsch   PetscFunctionBegin;
16127104ee2SJacob Faibussowitsch   if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized,1);
162e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
16327104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
164e5c89e4eSSatish Balay }
165e5c89e4eSSatish Balay 
16657171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []);
167e5c89e4eSSatish Balay 
168e5c89e4eSSatish Balay /*
169e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
170e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
171e5c89e4eSSatish Balay */
172367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
173e5c89e4eSSatish Balay 
174367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
175e5c89e4eSSatish Balay {
176e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
177e5c89e4eSSatish Balay 
178e5c89e4eSSatish Balay   PetscFunctionBegin;
179e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
180e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
18141e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
182e5c89e4eSSatish Balay   }
183e5c89e4eSSatish Balay 
184e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
185e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
186e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
187e5c89e4eSSatish Balay   }
188812af9f3SBarry Smith   PetscFunctionReturnVoid();
189e5c89e4eSSatish Balay }
190e5c89e4eSSatish Balay 
191e5c89e4eSSatish Balay /*
192e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
193e5c89e4eSSatish Balay sum of the second entry.
194b693b147SBarry Smith 
19576f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
196367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
197b693b147SBarry Smith there would be no place to store the both needed results.
198e5c89e4eSSatish Balay */
19976ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
200e5c89e4eSSatish Balay {
201e5c89e4eSSatish Balay   PetscFunctionBegin;
202d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
203d6e4c47cSJed Brown   {
204d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
2055f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm));
206d6e4c47cSJed Brown     *max = work.max;
207d6e4c47cSJed Brown     *sum = work.sum;
208d6e4c47cSJed Brown   }
209d6e4c47cSJed Brown #else
210d6e4c47cSJed Brown   {
211d6e4c47cSJed Brown     PetscMPIInt    size,rank;
212d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
2135f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_size(comm,&size));
2145f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(comm,&rank));
2155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(size,&work));
2165f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm));
2176ac3741eSJed Brown     *max = work[rank].max;
2186ac3741eSJed Brown     *sum = work[rank].sum;
2195f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(work));
220d6e4c47cSJed Brown   }
221d6e4c47cSJed Brown #endif
222e5c89e4eSSatish Balay   PetscFunctionReturn(0);
223e5c89e4eSSatish Balay }
224e5c89e4eSSatish Balay 
225e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
226e5c89e4eSSatish Balay 
227de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
22806a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
229e5c89e4eSSatish Balay 
230092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
231e5c89e4eSSatish Balay {
232e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
233e5c89e4eSSatish Balay 
234e5c89e4eSSatish Balay   PetscFunctionBegin;
2357c2de775SJed Brown   if (*datatype == MPIU_REAL) {
236e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
237a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2387c2de775SJed Brown   }
2397c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
240c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2417c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
242a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2437c2de775SJed Brown   }
2447c2de775SJed Brown #endif
2457c2de775SJed Brown   else {
2467c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
24741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
248e2e03761SBarry Smith   }
249812af9f3SBarry Smith   PetscFunctionReturnVoid();
250e5c89e4eSSatish Balay }
251e5c89e4eSSatish Balay #endif
252e5c89e4eSSatish Balay 
253570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
254d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
255d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
256d9822059SBarry Smith 
257092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
258d9822059SBarry Smith {
259d9822059SBarry Smith   PetscInt i,count = *cnt;
260d9822059SBarry Smith 
261d9822059SBarry Smith   PetscFunctionBegin;
2627c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2638c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
264a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2657c2de775SJed Brown   }
2667c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2677c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2687c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2697c2de775SJed Brown     for (i=0; i<count; i++) {
2707c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2717c2de775SJed Brown     }
2727c2de775SJed Brown   }
2737c2de775SJed Brown #endif
2747c2de775SJed Brown   else {
2757c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
27641e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2778c764dc5SJose Roman   }
278d9822059SBarry Smith   PetscFunctionReturnVoid();
279d9822059SBarry Smith }
280d9822059SBarry Smith 
281092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
282d9822059SBarry Smith {
283d9822059SBarry Smith   PetscInt    i,count = *cnt;
284d9822059SBarry Smith 
285d9822059SBarry Smith   PetscFunctionBegin;
2867c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2878c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
288a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2897c2de775SJed Brown   }
2907c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2917c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2927c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2937c2de775SJed Brown     for (i=0; i<count; i++) {
2947c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2957c2de775SJed Brown     }
2967c2de775SJed Brown   }
2977c2de775SJed Brown #endif
2987c2de775SJed Brown   else {
2998c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
30041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
3018c764dc5SJose Roman   }
302d9822059SBarry Smith   PetscFunctionReturnVoid();
303d9822059SBarry Smith }
304d9822059SBarry Smith #endif
305d9822059SBarry Smith 
306480cf27aSJed Brown /*
307480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
308480cf27aSJed Brown 
309ff0e51ddSBarry Smith    This is called by MPI, not by users. This is called by MPI_Comm_free() when the communicator that has this  data as an attribute is freed.
310480cf27aSJed Brown 
31112801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
312480cf27aSJed Brown 
313480cf27aSJed Brown */
31433779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
315480cf27aSJed Brown {
31605342407SJunchao Zhang   PetscCommCounter      *counter=(PetscCommCounter*)count_val;
31757f21012SBarry Smith   struct PetscCommStash *comms = counter->comms, *pcomm;
318480cf27aSJed Brown 
319480cf27aSJed Brown   PetscFunctionBegin;
3205f80ce2aSJacob Faibussowitsch   CHKERRMPI(PetscInfo(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm));
3215f80ce2aSJacob Faibussowitsch   CHKERRMPI(PetscFree(counter->iflags));
32257f21012SBarry Smith   while (comms) {
3235f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_free(&comms->comm));
32457f21012SBarry Smith     pcomm = comms;
32557f21012SBarry Smith     comms = comms->next;
3265f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(pcomm));
32757f21012SBarry Smith   }
3285f80ce2aSJacob Faibussowitsch   CHKERRMPI(PetscFree(counter));
329480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
330480cf27aSJed Brown }
331480cf27aSJed Brown 
332480cf27aSJed Brown /*
33347435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
334da3039f7SJed Brown   calls MPI_Comm_free().
335da3039f7SJed Brown 
336da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
337480cf27aSJed Brown 
338ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
339480cf27aSJed Brown 
34012801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
341480cf27aSJed Brown 
342480cf27aSJed Brown */
34333779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
344480cf27aSJed Brown {
34533779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
346480cf27aSJed Brown 
347480cf27aSJed Brown   PetscFunctionBegin;
34812801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
349ec4fadc2SJed Brown   icomm.ptr = attr_val;
35076bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
35133779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
35233779a13SJunchao Zhang     PetscMPIInt flg;
35333779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
3545f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg));
35533779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
35633779a13SJunchao Zhang     if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm's OuterComm attribute does not point to outer PETSc comm");
35733779a13SJunchao Zhang   }
3585f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval));
3595f80ce2aSJacob Faibussowitsch   CHKERRMPI(PetscInfo(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm));
360da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
361b89831e5SBarry Smith }
362da3039f7SJed Brown 
363da3039f7SJed Brown /*
36433779a13SJunchao Zhang  * This is invoked on the inner comm when Petsc_InnerComm_Attr_Delete_Fn calls MPI_Comm_delete_attr().  It should not be reached any other way.
365da3039f7SJed Brown  */
36633779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
367da3039f7SJed Brown {
368da3039f7SJed Brown   PetscFunctionBegin;
3695f80ce2aSJacob Faibussowitsch   CHKERRMPI(PetscInfo(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm));
370480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
371480cf27aSJed Brown }
372480cf27aSJed Brown 
37333779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
37442218b76SBarry Smith 
375951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3768cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3778cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3788cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
379e39fd77fSBarry Smith #endif
380e39fd77fSBarry Smith 
3810f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
38212801b39SBarry Smith 
383eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
384eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3856ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
38602c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
387dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
388e5c89e4eSSatish Balay 
389dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
390051e4cf2SJed Brown {
391051e4cf2SJed Brown   PetscFunctionBegin;
3925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferCreate(1,10000,&PetscCitationsList));
3935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCitationsRegister("@TechReport{petsc-user-ref,\n  Author = {Satish Balay and Shrirang Abhyankar and Mark F. Adams and Jed Brown \n            and Peter Brune and Kris Buschelman and Lisandro Dalcin and\n            Victor Eijkhout and William D. Gropp and Dmitry Karpeyev and\n            Dinesh Kaushik and Matthew G. Knepley and Dave A. May and Lois Curfman McInnes\n            and Richard Tran Mills and Todd Munson and Karl Rupp and Patrick Sanan\n            and Barry F. Smith and Stefano Zampini and Hong Zhang and Hong Zhang},\n  Title = {{PETS}c Users Manual},\n  Number = {ANL-95/11 - Revision 3.11},\n  Institution = {Argonne National Laboratory},\n  Year = {2019}\n}\n",NULL));
3945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCitationsRegister("@InProceedings{petsc-efficient,\n  Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n  Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n  Booktitle = {Modern Software Tools in Scientific Computing},\n  Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n  Pages = {163--202},\n  Publisher = {Birkh{\\\"{a}}user Press},\n  Year = {1997}\n}\n",NULL));
395051e4cf2SJed Brown   PetscFunctionReturn(0);
396051e4cf2SJed Brown }
397e5c89e4eSSatish Balay 
3982d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3992d747510SLisandro Dalcin 
4002d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4012d747510SLisandro Dalcin {
4022d747510SLisandro Dalcin   PetscFunctionBegin;
4035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(programname,name,sizeof(programname)));
4042d747510SLisandro Dalcin   PetscFunctionReturn(0);
4052d747510SLisandro Dalcin }
4062d747510SLisandro Dalcin 
4072d747510SLisandro Dalcin /*@C
4082d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4092d747510SLisandro Dalcin 
4102d747510SLisandro Dalcin     Not Collective
4112d747510SLisandro Dalcin 
4122d747510SLisandro Dalcin     Input Parameter:
4132d747510SLisandro Dalcin .   len - length of the string name
4142d747510SLisandro Dalcin 
4152d747510SLisandro Dalcin     Output Parameter:
4162d747510SLisandro Dalcin .   name - the name of the running program
4172d747510SLisandro Dalcin 
4182d747510SLisandro Dalcin    Level: advanced
4192d747510SLisandro Dalcin 
4202d747510SLisandro Dalcin     Notes:
4212d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4222d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4232d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4242d747510SLisandro Dalcin @*/
4252d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4262d747510SLisandro Dalcin {
4272d747510SLisandro Dalcin   PetscFunctionBegin;
4285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(name,programname,len));
4292d747510SLisandro Dalcin   PetscFunctionReturn(0);
4302d747510SLisandro Dalcin }
4312d747510SLisandro Dalcin 
432e5c89e4eSSatish Balay /*@C
433e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
434e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
435e5c89e4eSSatish Balay 
436e5c89e4eSSatish Balay    Not Collective
437e5c89e4eSSatish Balay 
438e5c89e4eSSatish Balay    Output Parameters:
439e5c89e4eSSatish Balay +  argc - count of number of command line arguments
440e5c89e4eSSatish Balay -  args - the command line arguments
441e5c89e4eSSatish Balay 
442e5c89e4eSSatish Balay    Level: intermediate
443e5c89e4eSSatish Balay 
444e5c89e4eSSatish Balay    Notes:
445e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
446e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
447e5c89e4eSSatish Balay 
448f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
449f177e3b1SBarry Smith 
450793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
451e5c89e4eSSatish Balay 
452e5c89e4eSSatish Balay @*/
4537087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
454e5c89e4eSSatish Balay {
455e5c89e4eSSatish Balay   PetscFunctionBegin;
4562c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!PetscInitializeCalled && PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
457e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
458e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
459e5c89e4eSSatish Balay   PetscFunctionReturn(0);
460e5c89e4eSSatish Balay }
461e5c89e4eSSatish Balay 
462793721a6SBarry Smith /*@C
463793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
464793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
465793721a6SBarry Smith 
466793721a6SBarry Smith    Not Collective
467793721a6SBarry Smith 
468793721a6SBarry Smith    Output Parameters:
469793721a6SBarry Smith .  args - the command line arguments
470793721a6SBarry Smith 
471793721a6SBarry Smith    Level: intermediate
472793721a6SBarry Smith 
473793721a6SBarry Smith    Notes:
474793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
475793721a6SBarry Smith 
476793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
477793721a6SBarry Smith 
478793721a6SBarry Smith @*/
4797087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
480793721a6SBarry Smith {
481793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
482793721a6SBarry Smith 
483793721a6SBarry Smith   PetscFunctionBegin;
4842c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!PetscInitializeCalled && PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4852d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
4865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(argc,args));
487793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
4885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]));
489793721a6SBarry Smith   }
4902d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
491793721a6SBarry Smith   PetscFunctionReturn(0);
492793721a6SBarry Smith }
493793721a6SBarry Smith 
494793721a6SBarry Smith /*@C
495793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
496793721a6SBarry Smith 
497793721a6SBarry Smith    Not Collective
498793721a6SBarry Smith 
499793721a6SBarry Smith    Output Parameters:
500793721a6SBarry Smith .  args - the command line arguments
501793721a6SBarry Smith 
502793721a6SBarry Smith    Level: intermediate
503793721a6SBarry Smith 
504793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
505793721a6SBarry Smith 
506793721a6SBarry Smith @*/
5077087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
508793721a6SBarry Smith {
509793721a6SBarry Smith   PetscInt       i = 0;
510793721a6SBarry Smith 
511793721a6SBarry Smith   PetscFunctionBegin;
512a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
513793721a6SBarry Smith   while (args[i]) {
5145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(args[i]));
515793721a6SBarry Smith     i++;
516793721a6SBarry Smith   }
5175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(args));
518793721a6SBarry Smith   PetscFunctionReturn(0);
519793721a6SBarry Smith }
520793721a6SBarry Smith 
52127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
52230befbd2SBarry Smith #include <petscconfiginfo.h>
52330befbd2SBarry Smith 
52495c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
52511525c0dSBarry Smith {
52627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
52711525c0dSBarry Smith   if (!PetscGlobalRank) {
52830befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
52911525c0dSBarry Smith     int            port;
530ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
53111525c0dSBarry Smith     size_t         applinelen,introlen;
53211525c0dSBarry Smith     PetscErrorCode ierr;
533ffbd1cfbSBarry Smith     char           sawsurl[256];
53411525c0dSBarry Smith 
5355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsHasName(NULL,NULL,"-saws_log",&flg));
53611525c0dSBarry Smith     if (flg) {
53711525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
53811525c0dSBarry Smith 
5395f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL));
54011525c0dSBarry Smith       if (sawslog[0]) {
54111525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
54211525c0dSBarry Smith       } else {
54311525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
54411525c0dSBarry Smith       }
54511525c0dSBarry Smith     }
5465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg));
54711525c0dSBarry Smith     if (flg) {
54811525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
54911525c0dSBarry Smith     }
5505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL));
551ffbd1cfbSBarry Smith     if (selectport) {
552ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
553ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
554ffbd1cfbSBarry Smith     } else {
5555f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg));
55611525c0dSBarry Smith       if (flg) {
55711525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
55811525c0dSBarry Smith       }
559ffbd1cfbSBarry Smith     }
5605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg));
56111525c0dSBarry Smith     if (flg) {
5622f613bf5SBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
5635f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcmp(root,".",&rootlocal));
5649c1e0ce8SBarry Smith     } else {
5655f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsHasName(NULL,NULL,"-saws_options",&flg));
5669c1e0ce8SBarry Smith       if (flg) {
5675f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root)));
5682f613bf5SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
5699c1e0ce8SBarry Smith       }
57011525c0dSBarry Smith     }
5715f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2));
57211525c0dSBarry Smith     if (flg2) {
57311525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
574*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
5755f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root));
5765f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscTestDirectory(jsdir,'r',&flg));
577*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
5782f613bf5SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());
57911525c0dSBarry Smith     }
5805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscGetProgramName(programname,sizeof(programname)));
5815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrlen(help,&applinelen));
58211525c0dSBarry Smith     introlen   = 4096 + applinelen;
58330a8c9c0SSurtai Han     applinelen += 1024;
5845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc(applinelen,&appline));
5855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc(introlen,&intro));
58611525c0dSBarry Smith 
58711525c0dSBarry Smith     if (rootlocal) {
5885f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSNPrintf(appline,applinelen,"%s.c.html",programname));
5895f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscTestFile(appline,'r',&rootlocal));
59011525c0dSBarry Smith     }
5915f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetAll(NULL,&options));
59211525c0dSBarry Smith     if (rootlocal && help) {
5935f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSNPrintf(appline,applinelen,"<center> Running <a href=\"%s.c.html\">%s</a> %s</center><br><center><pre>%s</pre></center><br>\n",programname,programname,options,help));
59411525c0dSBarry Smith     } else if (help) {
5955f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help));
59611525c0dSBarry Smith     } else {
5975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options));
59811525c0dSBarry Smith     }
5995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(options));
6005f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscGetVersion(version,sizeof(version)));
6015f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSNPrintf(intro,introlen,"<body>\n"
602a17b96a8SKyle Gerard Felker                           "<center><h2> <a href=\"https://petsc.org/\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n"
603df62c222SBarry Smith                           "<center>This is the default PETSc application dashboard, from it you can access any published PETSc objects or logging data</center><br><center>%s configured with %s</center><br>\n"
6045f80ce2aSJacob Faibussowitsch                           "%s",version,petscconfigureoptions,appline));
60543da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
6065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(intro));
6075f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(appline));
608ffbd1cfbSBarry Smith     if (selectport) {
609aa573868SBarry Smith       PetscBool silent;
6107d812c46SBarry Smith 
6117d812c46SBarry Smith       ierr = SAWs_Initialize();
6127d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6137d812c46SBarry Smith       while (ierr) {
6147d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6157d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6167d812c46SBarry Smith         ierr = SAWs_Initialize();
6177d812c46SBarry Smith       }
6187d812c46SBarry Smith 
6195f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL));
620aa573868SBarry Smith       if (!silent) {
621ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
6225f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl));
623ffbd1cfbSBarry Smith       }
6247d812c46SBarry Smith     } else {
6257d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
626aa573868SBarry Smith     }
6275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCitationsRegister("@TechReport{ saws,\n"
6280af79b04SBarry Smith                                    "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6290af79b04SBarry Smith                                    "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6300af79b04SBarry Smith                                    "  Institution = {Argonne National Laboratory},\n"
6315f80ce2aSJacob Faibussowitsch                                    "  Year   = 2013\n}\n",NULL));
63211525c0dSBarry Smith   }
63311525c0dSBarry Smith   PetscFunctionReturn(0);
63411525c0dSBarry Smith }
63511525c0dSBarry Smith #endif
63611525c0dSBarry Smith 
6374dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6384dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6394dfee713SSatish Balay {
6404dfee713SSatish Balay   PetscFunctionBegin;
6414dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6424dfee713SSatish Balay     /* see MPI.py for details on this bug */
6434dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6444dfee713SSatish Balay #endif
6454dfee713SSatish Balay   PetscFunctionReturn(0);
6464dfee713SSatish Balay }
6474dfee713SSatish Balay 
648a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
649a56f64adSBarry Smith #include <adios.h>
65022580e64SBarry Smith #include <adios_read.h>
6517b56e58cSSatish Balay int64_t Petsc_adios_group;
652a56f64adSBarry Smith #endif
653a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
654cd1b99a6SBarry Smith #include <omp.h>
655f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
656cd1b99a6SBarry Smith #endif
657a56f64adSBarry Smith 
658a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_DEVICE)
659a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
660a4af0ceeSJacob Faibussowitsch #  if PetscDefined(HAVE_CUDA)
661a4af0ceeSJacob Faibussowitsch // REMOVE ME
662a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
663a4af0ceeSJacob Faibussowitsch #  endif
664a4af0ceeSJacob Faibussowitsch #  if PetscDefined(HAVE_HIP)
665a4af0ceeSJacob Faibussowitsch // REMOVE ME
666a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
667a4af0ceeSJacob Faibussowitsch #  endif
668a4af0ceeSJacob Faibussowitsch #endif
669a4af0ceeSJacob Faibussowitsch 
67027104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
671bc8a24ecSBarry Smith #include <dlfcn.h>
672bc8a24ecSBarry Smith #endif
673a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
674a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
675a4af0ceeSJacob Faibussowitsch #endif
676a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
677a4af0ceeSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViennaCLInit();
678a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE;
679a4af0ceeSJacob Faibussowitsch #endif
680bc8a24ecSBarry Smith 
68185649d77SJunchao Zhang /*
68285649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
68385649d77SJunchao Zhang 
68485649d77SJunchao Zhang   prog:     program name
68502101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
68685649d77SJunchao Zhang   help:     program help message
68702101c96SSatish Balay   ftn:      is it called from Fortran initilization (petscinitializef_)?
68885649d77SJunchao Zhang   readarguments,len: used when fortran is true
68985649d77SJunchao Zhang */
69002101c96SSatish Balay PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char* prog,const char* file,const char *help,PetscBool ftn,PetscBool readarguments,PetscInt len)
69185649d77SJunchao Zhang {
69285649d77SJunchao Zhang   PetscMPIInt size;
69385649d77SJunchao Zhang   PetscBool   flg = PETSC_TRUE;
69485649d77SJunchao Zhang   char        hostname[256];
69585649d77SJunchao Zhang 
69627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
69727104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
69885649d77SJunchao Zhang   /*
69985649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
70085649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
70185649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
70285649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
70385649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
70485649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
70585649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
70685649d77SJunchao Zhang       listed above and since that time are compatible.
70785649d77SJunchao Zhang 
70885649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
70985649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
71085649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
71185649d77SJunchao Zhang       to perform the checking.
71285649d77SJunchao Zhang 
71385649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
71485649d77SJunchao Zhang 
71585649d77SJunchao Zhang       Questions:
71685649d77SJunchao Zhang 
71785649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
71885649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
71985649d77SJunchao Zhang   */
72085649d77SJunchao Zhang 
72185649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
72285649d77SJunchao Zhang   {
72385649d77SJunchao Zhang     char           mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
72485649d77SJunchao Zhang     PetscMPIInt    mpilibraryversionlength;
7255f80ce2aSJacob Faibussowitsch     PetscErrorCode ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);
72627104ee2SJacob Faibussowitsch     if (ierr) PetscFunctionReturn(ierr);
72785649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
72885649d77SJunchao Zhang #if defined(MPICH_VERSION)
72985649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000
73085649d77SJunchao Zhang     {
73185649d77SJunchao Zhang       char      *ver,*lf;
73285649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
73327104ee2SJacob Faibussowitsch       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);
73427104ee2SJacob Faibussowitsch       if (ierr) PetscFunctionReturn(ierr);
73527104ee2SJacob Faibussowitsch       else if (ver) {
73627104ee2SJacob Faibussowitsch         ierr = PetscStrchr(ver,'\n',&lf);
73727104ee2SJacob Faibussowitsch         if (ierr) PetscFunctionReturn(ierr);
73827104ee2SJacob Faibussowitsch         else if (lf) {
73985649d77SJunchao Zhang           *lf = 0;
74027104ee2SJacob Faibussowitsch           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) PetscFunctionReturn(ierr);
74185649d77SJunchao Zhang         }
74285649d77SJunchao Zhang       }
74385649d77SJunchao Zhang       if (!flg) {
7447d3de750SJacob Faibussowitsch         PetscInfo(NULL,"PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n",mpilibraryversion,MPICH_VESION);
74585649d77SJunchao Zhang         flg = PETSC_TRUE;
74685649d77SJunchao Zhang       }
74785649d77SJunchao Zhang     }
74885649d77SJunchao Zhang #endif
74985649d77SJunchao Zhang     /* check for OpenMPI version, it is not part of the MPI ABI initiative (is it part of another initiative that needs to be handled?) */
75085649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION)
75185649d77SJunchao Zhang     {
75285649d77SJunchao Zhang       char *ver,bs[MPI_MAX_LIBRARY_VERSION_STRING],*bsf;
75385649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
75485649d77SJunchao Zhang #define PSTRSZ 2
75585649d77SJunchao Zhang       char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI","FUJITSU MPI"};
75685649d77SJunchao Zhang       char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v","Library "};
75785649d77SJunchao Zhang       int i;
75885649d77SJunchao Zhang       for (i=0; i<PSTRSZ; i++) {
75927104ee2SJacob Faibussowitsch         ierr = PetscStrstr(mpilibraryversion,ompistr1[i],&ver);
76027104ee2SJacob Faibussowitsch         if (ierr) PetscFunctionReturn(ierr);
76127104ee2SJacob Faibussowitsch         else if (ver) {
76285649d77SJunchao Zhang           PetscSNPrintf(bs,MPI_MAX_LIBRARY_VERSION_STRING,"%s%d.%d",ompistr2[i],OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
76327104ee2SJacob Faibussowitsch           ierr = PetscStrstr(ver,bs,&bsf);
76427104ee2SJacob Faibussowitsch           if (ierr) PetscFunctionReturn(ierr);
76527104ee2SJacob Faibussowitsch           else if (bsf) flg = PETSC_TRUE;
76685649d77SJunchao Zhang           break;
76785649d77SJunchao Zhang         }
76885649d77SJunchao Zhang       }
76985649d77SJunchao Zhang       if (!flg) {
7707d3de750SJacob Faibussowitsch         PetscInfo(NULL,"PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n",mpilibraryversion,OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
77185649d77SJunchao Zhang         flg = PETSC_TRUE;
77285649d77SJunchao Zhang       }
77385649d77SJunchao Zhang     }
77485649d77SJunchao Zhang #endif
77585649d77SJunchao Zhang   }
77685649d77SJunchao Zhang #endif
77785649d77SJunchao Zhang 
77885649d77SJunchao Zhang #if defined(PETSC_HAVE_DLSYM)
77985649d77SJunchao Zhang   /* These symbols are currently in the OpenMPI and MPICH libraries; they may not always be, in that case the test will simply not detect the problem */
78027104ee2SJacob Faibussowitsch   if (PetscUnlikely(dlsym(RTLD_DEFAULT,"ompi_mpi_init") && dlsym(RTLD_DEFAULT,"MPID_Abort"))) {
78185649d77SJunchao Zhang     fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
7825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStackView(stderr));
78327104ee2SJacob Faibussowitsch     PetscFunctionReturn(PETSC_ERR_MPI_LIB_INCOMP);
78485649d77SJunchao Zhang   }
78585649d77SJunchao Zhang #endif
78685649d77SJunchao Zhang 
78785649d77SJunchao Zhang   /* these must be initialized in a routine, not as a constant declaration*/
78885649d77SJunchao Zhang   PETSC_STDOUT = stdout;
78985649d77SJunchao Zhang   PETSC_STDERR = stderr;
79085649d77SJunchao Zhang 
79185649d77SJunchao Zhang   /*CHKERRQ can be used from now */
79285649d77SJunchao Zhang   PetscErrorHandlingInitialized = PETSC_TRUE;
79385649d77SJunchao Zhang 
79485649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
79585649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
79685649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
79785649d77SJunchao Zhang #endif
79885649d77SJunchao Zhang 
7995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsCreateDefault());
80085649d77SJunchao Zhang 
80185649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
80285649d77SJunchao Zhang 
8035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSetProgramName(prog));
8045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen));
8055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout));
8065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr));
8075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockCreate(&PetscCommSpinLock));
80885649d77SJunchao Zhang 
80985649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
8105f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN));
81185649d77SJunchao Zhang 
81285649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
8135f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS));
8145f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE));
81585649d77SJunchao Zhang   }
81685649d77SJunchao Zhang 
81785649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
8185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscErrorPrintfInitialize());
81985649d77SJunchao Zhang 
8205f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank));
8215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize));
82285649d77SJunchao Zhang 
82385649d77SJunchao Zhang   MPIU_BOOL = MPI_INT;
82485649d77SJunchao Zhang   MPIU_ENUM = MPI_INT;
82585649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
82685649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
82785649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
82885649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
82985649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
83085649d77SJunchao Zhang #endif
83127104ee2SJacob Faibussowitsch   else {
83227104ee2SJacob Faibussowitsch     (*PetscErrorPrintf)("PetscInitialize_Common: Could not find MPI type for size_t\n");
83327104ee2SJacob Faibussowitsch     PetscFunctionReturn(PETSC_ERR_SUP_SYS);
83427104ee2SJacob Faibussowitsch   }
83585649d77SJunchao Zhang 
83685649d77SJunchao Zhang   /*
83785649d77SJunchao Zhang      Initialized the global complex variable; this is because with
83885649d77SJunchao Zhang      shared libraries the constructors for global variables
83985649d77SJunchao Zhang      are not called; at least on IRIX.
84085649d77SJunchao Zhang   */
84185649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
84285649d77SJunchao Zhang   {
84385649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
84485649d77SJunchao Zhang     PetscComplex ic(0.0,1.0);
84585649d77SJunchao Zhang     PETSC_i = ic;
84685649d77SJunchao Zhang #else
84785649d77SJunchao Zhang     PETSC_i = _Complex_I;
84885649d77SJunchao Zhang #endif
84985649d77SJunchao Zhang   }
85085649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
85185649d77SJunchao Zhang 
85285649d77SJunchao Zhang   /*
85385649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
85485649d77SJunchao Zhang      half of the entries and maxes the second half.
85585649d77SJunchao Zhang   */
8565f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP));
85785649d77SJunchao Zhang 
85885649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128)
8595f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128));
8605f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&MPIU___FLOAT128));
86185649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
8625f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128));
8635f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&MPIU___COMPLEX128));
86485649d77SJunchao Zhang #endif
8655f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_create(PetscMax_Local,1,&MPIU_MAX));
8665f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_create(PetscMin_Local,1,&MPIU_MIN));
86785649d77SJunchao Zhang #elif defined(PETSC_USE_REAL___FP16)
8685f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16));
8695f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&MPIU___FP16));
8705f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_create(PetscMax_Local,1,&MPIU_MAX));
8715f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_create(PetscMin_Local,1,&MPIU_MIN));
87285649d77SJunchao Zhang #endif
87385649d77SJunchao Zhang 
87485649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
8755f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_create(PetscSum_Local,1,&MPIU_SUM));
87685649d77SJunchao Zhang #endif
87785649d77SJunchao Zhang 
8785f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR));
8795f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&MPIU_2SCALAR));
88085649d77SJunchao Zhang 
88185649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
88285649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
88385649d77SJunchao Zhang   {
88485649d77SJunchao Zhang     struct PetscRealInt { PetscReal v; PetscInt i; };
88585649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
88685649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscRealInt,v),offsetof(struct PetscRealInt,i)};
88785649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_REAL,MPIU_INT}, tmpStruct;
88885649d77SJunchao Zhang 
8895f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct));
8905f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscRealInt),&MPIU_REAL_INT));
8915f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_free(&tmpStruct));
8925f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_commit(&MPIU_REAL_INT));
89385649d77SJunchao Zhang   }
89485649d77SJunchao Zhang   {
89585649d77SJunchao Zhang     struct PetscScalarInt { PetscScalar v; PetscInt i; };
89685649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
89785649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscScalarInt,v),offsetof(struct PetscScalarInt,i)};
89885649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_SCALAR,MPIU_INT}, tmpStruct;
89985649d77SJunchao Zhang 
9005f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct));
9015f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscScalarInt),&MPIU_SCALAR_INT));
9025f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_free(&tmpStruct));
9035f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
90485649d77SJunchao Zhang   }
90585649d77SJunchao Zhang #endif
90685649d77SJunchao Zhang 
90785649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9085f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT));
9095f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&MPIU_2INT));
91085649d77SJunchao Zhang #endif
9115f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(4,MPI_INT,&MPI_4INT));
9125f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&MPI_4INT));
9135f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(4,MPIU_INT,&MPIU_4INT));
9145f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&MPIU_4INT));
91585649d77SJunchao Zhang 
91685649d77SJunchao Zhang   /*
91785649d77SJunchao Zhang      Attributes to be set on PETSc communicators
91885649d77SJunchao Zhang   */
9195f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0));
9205f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0));
9215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0));
9225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0));
92385649d77SJunchao Zhang 
92485649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
9255f80ce2aSJacob Faibussowitsch   if (ftn) CHKERRQ(PetscInitFortran_Private(readarguments,file,len));
92685649d77SJunchao Zhang   else
92785649d77SJunchao Zhang #endif
9285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInsert(NULL,&PetscGlobalArgc,&PetscGlobalArgs,file));
92985649d77SJunchao Zhang 
93085649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscErrorPrintfInitialize());
93285649d77SJunchao Zhang 
93385649d77SJunchao Zhang   /*
93485649d77SJunchao Zhang      Check system options and print help
93585649d77SJunchao Zhang   */
9365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsCheckInitial_Private(help));
93785649d77SJunchao Zhang 
938a4af0ceeSJacob Faibussowitsch   /*
939a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
940a4af0ceeSJacob Faibussowitsch 
941a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
942a4af0ceeSJacob Faibussowitsch    1. MPI initialized
943a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
944a4af0ceeSJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up its own SIGSEV handler via
945a4af0ceeSJacob Faibussowitsch       the push/pop interface.
946a4af0ceeSJacob Faibussowitsch   */
947a2158755SJunchao Zhang #if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP) || PetscDefined(HAVE_SYCL))
9485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
949a4af0ceeSJacob Faibussowitsch #endif
950a4af0ceeSJacob Faibussowitsch 
951a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
952a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
9535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-log_summary",&flg));
9545f80ce2aSJacob Faibussowitsch   if (!flg) CHKERRQ(PetscOptionsHasName(NULL,NULL,"-log_view",&flg));
9555f80ce2aSJacob Faibussowitsch   if (!flg) CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-viennacl_synchronize",&flg,NULL));
956a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
9575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViennaCLInit());
958a4af0ceeSJacob Faibussowitsch #endif
959a4af0ceeSJacob Faibussowitsch 
960a4af0ceeSJacob Faibussowitsch   /*
961a4af0ceeSJacob Faibussowitsch      Creates the logging data structures; this is enabled even if logging is not turned on
962a4af0ceeSJacob Faibussowitsch      This is the last thing we do before returning to the user code to prevent having the
963a4af0ceeSJacob Faibussowitsch      logging numbers contaminated by any startup time associated with MPI
964a4af0ceeSJacob Faibussowitsch   */
965a4af0ceeSJacob Faibussowitsch #if defined(PETSC_USE_LOG)
9665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogInitialize());
967a4af0ceeSJacob Faibussowitsch #endif
968a4af0ceeSJacob Faibussowitsch 
9695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCitationsInitialize());
97085649d77SJunchao Zhang 
97185649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
9725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInitializeSAWs(ftn ? NULL : help));
97327104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
9745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-stack_view",&flg));
9755f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PetscStackViewSAWs());
97685649d77SJunchao Zhang #endif
97785649d77SJunchao Zhang 
97885649d77SJunchao Zhang   /*
97985649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
98085649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
98185649d77SJunchao Zhang   */
9825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInitialize_DynamicLibraries());
98385649d77SJunchao Zhang 
9845f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
9855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL,"PETSc successfully started: number of processors = %d\n",size));
9865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGetHostName(hostname,256));
9875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL,"Running on machine: %s\n",hostname));
98885649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
98985649d77SJunchao Zhang   {
99085649d77SJunchao Zhang     PetscBool       omp_view_flag;
99185649d77SJunchao Zhang     char           *threads = getenv("OMP_NUM_THREADS");
9925f80ce2aSJacob Faibussowitsch     PetscErrorCode  ierr;
99385649d77SJunchao Zhang 
99485649d77SJunchao Zhang     if (threads) {
9955f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(NULL,"Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n",threads));
99685649d77SJunchao Zhang       (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
99785649d77SJunchao Zhang     } else {
9982f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt) omp_get_max_threads();
9995f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(NULL,"Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n",PetscNumOMPThreads));
100085649d77SJunchao Zhang     }
100185649d77SJunchao Zhang     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
10025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-omp_num_threads","Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS","None",PetscNumOMPThreads,&PetscNumOMPThreads,&flg));
10035f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag));
100485649d77SJunchao Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
100585649d77SJunchao Zhang     if (flg) {
10065f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(NULL,"Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n",PetscNumOMPThreads));
100785649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
100885649d77SJunchao Zhang     }
100985649d77SJunchao Zhang     if (omp_view_flag) {
10105f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %" PetscInt_FMT "\n",PetscNumOMPThreads));
101185649d77SJunchao Zhang     }
101285649d77SJunchao Zhang   }
101385649d77SJunchao Zhang #endif
101485649d77SJunchao Zhang 
101585649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
101685649d77SJunchao Zhang   /*
101785649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
101885649d77SJunchao Zhang 
101985649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
102085649d77SJunchao Zhang   */
10215f80ce2aSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) CHKERRMPI(MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL));
102285649d77SJunchao Zhang #endif
102385649d77SJunchao Zhang 
102485649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFPTCreate(10000));
102685649d77SJunchao Zhang #endif
102785649d77SJunchao Zhang 
102885649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
102985649d77SJunchao Zhang   {
103085649d77SJunchao Zhang     PetscViewer viewer;
10315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg));
103285649d77SJunchao Zhang     if (flg) {
10335f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscProcessPlacementView(viewer));
10345f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&viewer));
103585649d77SJunchao Zhang     }
103685649d77SJunchao Zhang   }
103785649d77SJunchao Zhang #endif
103885649d77SJunchao Zhang 
103985649d77SJunchao Zhang   flg  = PETSC_TRUE;
10405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL));
10415f80ce2aSJacob Faibussowitsch   if (!flg) CHKERRQ(PetscOptionsPushGetViewerOff(PETSC_TRUE));
104285649d77SJunchao Zhang 
104385649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10445f80ce2aSJacob Faibussowitsch   CHKERRQ(adios_init_noxml(PETSC_COMM_WORLD));
10455f80ce2aSJacob Faibussowitsch   CHKERRQ(adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default));
10465f80ce2aSJacob Faibussowitsch   CHKERRQ(adios_select_method(Petsc_adios_group,"MPI","",""));
10475f80ce2aSJacob Faibussowitsch   CHKERRQ(adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,""));
104885649d77SJunchao Zhang #endif
104985649d77SJunchao Zhang 
105085649d77SJunchao Zhang #if defined(__VALGRIND_H)
105185649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND? PETSC_TRUE: PETSC_FALSE;
105285649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
10535f80ce2aSJacob Faibussowitsch   if (PETSC_RUNNING_ON_VALGRIND) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"WARNING: Running valgrind with the MacOS native BLAS and LAPACK can fail. If it fails suggest configuring with --download-fblaslapack or --download-f2cblaslapack"));
105485649d77SJunchao Zhang #endif
105585649d77SJunchao Zhang #endif
105685649d77SJunchao Zhang   /*
105785649d77SJunchao Zhang       Set flag that we are completely initialized
105885649d77SJunchao Zhang   */
105985649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
106085649d77SJunchao Zhang 
10615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-python",&flg));
10625f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PetscPythonInitialize(NULL,NULL));
106327104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
106485649d77SJunchao Zhang }
106585649d77SJunchao Zhang 
1066e5c89e4eSSatish Balay /*@C
1067e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1068e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
1069e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1070e5c89e4eSSatish Balay    your program -- usually the very first line!
1071e5c89e4eSSatish Balay 
1072e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
1073e5c89e4eSSatish Balay 
1074e5c89e4eSSatish Balay    Input Parameters:
1075e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1076e5c89e4eSSatish Balay .  args - the command line arguments
1077be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1078be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1079be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1080c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
10810298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1082e5c89e4eSSatish Balay 
108305827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
108405827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
108505827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
108605827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
108705827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1088e5c89e4eSSatish Balay 
1089e5c89e4eSSatish Balay    Options Database Keys:
10907ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
10917ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1092e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
10938a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
10948a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
10958a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
10968a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
10978a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1098bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1099e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1100e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1101e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
110279dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
110379dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
110479dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
1105aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
110692f119d6SBarry Smith .  -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds, ignored in optimized build. May want to set in PETSC_OPTIONS environmental variable
110792f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
110892f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1109608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
111092f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1111e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1112e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1113e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1114e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1115e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
11160841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1117e5c89e4eSSatish Balay 
1118c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1119c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1120c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1121c5b5d8d5SVaclav Hapla -  -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet()
1122c5b5d8d5SVaclav Hapla 
1123c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1124c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1125c5b5d8d5SVaclav Hapla    They can be used also in option files.
1126c5b5d8d5SVaclav Hapla 
1127c5b5d8d5SVaclav Hapla    See PetscOptionsMonitorSet() to do monitoring programmatically.
1128c5b5d8d5SVaclav Hapla 
1129e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1130a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1131fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
1132217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1133217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1134495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1135e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
11369a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
113779dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
11389a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1139495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1140f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1141495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
1142495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
1143c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
114487c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
1145c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1146495fc317SBarry Smith 
1147609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
1148e5c89e4eSSatish Balay 
1149ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1150ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1151ffbd1cfbSBarry Smith .  -saws_port_auto_select - have SAWs select a new unique port number where it publishes the data, the URL is printed to the screen
1152ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1153ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1154ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1155ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1156ffbd1cfbSBarry Smith 
1157e5c89e4eSSatish Balay    Environmental Variables:
1158e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
1159e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
1160e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
11614a971ea4SToby Isaac .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
11624a971ea4SToby Isaac .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1163e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
1164e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
1165e5c89e4eSSatish Balay 
1166e5c89e4eSSatish Balay    Level: beginner
1167e5c89e4eSSatish Balay 
1168e5c89e4eSSatish Balay    Notes:
1169e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
1170e5c89e4eSSatish Balay    it before PetscInitialize().
1171e5c89e4eSSatish Balay 
1172e5c89e4eSSatish Balay    Fortran Version:
1173e5c89e4eSSatish Balay    In Fortran this routine has the format
1174e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
1175e5c89e4eSSatish Balay 
1176e5c89e4eSSatish Balay +  ierr - error return code
1177c5b5d8d5SVaclav Hapla -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
1178c5b5d8d5SVaclav Hapla           Use PETSC_NULL_CHARACTER to not check for code specific file.
1179c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
1180e5c89e4eSSatish Balay 
1181e5c89e4eSSatish Balay    Important Fortran Note:
11820eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
11830298fd71SBarry Smith    null character string; you CANNOT just use NULL as
1184a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
1185e5c89e4eSSatish Balay 
118601cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
118701cb0274SBarry Smith    calling PetscInitialize().
1188e5c89e4eSSatish Balay 
118901cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
1190e5c89e4eSSatish Balay 
1191e5c89e4eSSatish Balay @*/
11927087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
1193e5c89e4eSSatish Balay {
119485649d77SJunchao Zhang   PetscMPIInt    flag;
119585649d77SJunchao Zhang   const char     *prog = "Unknown Name";
1196e5c89e4eSSatish Balay 
119727104ee2SJacob Faibussowitsch   PetscFunctionBegin;
119827104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
11995f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Initialized(&flag));
1200e5c89e4eSSatish Balay   if (!flag) {
12012c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PETSC_COMM_WORLD != MPI_COMM_NULL,PETSC_COMM_SELF,PETSC_ERR_SUP,"You cannot set PETSC_COMM_WORLD if you have not initialized MPI first");
12025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPreMPIInit_Private());
12035e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12045e765c61SJed Brown     {
12055e765c61SJed Brown       PetscMPIInt provided;
12065f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided));
12075e765c61SJed Brown     }
12085e765c61SJed Brown #else
12095f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Init(argc,args));
12105e765c61SJed Brown #endif
1211e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1212e5c89e4eSSatish Balay   }
12134dfee713SSatish Balay 
121485649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1215e5c89e4eSSatish Balay   if (argc && args) {
1216e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1217e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1218e5c89e4eSSatish Balay   }
12195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInitialize_Common(prog,file,help,PETSC_FALSE/*C*/,PETSC_FALSE,0));
122027104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
1221e5c89e4eSSatish Balay }
1222e5c89e4eSSatish Balay 
1223a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
122495c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1225ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1226ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
122795c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
12284097062eSBarry Smith #endif
1229e5c89e4eSSatish Balay 
1230008a6e76SBarry Smith /*
1231008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1232008a6e76SBarry Smith */
1233008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1234008a6e76SBarry Smith {
1235008a6e76SBarry Smith   PetscFunctionBegin;
1236008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
12375f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU___FLOAT128));
1238008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
12395f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU___COMPLEX128));
1240008a6e76SBarry Smith #endif
12415f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_free(&MPIU_MAX));
12425f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_free(&MPIU_MIN));
1243008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
12445f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU___FP16));
12455f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_free(&MPIU_MAX));
12465f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_free(&MPIU_MIN));
1247008a6e76SBarry Smith #endif
1248008a6e76SBarry Smith 
1249de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
12505f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_free(&MPIU_SUM));
1251008a6e76SBarry Smith #endif
1252008a6e76SBarry Smith 
12535f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU_2SCALAR));
12545f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU_REAL_INT));
12555f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU_SCALAR_INT));
125640df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
12575f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU_2INT));
1258008a6e76SBarry Smith #endif
12595f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPI_4INT));
12605f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&MPIU_4INT));
12615f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
1262008a6e76SBarry Smith   PetscFunctionReturn(0);
1263008a6e76SBarry Smith }
1264008a6e76SBarry Smith 
1265a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
1266a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1267a4af0ceeSJacob Faibussowitsch #endif
1268a4af0ceeSJacob Faibussowitsch 
1269e5c89e4eSSatish Balay /*@C
1270e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1271e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1272e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1273e5c89e4eSSatish Balay 
1274e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1275e5c89e4eSSatish Balay 
1276e5c89e4eSSatish Balay    Options Database Keys:
127726a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1278e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
12797eb1d149SBarry Smith .  -objects_dump [all] - Prints list of objects allocated by the user that have not been freed, the option all cause all outstanding objects to be listed
1280e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
128192f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1282e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
128392f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1284e5c89e4eSSatish Balay 
1285e5c89e4eSSatish Balay    Level: beginner
1286e5c89e4eSSatish Balay 
1287e5c89e4eSSatish Balay    Note:
1288e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1289e5c89e4eSSatish Balay 
129088c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1291e5c89e4eSSatish Balay @*/
12927087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1293e5c89e4eSSatish Balay {
12944bb5149bSJed Brown   PetscMPIInt    rank;
1295a8d2bbe5SBarry Smith   PetscInt       nopt;
12962bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1297dff31646SBarry Smith   PetscBool      flg;
129810463e74SBarry Smith #if defined(PETSC_USE_LOG)
129910463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
130010463e74SBarry Smith #endif
1301e5c89e4eSSatish Balay 
13023cbbc5ffSBarry Smith   PetscFunctionBegin;
130327104ee2SJacob Faibussowitsch   if (PetscUnlikely(!PetscInitializeCalled)) {
13041724198aSStefano Zampini     fprintf(PETSC_STDOUT,"PetscInitialize() must be called before PetscFinalize()\n");
13055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStackView(PETSC_STDOUT));
13061724198aSStefano Zampini     PetscStackClearTop;
130727104ee2SJacob Faibussowitsch     return PETSC_ERR_ARG_WRONGSTATE;
130827104ee2SJacob Faibussowitsch   }
13095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL,"PetscFinalize() called\n"));
1310b022a5c1SBarry Smith 
13115f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1312a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13135f80ce2aSJacob Faibussowitsch   CHKERRQ(adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE));
13145f80ce2aSJacob Faibussowitsch   CHKERRQ(adios_finalize(rank));
1315a56f64adSBarry Smith #endif
13165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-citations",&flg));
1317dff31646SBarry Smith   if (flg) {
13181f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
13191f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
13201f817a21SBarry Smith 
13215f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL));
13221f817a21SBarry Smith     if (filename[0]) {
13235f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd));
13241f817a21SBarry Smith     }
13255f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSegBufferGet(PetscCitationsList,1,&cits));
1326dff31646SBarry Smith     cits[0] = 0;
13275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSegBufferExtractAlloc(PetscCitationsList,&cits));
13285f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n"));
13295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n"));
13305f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits));
13315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n"));
13325f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFClose(PETSC_COMM_WORLD,fd));
13335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(cits));
1334dff31646SBarry Smith   }
13355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferDestroy(&PetscCitationsList));
1336dff31646SBarry Smith 
1337c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
133804102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
133904102261SBarry Smith   {
134004102261SBarry Smith     PetscInt nmax = 2;
134104102261SBarry Smith     char     **buffs;
13425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(2,&buffs));
13435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1));
134404102261SBarry Smith     if (flg1) {
1345*28b400f6SJacob Faibussowitsch       PetscCheck(nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
134604102261SBarry Smith       if (nmax == 1) {
13475f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(128,&buffs[1]));
13485f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscGetProgramName(buffs[1],32));
13495f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscStrcat(buffs[1]," has completed"));
135004102261SBarry Smith       }
13515f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL));
13525f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(buffs[0]));
13535f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(buffs[1]));
135404102261SBarry Smith     }
13555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(buffs));
135604102261SBarry Smith   }
1357763ec7b1SBarry Smith   {
1358763ec7b1SBarry Smith     PetscInt nmax = 2;
1359763ec7b1SBarry Smith     char     **buffs;
13605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(2,&buffs));
13615f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1));
1362763ec7b1SBarry Smith     if (flg1) {
1363*28b400f6SJacob Faibussowitsch       PetscCheck(nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1364763ec7b1SBarry Smith       if (nmax == 1) {
13655f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(128,&buffs[1]));
13665f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscGetProgramName(buffs[1],32));
13675f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscStrcat(buffs[1]," has completed"));
1368763ec7b1SBarry Smith       }
13695f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL));
13705f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(buffs[0]));
13715f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(buffs[1]));
1372763ec7b1SBarry Smith     }
13735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(buffs));
1374763ec7b1SBarry Smith   }
137504102261SBarry Smith #endif
137604102261SBarry Smith 
13772d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
13785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFPTDestroy());
13792d53ad75SBarry Smith #endif
13802d53ad75SBarry Smith 
1381e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1382dff31646SBarry Smith   flg = PETSC_FALSE;
13835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL));
1384d5649816SBarry Smith   if (flg) {
13855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsSAWsDestroy());
1386d5649816SBarry Smith   }
1387d5649816SBarry Smith #endif
1388d5649816SBarry Smith 
1389681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1390681455b2SBarry Smith   flg1 = PETSC_FALSE;
13915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL));
1392681455b2SBarry Smith   if (flg1) {
1393681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
13945f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL));
1395681455b2SBarry Smith   }
1396681455b2SBarry Smith #endif
1397681455b2SBarry Smith 
139867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
13995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL));
1400e5c89e4eSSatish Balay   if (!flg2) {
140190d69ab7SBarry Smith     flg2 = PETSC_FALSE;
14025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL));
1403e5c89e4eSSatish Balay   }
1404e5c89e4eSSatish Balay   if (flg2) {
14055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n"));
1406e5c89e4eSSatish Balay   }
140767584ceeSBarry Smith #endif
1408e5c89e4eSSatish Balay 
1409e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
141090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL));
1412e5c89e4eSSatish Balay   if (flg1) {
1413e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
14145f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD));
14155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops));
1416e5c89e4eSSatish Balay   }
1417e5c89e4eSSatish Balay #endif
1418e5c89e4eSSatish Balay 
1419e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1420e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1421e5c89e4eSSatish Balay   mname[0] = 0;
14225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1));
1423e5c89e4eSSatish Balay   if (flg1) {
14245f80ce2aSJacob Faibussowitsch     if (mname[0]) CHKERRQ(PetscLogMPEDump(mname));
14255f80ce2aSJacob Faibussowitsch     else          CHKERRQ(PetscLogMPEDump(0));
1426e5c89e4eSSatish Balay   }
1427e5c89e4eSSatish Balay #endif
1428356e5837SBarry Smith #endif
1429a297a907SKarl Rupp 
1430dd710f27SBarry Smith   /*
1431dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1432dd710f27SBarry Smith   */
14335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectRegisterDestroyAll());
1434dd710f27SBarry Smith 
1435356e5837SBarry Smith #if defined(PETSC_USE_LOG)
14365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsPushGetViewerOff(PETSC_FALSE));
14375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogViewFromOptions());
14385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsPopGetViewerOff());
143987c3beb6SLisandro Dalcin 
1440356e5837SBarry Smith   mname[0] = 0;
14415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1));
1442e5c89e4eSSatish Balay   if (flg1) {
144391eabc43SBarry Smith     PetscViewer viewer;
14445f80ce2aSJacob Faibussowitsch     CHKERRQ((*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n"));
144591eabc43SBarry Smith     if (mname[0]) {
14465f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer));
14475f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogView(viewer));
14485f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&viewer));
144933f85c2fSBarry Smith     } else {
145033f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
14515f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT));
14525f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogView(viewer));
14535f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPopFormat(viewer));
145433f85c2fSBarry Smith     }
1455e5c89e4eSSatish Balay   }
1456a297a907SKarl Rupp 
1457dd710f27SBarry Smith   /*
1458dd710f27SBarry Smith      Free any objects created by the last block of code.
1459dd710f27SBarry Smith   */
14605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectRegisterDestroyAll());
1461dd710f27SBarry Smith 
1462dd710f27SBarry Smith   mname[0] = 0;
14635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1));
14645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2));
14655f80ce2aSJacob Faibussowitsch   if (flg1 || flg2) CHKERRQ(PetscLogDump(mname));
1466e5c89e4eSSatish Balay #endif
146710463e74SBarry Smith 
146890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL));
14705f80ce2aSJacob Faibussowitsch   if (!flg1) CHKERRQ(PetscPopSignalHandler());
147190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL));
1473e5c89e4eSSatish Balay   if (flg1) {
14745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMPIDump(stdout));
1475e5c89e4eSSatish Balay   }
147690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
147790d69ab7SBarry Smith   flg2 = PETSC_FALSE;
14788bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
14795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1));
14805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1));
14815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL));
1482e4c476e2SSatish Balay 
1483e5c89e4eSSatish Balay   if (flg2) {
1484be56827dSJed Brown     PetscViewer viewer;
14855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
14865f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerSetType(viewer,PETSCVIEWERASCII));
14875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsView(NULL,viewer));
14885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
1489e5c89e4eSSatish Balay   }
1490e5c89e4eSSatish Balay 
1491e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
14925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-nox",&flg1));
14935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1));
1494e5c89e4eSSatish Balay 
149533fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
14965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1));
14979245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1498e5c89e4eSSatish Balay   if (flg3) {
14993de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1500be56827dSJed Brown       PetscViewer viewer;
15015f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
15025f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerSetType(viewer,PETSCVIEWERASCII));
15035f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsView(NULL,viewer));
15045f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&viewer));
1505e5c89e4eSSatish Balay     }
15065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsAllUsed(NULL,&nopt));
15073de2bfdfSBarry Smith     if (nopt) {
15085f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n"));
15095f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n"));
15103de2bfdfSBarry Smith       if (nopt == 1) {
15115f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n"));
1512e5c89e4eSSatish Balay       } else {
15135f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"There are %" PetscInt_FMT " unused database options. They are:\n",nopt));
1514e5c89e4eSSatish Balay       }
15153de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15165f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n"));
1517df12ba86SBarry Smith     }
15185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsLeft(NULL));
1519e5c89e4eSSatish Balay   }
1520e5c89e4eSSatish Balay 
1521e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1522d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15235f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStackSAWsViewOff());
152416ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1525d45a07a7SBarry Smith   }
1526ec957eceSBarry Smith #endif
1527ec957eceSBarry Smith 
15284097062eSBarry Smith #if defined(PETSC_USE_LOG)
152910463e74SBarry Smith   /*
1530dbc8283eSBarry Smith        List all objects the user may have forgot to free
15312eff7a51SBarry Smith   */
153205df10baSBarry Smith   if (PetscObjectsLog) {
15335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1));
1534a64a8e02SBarry Smith     if (flg1) {
1535a64a8e02SBarry Smith       MPI_Comm local_comm;
15367eb1d149SBarry Smith       char     string[64];
1537a64a8e02SBarry Smith 
15385f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL));
15395f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_dup(MPI_COMM_WORLD,&local_comm));
15405f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSequentialPhaseBegin_Private(local_comm,1));
15415f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
15425f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSequentialPhaseEnd_Private(local_comm,1));
15435f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_free(&local_comm));
15440a1571b3SBarry Smith     }
154505df10baSBarry Smith   }
15464097062eSBarry Smith #endif
15474097062eSBarry Smith 
15484097062eSBarry Smith #if defined(PETSC_USE_LOG)
1549dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1550dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
15515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(PetscObjects));
15524097062eSBarry Smith #endif
15532eff7a51SBarry Smith 
155433f85c2fSBarry Smith   /*
155533f85c2fSBarry Smith      Destroy any packages that registered a finalize
155633f85c2fSBarry Smith   */
15575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRegisterFinalizeAll());
155833f85c2fSBarry Smith 
1559101409b8SToby Isaac #if defined(PETSC_USE_LOG)
15605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFinalize());
1561101409b8SToby Isaac #endif
1562101409b8SToby Isaac 
156333f85c2fSBarry Smith   /*
156448dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
156548dd1dffSBarry Smith 
15665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListPrintAll());
156748dd1dffSBarry Smith   */
156837e93019SBarry Smith 
15694028d114SSatish Balay   if (petsc_history) {
15705f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCloseHistoryFile(&petsc_history));
157102c9f0b5SLisandro Dalcin     petsc_history = NULL;
1572e5c89e4eSSatish Balay   }
15735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
15745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfoDestroy());
1575e5c89e4eSSatish Balay 
157667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
157792f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1578e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
157992f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1580e5c89e4eSSatish Balay     FILE *fd;
1581ed9cf6e9SBarry Smith     int  err;
1582e5c89e4eSSatish Balay 
1583dc92acbaSJed Brown     flg2 = PETSC_FALSE;
158492f119d6SBarry Smith     flg3 = PETSC_FALSE;
15855f80ce2aSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL));
15865f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL));
158792f119d6SBarry Smith     fname[0] = 0;
15885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1));
1589e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1590e5c89e4eSSatish Balay 
1591589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1592*28b400f6SJacob Faibussowitsch       fd   = fopen(sname,"w"); PetscCheck(fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
15935f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMallocDump(fd));
1594ed9cf6e9SBarry Smith       err  = fclose(fd);
1595*28b400f6SJacob Faibussowitsch       PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
159692f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1597e5c89e4eSSatish Balay       MPI_Comm local_comm;
1598e5c89e4eSSatish Balay 
15995f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_dup(MPI_COMM_WORLD,&local_comm));
16005f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSequentialPhaseBegin_Private(local_comm,1));
16015f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMallocDump(stdout));
16025f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSequentialPhaseEnd_Private(local_comm,1));
16035f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_free(&local_comm));
1604e5c89e4eSSatish Balay     }
1605e5c89e4eSSatish Balay     fname[0] = 0;
16065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1));
1607e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1608e5c89e4eSSatish Balay 
1609589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1610*28b400f6SJacob Faibussowitsch       fd   = fopen(sname,"w"); PetscCheck(fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
16115f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMallocView(fd));
1612ed9cf6e9SBarry Smith       err  = fclose(fd);
1613*28b400f6SJacob Faibussowitsch       PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
161492f119d6SBarry Smith     } else if (flg1) {
161592f119d6SBarry Smith       MPI_Comm local_comm;
161692f119d6SBarry Smith 
16175f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_dup(MPI_COMM_WORLD,&local_comm));
16185f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSequentialPhaseBegin_Private(local_comm,1));
16195f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMallocView(stdout));
16205f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSequentialPhaseEnd_Private(local_comm,1));
16215f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_free(&local_comm));
1622e5c89e4eSSatish Balay     }
1623e5c89e4eSSatish Balay   }
162467584ceeSBarry Smith #endif
162520e2c332SMatthew G. Knepley 
16265486ca60SMatthew G. Knepley   /*
16275486ca60SMatthew G. Knepley      Close any open dynamic libraries
16285486ca60SMatthew G. Knepley   */
16295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFinalize_DynamicLibraries());
16305486ca60SMatthew G. Knepley 
1631e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsDestroyDefault());
1633e5c89e4eSSatish Balay 
1634e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
163502c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1636e5c89e4eSSatish Balay 
1637c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1638c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
16395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscKokkosFinalize_Private());
1640c2b86a48SJunchao Zhang     PetscBeganKokkos = PETSC_FALSE;
164145639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1642c2b86a48SJunchao Zhang   }
1643c2b86a48SJunchao Zhang #endif
1644c2b86a48SJunchao Zhang 
164571438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
164671438e86SJunchao Zhang   if (PetscBeganNvshmem) {
16475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscNvshmemFinalize());
164871438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
164971438e86SJunchao Zhang   }
165071438e86SJunchao Zhang #endif
1651a0e72f99SJunchao Zhang 
16525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFreeMPIResources());
1653e5c89e4eSSatish Balay 
1654dbc8283eSBarry Smith   /*
1655efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1656efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1657efb80d3cSBarry Smith 
1658efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1659efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1660dbc8283eSBarry Smith  */
1661b770b1f6SSatish Balay   {
1662dbc8283eSBarry Smith     PetscCommCounter *counter;
1663dbc8283eSBarry Smith     PetscMPIInt      flg;
1664dbc8283eSBarry Smith     MPI_Comm         icomm;
1665265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
16665f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg));
1667dbc8283eSBarry Smith     if (flg) {
1668265f3f35SJed Brown       icomm = ucomm.comm;
16695f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg));
1670*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1671dbc8283eSBarry Smith 
16725f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval));
16735f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval));
16745f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_free(&icomm));
1675dbc8283eSBarry Smith     }
16765f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg));
1677dbc8283eSBarry Smith     if (flg) {
1678265f3f35SJed Brown       icomm = ucomm.comm;
16795f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg));
1680*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1681dbc8283eSBarry Smith 
16825f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval));
16835f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval));
16845f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_free(&icomm));
1685dbc8283eSBarry Smith     }
1686b770b1f6SSatish Balay   }
1687dbc8283eSBarry Smith 
16885f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
16895f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
16905f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
16915f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
1692480cf27aSJed Brown 
16935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
16945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
16955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
16965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpinlockDestroy(&PetscCommSpinLock));
1697ef19f930SBarry Smith 
1698e5c89e4eSSatish Balay   if (PetscBeganMPI) {
169999b1327fSBarry Smith     PetscMPIInt flag;
17005f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Finalized(&flag));
17012c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
17025f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Finalize());
1703e5c89e4eSSatish Balay   }
1704e5c89e4eSSatish Balay /*
1705e5c89e4eSSatish Balay 
1706e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1707e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1708e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1709e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1710e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17110e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1712e5c89e4eSSatish Balay    memory was not freed.
1713e5c89e4eSSatish Balay 
1714e5c89e4eSSatish Balay */
17155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMallocClear());
17165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStackReset());
1717a297a907SKarl Rupp 
17188ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_FALSE;
1719e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1720e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
172156883f60SBarry Smith #if defined(PETSC_USE_GCOV)
172256883f60SBarry Smith   /*
172356883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
172456883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
172556883f60SBarry Smith    */
172656883f60SBarry Smith   __gcov_flush();
172756883f60SBarry Smith #endif
17281724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
17291724198aSStefano Zampini   PetscStackClearTop;
173027104ee2SJacob Faibussowitsch   return 0;
1731e5c89e4eSSatish Balay }
1732e5c89e4eSSatish Balay 
173343db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
17348cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
173543db4dbbSBarry Smith {
173643db4dbbSBarry Smith   if (*a == *b) return 1;
173743db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
173843db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
173943db4dbbSBarry Smith   return 0;
174043db4dbbSBarry Smith }
1741a70650f6SBarry Smith #endif
174243db4dbbSBarry Smith 
174343db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
17448cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
174543db4dbbSBarry Smith {
174643db4dbbSBarry Smith   if (*a == *b) return 1;
174743db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
174843db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
174943db4dbbSBarry Smith   return 0;
175043db4dbbSBarry Smith }
175143db4dbbSBarry Smith #endif
1752