xref: /petsc/src/sys/objects/pinit.c (revision 613bf2b20e6d15d382f10220b3a3e6e72f647221)
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
18aaf3846cSSatish Balay #if defined(PETSC_HAVE___GCOV_DUMP)
19aaf3846cSSatish Balay #define __gcov_flush(x) __gcov_dump(x)
20aaf3846cSSatish Balay #endif
2156883f60SBarry Smith void  __gcov_flush(void);
2256883f60SBarry Smith EXTERN_C_END
2356883f60SBarry Smith #endif
248101f56cSMatthew Knepley 
252d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
2695c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
272d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
282d53ad75SBarry Smith #endif
292d53ad75SBarry Smith 
3027104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
3116ad0300SBarry Smith #include <petscviewersaws.h>
32a6790183SBarry Smith #endif
33eb02082bSJunchao Zhang 
34e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
35e5c89e4eSSatish Balay 
3695c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
37e5c89e4eSSatish Balay 
3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
4295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
430069ddf5SShri Abhyankar 
446de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
45e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
476de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
486de5d289SStefano Zampini #else
496de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
506de5d289SStefano Zampini #endif
51e5c89e4eSSatish Balay 
52480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
53480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
54480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
555f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
56480cf27aSJed Brown 
57e5c89e4eSSatish Balay /*
58e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
59e5c89e4eSSatish Balay */
6002c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
6102c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
62e5c89e4eSSatish Balay 
63ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
64ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
650f8e0872SSatish Balay 
66a2f94806SJed Brown PetscInt PetscHotRegionDepth;
67a2f94806SJed Brown 
686edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
696edef35eSSatish Balay 
70b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
71b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
72b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
73b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
74b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
75b22622e2STadeu Manoel #endif
76b22622e2STadeu Manoel 
77e5c89e4eSSatish Balay /*
78945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
7972a42c3cSBarry Smith 
8072a42c3cSBarry Smith    Collective
8172a42c3cSBarry Smith 
8272a42c3cSBarry Smith    Level: advanced
8372a42c3cSBarry Smith 
8495452b02SPatrick Sanan     Notes:
85a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
860f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
87a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
880f11a792SBarry Smith 
89a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
901ea5a559SBarry Smith 
91db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()`
920f11a792SBarry Smith */
93945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
9472a42c3cSBarry Smith {
9572a42c3cSBarry Smith   int    myargc = argc;
9672a42c3cSBarry Smith   char **myargs = args;
9772a42c3cSBarry Smith 
9872a42c3cSBarry Smith   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&myargc,&myargs,filename,help));
1009566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
101df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
10227104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
10372a42c3cSBarry Smith }
10472a42c3cSBarry Smith 
105f0865b08SBarry Smith /*
106a56f64adSBarry Smith       Used by Julia interface to get communicator
107f0865b08SBarry Smith */
108945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
109f0865b08SBarry Smith {
110f0865b08SBarry Smith   PetscFunctionBegin;
11127104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidPointer(comm,1);
112f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
113f0865b08SBarry Smith   PetscFunctionReturn(0);
114f0865b08SBarry Smith }
115f0865b08SBarry Smith 
116e5c89e4eSSatish Balay /*@C
117e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
118e5c89e4eSSatish Balay         the command line arguments.
119e5c89e4eSSatish Balay 
120e5c89e4eSSatish Balay    Collective
121e5c89e4eSSatish Balay 
122e5c89e4eSSatish Balay    Level: advanced
123e5c89e4eSSatish Balay 
124db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`
125e5c89e4eSSatish Balay @*/
1267087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
127e5c89e4eSSatish Balay {
128e5c89e4eSSatish Balay   int    argc = 0;
12902c9f0b5SLisandro Dalcin   char **args = NULL;
130e5c89e4eSSatish Balay 
131e5c89e4eSSatish Balay   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,NULL,NULL));
13339a651e2SJacob Faibussowitsch   PetscFunctionReturn(0);
134e5c89e4eSSatish Balay }
135e5c89e4eSSatish Balay 
136e5c89e4eSSatish Balay /*@
137e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
138e5c89e4eSSatish Balay 
13993b6d2d1SJed Brown    Level: beginner
140e5c89e4eSSatish Balay 
141db781477SPatrick Sanan .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 
156db781477SPatrick Sanan .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;
2059566063dSJacob Faibussowitsch     PetscCallMPI(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;
2139566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm,&size));
2149566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm,&rank));
2159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&work));
2161c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm));
2176ac3741eSJed Brown     *max = work[rank].max;
2186ac3741eSJed Brown     *sum = work[rank].sum;
2199566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
220d6e4c47cSJed Brown   }
221d6e4c47cSJed Brown #endif
222e5c89e4eSSatish Balay   PetscFunctionReturn(0);
223e5c89e4eSSatish Balay }
224e5c89e4eSSatish Balay 
225e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
226e5c89e4eSSatish Balay 
227*613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
228*613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
229*613bf2b2SPierre Jolivet #include <quadmath.h>
230*613bf2b2SPierre Jolivet MPI_Op MPIU_SUM___FLOAT128 = 0;
231*613bf2b2SPierre Jolivet #endif
232de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
23306a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
234*613bf2b2SPierre Jolivet #endif
235e5c89e4eSSatish Balay 
236092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
237e5c89e4eSSatish Balay {
238e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
239e5c89e4eSSatish Balay 
240e5c89e4eSSatish Balay   PetscFunctionBegin;
2417c2de775SJed Brown   if (*datatype == MPIU_REAL) {
242e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
243a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2447c2de775SJed Brown   }
2457c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
246c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2477c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
248a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2497c2de775SJed Brown   }
2507c2de775SJed Brown #endif
251*613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
252*613bf2b2SPierre Jolivet   else if (*datatype == MPIU___FLOAT128) {
253*613bf2b2SPierre Jolivet     __float128 *xin = (__float128*)in,*xout = (__float128*)out;
254*613bf2b2SPierre Jolivet     for (i=0; i<count; i++) xout[i] += xin[i];
255*613bf2b2SPierre Jolivet   }
256*613bf2b2SPierre Jolivet   else if (*datatype == MPIU___COMPLEX128) {
257*613bf2b2SPierre Jolivet     __complex128 *xin = (__complex128*)in,*xout = (__complex128*)out;
258*613bf2b2SPierre Jolivet     for (i=0; i<count; i++) xout[i] += xin[i];
259*613bf2b2SPierre Jolivet   }
260*613bf2b2SPierre Jolivet #endif
2617c2de775SJed Brown   else {
262*613bf2b2SPierre Jolivet #if !defined(PETSC_HAVE_REAL___FLOAT128)
2637c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
264*613bf2b2SPierre Jolivet #else
265*613bf2b2SPierre Jolivet     (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types");
266*613bf2b2SPierre Jolivet #endif
26741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
268e2e03761SBarry Smith   }
269812af9f3SBarry Smith   PetscFunctionReturnVoid();
270e5c89e4eSSatish Balay }
271e5c89e4eSSatish Balay #endif
272e5c89e4eSSatish Balay 
273570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
274d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
275d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
276d9822059SBarry Smith 
277092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
278d9822059SBarry Smith {
279d9822059SBarry Smith   PetscInt i,count = *cnt;
280d9822059SBarry Smith 
281d9822059SBarry Smith   PetscFunctionBegin;
2827c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2838c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
284a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2857c2de775SJed Brown   }
2867c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2877c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2887c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2897c2de775SJed Brown     for (i=0; i<count; i++) {
2907c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2917c2de775SJed Brown     }
2927c2de775SJed Brown   }
2937c2de775SJed Brown #endif
2947c2de775SJed Brown   else {
2957c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
29641e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2978c764dc5SJose Roman   }
298d9822059SBarry Smith   PetscFunctionReturnVoid();
299d9822059SBarry Smith }
300d9822059SBarry Smith 
301092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
302d9822059SBarry Smith {
303d9822059SBarry Smith   PetscInt    i,count = *cnt;
304d9822059SBarry Smith 
305d9822059SBarry Smith   PetscFunctionBegin;
3067c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3078c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
308a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
3097c2de775SJed Brown   }
3107c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
3117c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3127c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
3137c2de775SJed Brown     for (i=0; i<count; i++) {
3147c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3157c2de775SJed Brown     }
3167c2de775SJed Brown   }
3177c2de775SJed Brown #endif
3187c2de775SJed Brown   else {
3198c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
32041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
3218c764dc5SJose Roman   }
322d9822059SBarry Smith   PetscFunctionReturnVoid();
323d9822059SBarry Smith }
324d9822059SBarry Smith #endif
325d9822059SBarry Smith 
326480cf27aSJed Brown /*
327480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
328480cf27aSJed Brown 
329ff0e51ddSBarry 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.
330480cf27aSJed Brown 
33112801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
332480cf27aSJed Brown 
333480cf27aSJed Brown */
33433779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
335480cf27aSJed Brown {
33605342407SJunchao Zhang   PetscCommCounter      *counter=(PetscCommCounter*)count_val;
33757f21012SBarry Smith   struct PetscCommStash *comms = counter->comms, *pcomm;
338480cf27aSJed Brown 
339480cf27aSJed Brown   PetscFunctionBegin;
3409566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm));
3419566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter->iflags));
34257f21012SBarry Smith   while (comms) {
3439566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&comms->comm));
34457f21012SBarry Smith     pcomm = comms;
34557f21012SBarry Smith     comms = comms->next;
3469566063dSJacob Faibussowitsch     PetscCall(PetscFree(pcomm));
34757f21012SBarry Smith   }
3489566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter));
349480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
350480cf27aSJed Brown }
351480cf27aSJed Brown 
352480cf27aSJed Brown /*
35347435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
354da3039f7SJed Brown   calls MPI_Comm_free().
355da3039f7SJed Brown 
356da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
357480cf27aSJed Brown 
358ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
359480cf27aSJed Brown 
36012801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
361480cf27aSJed Brown 
362480cf27aSJed Brown */
36333779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
364480cf27aSJed Brown {
36533779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
366480cf27aSJed Brown 
367480cf27aSJed Brown   PetscFunctionBegin;
36812801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
369ec4fadc2SJed Brown   icomm.ptr = attr_val;
37076bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
37133779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
37233779a13SJunchao Zhang     PetscMPIInt flg;
37333779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
3749566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg));
37533779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
37633779a13SJunchao 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");
37733779a13SJunchao Zhang   }
3789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval));
3799566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm));
380da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
381b89831e5SBarry Smith }
382da3039f7SJed Brown 
383da3039f7SJed Brown /*
38433779a13SJunchao 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.
385da3039f7SJed Brown  */
38633779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
387da3039f7SJed Brown {
388da3039f7SJed Brown   PetscFunctionBegin;
3899566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm));
390480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
391480cf27aSJed Brown }
392480cf27aSJed Brown 
39333779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
39442218b76SBarry Smith 
395951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3968cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3978cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3988cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
399e39fd77fSBarry Smith #endif
400e39fd77fSBarry Smith 
4010f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
40212801b39SBarry Smith 
403eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
404eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4056ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
40602c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
407dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
408e5c89e4eSSatish Balay 
409dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
410051e4cf2SJed Brown {
411051e4cf2SJed Brown   PetscFunctionBegin;
4129566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(1,10000,&PetscCitationsList));
41330c35bf2SSatish Balay 
41430c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\
41530c35bf2SSatish Balay   Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\
41630c35bf2SSatish Balay     and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\
41730c35bf2SSatish Balay     and Victor Eijkhout and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\
41830c35bf2SSatish Balay     and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\
41930c35bf2SSatish Balay     and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\
42030c35bf2SSatish Balay     and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\
42130c35bf2SSatish Balay     and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\
42230c35bf2SSatish Balay   Title = {{PETSc/TAO} Users Manual},\n\
42330c35bf2SSatish Balay   Number = {ANL-21/39 - Revision 3.17},\n\
42430c35bf2SSatish Balay   Institution = {Argonne National Laboratory},\n\
42530c35bf2SSatish Balay   Year = {2022}\n}\n",NULL));
42630c35bf2SSatish Balay 
42730c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\
42830c35bf2SSatish Balay   Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\
42930c35bf2SSatish Balay   Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\
43030c35bf2SSatish Balay   Booktitle = {Modern Software Tools in Scientific Computing},\n\
43130c35bf2SSatish Balay   Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\
43230c35bf2SSatish Balay   Pages = {163--202},\n\
43330c35bf2SSatish Balay   Publisher = {Birkh{\\\"{a}}user Press},\n\
43430c35bf2SSatish Balay   Year = {1997}\n}\n",NULL));
43530c35bf2SSatish Balay 
436051e4cf2SJed Brown   PetscFunctionReturn(0);
437051e4cf2SJed Brown }
438e5c89e4eSSatish Balay 
4392d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4402d747510SLisandro Dalcin 
4412d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4422d747510SLisandro Dalcin {
4432d747510SLisandro Dalcin   PetscFunctionBegin;
4449566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(programname,name,sizeof(programname)));
4452d747510SLisandro Dalcin   PetscFunctionReturn(0);
4462d747510SLisandro Dalcin }
4472d747510SLisandro Dalcin 
4482d747510SLisandro Dalcin /*@C
4492d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4502d747510SLisandro Dalcin 
4512d747510SLisandro Dalcin     Not Collective
4522d747510SLisandro Dalcin 
4532d747510SLisandro Dalcin     Input Parameter:
4542d747510SLisandro Dalcin .   len - length of the string name
4552d747510SLisandro Dalcin 
4562d747510SLisandro Dalcin     Output Parameter:
4572d747510SLisandro Dalcin .   name - the name of the running program
4582d747510SLisandro Dalcin 
4592d747510SLisandro Dalcin    Level: advanced
4602d747510SLisandro Dalcin 
4612d747510SLisandro Dalcin     Notes:
4622d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4632d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4642d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4652d747510SLisandro Dalcin @*/
4662d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4672d747510SLisandro Dalcin {
4682d747510SLisandro Dalcin   PetscFunctionBegin;
4699566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(name,programname,len));
4702d747510SLisandro Dalcin   PetscFunctionReturn(0);
4712d747510SLisandro Dalcin }
4722d747510SLisandro Dalcin 
473e5c89e4eSSatish Balay /*@C
474e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
475e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
476e5c89e4eSSatish Balay 
477e5c89e4eSSatish Balay    Not Collective
478e5c89e4eSSatish Balay 
479e5c89e4eSSatish Balay    Output Parameters:
480e5c89e4eSSatish Balay +  argc - count of number of command line arguments
481e5c89e4eSSatish Balay -  args - the command line arguments
482e5c89e4eSSatish Balay 
483e5c89e4eSSatish Balay    Level: intermediate
484e5c89e4eSSatish Balay 
485e5c89e4eSSatish Balay    Notes:
486e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
487e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
488e5c89e4eSSatish Balay 
489f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
490f177e3b1SBarry Smith 
491db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`
492e5c89e4eSSatish Balay 
493e5c89e4eSSatish Balay @*/
4947087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
495e5c89e4eSSatish Balay {
496e5c89e4eSSatish Balay   PetscFunctionBegin;
49739a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
498e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
499e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
500e5c89e4eSSatish Balay   PetscFunctionReturn(0);
501e5c89e4eSSatish Balay }
502e5c89e4eSSatish Balay 
503793721a6SBarry Smith /*@C
504793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
505793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
506793721a6SBarry Smith 
507793721a6SBarry Smith    Not Collective
508793721a6SBarry Smith 
509793721a6SBarry Smith    Output Parameters:
510793721a6SBarry Smith .  args - the command line arguments
511793721a6SBarry Smith 
512793721a6SBarry Smith    Level: intermediate
513793721a6SBarry Smith 
514793721a6SBarry Smith    Notes:
515793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
516793721a6SBarry Smith 
517db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`
518793721a6SBarry Smith 
519793721a6SBarry Smith @*/
5207087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
521793721a6SBarry Smith {
522793721a6SBarry Smith   PetscInt i,argc = PetscGlobalArgc;
523793721a6SBarry Smith 
524793721a6SBarry Smith   PetscFunctionBegin;
52539a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
5262d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
5279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(argc,args));
5289566063dSJacob Faibussowitsch   for (i=0; i<argc-1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]));
5292d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
530793721a6SBarry Smith   PetscFunctionReturn(0);
531793721a6SBarry Smith }
532793721a6SBarry Smith 
533793721a6SBarry Smith /*@C
534793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
535793721a6SBarry Smith 
536793721a6SBarry Smith    Not Collective
537793721a6SBarry Smith 
538793721a6SBarry Smith    Output Parameters:
539793721a6SBarry Smith .  args - the command line arguments
540793721a6SBarry Smith 
541793721a6SBarry Smith    Level: intermediate
542793721a6SBarry Smith 
543db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()`
544793721a6SBarry Smith 
545793721a6SBarry Smith @*/
5467087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
547793721a6SBarry Smith {
54839a651e2SJacob Faibussowitsch   PetscFunctionBegin;
54939a651e2SJacob Faibussowitsch   if (args) {
550793721a6SBarry Smith     PetscInt i = 0;
551793721a6SBarry Smith 
5529566063dSJacob Faibussowitsch     while (args[i]) PetscCall(PetscFree(args[i++]));
5539566063dSJacob Faibussowitsch     PetscCall(PetscFree(args));
55439a651e2SJacob Faibussowitsch   }
555793721a6SBarry Smith   PetscFunctionReturn(0);
556793721a6SBarry Smith }
557793721a6SBarry Smith 
55827104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
55930befbd2SBarry Smith #include <petscconfiginfo.h>
56030befbd2SBarry Smith 
56195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
56211525c0dSBarry Smith {
56327104ee2SJacob Faibussowitsch   PetscFunctionBegin;
56411525c0dSBarry Smith   if (!PetscGlobalRank) {
56530befbd2SBarry Smith     char      cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
56611525c0dSBarry Smith     int       port;
567ffbd1cfbSBarry Smith     PetscBool flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
56811525c0dSBarry Smith     size_t    applinelen,introlen;
569ffbd1cfbSBarry Smith     char      sawsurl[256];
57011525c0dSBarry Smith 
5719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL,NULL,"-saws_log",&flg));
57211525c0dSBarry Smith     if (flg) {
57311525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
57411525c0dSBarry Smith 
5759566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL));
57611525c0dSBarry Smith       if (sawslog[0]) {
577792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
57811525c0dSBarry Smith       } else {
579792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile,(NULL));
58011525c0dSBarry Smith       }
58111525c0dSBarry Smith     }
5829566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg));
58311525c0dSBarry Smith     if (flg) {
584792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Use_HTTPS,(cert));
58511525c0dSBarry Smith     }
5869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL));
587ffbd1cfbSBarry Smith     if (selectport) {
588792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_Available_Port,(&port));
589792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port,(port));
590ffbd1cfbSBarry Smith     } else {
5919566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg));
59211525c0dSBarry Smith       if (flg) {
593792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port,(port));
59411525c0dSBarry Smith       }
595ffbd1cfbSBarry Smith     }
5969566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg));
59711525c0dSBarry Smith     if (flg) {
598792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Document_Root,(root));
5999566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(root,".",&rootlocal));
6009c1e0ce8SBarry Smith     } else {
6019566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL,NULL,"-saws_options",&flg));
6029c1e0ce8SBarry Smith       if (flg) {
6039566063dSJacob Faibussowitsch         PetscCall(PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root)));
604792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Document_Root,(root));
6059c1e0ce8SBarry Smith       }
60611525c0dSBarry Smith     }
6079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2));
60811525c0dSBarry Smith     if (flg2) {
60911525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
61028b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
6119566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root));
6129566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(jsdir,'r',&flg));
61328b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
614792fecdfSBarry Smith       PetscCallSAWs(SAWs_Push_Local_Header,());
61511525c0dSBarry Smith     }
6169566063dSJacob Faibussowitsch     PetscCall(PetscGetProgramName(programname,sizeof(programname)));
6179566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(help,&applinelen));
61811525c0dSBarry Smith     introlen   = 4096 + applinelen;
61930a8c9c0SSurtai Han     applinelen += 1024;
6209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(applinelen,&appline));
6219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(introlen,&intro));
62211525c0dSBarry Smith 
62311525c0dSBarry Smith     if (rootlocal) {
6249566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline,applinelen,"%s.c.html",programname));
6259566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(appline,'r',&rootlocal));
62611525c0dSBarry Smith     }
6279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetAll(NULL,&options));
62811525c0dSBarry Smith     if (rootlocal && help) {
6299566063dSJacob Faibussowitsch       PetscCall(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));
63011525c0dSBarry Smith     } else if (help) {
6319566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help));
63211525c0dSBarry Smith     } else {
6339566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options));
63411525c0dSBarry Smith     }
6359566063dSJacob Faibussowitsch     PetscCall(PetscFree(options));
6369566063dSJacob Faibussowitsch     PetscCall(PetscGetVersion(version,sizeof(version)));
6379566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(intro,introlen,"<body>\n"
638a17b96a8SKyle 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"
639df62c222SBarry 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"
6405f80ce2aSJacob Faibussowitsch                           "%s",version,petscconfigureoptions,appline));
641792fecdfSBarry Smith     PetscCallSAWs(SAWs_Push_Body,("index.html",0,intro));
6429566063dSJacob Faibussowitsch     PetscCall(PetscFree(intro));
6439566063dSJacob Faibussowitsch     PetscCall(PetscFree(appline));
644ffbd1cfbSBarry Smith     if (selectport) {
645aa573868SBarry Smith       PetscBool silent;
6467d812c46SBarry Smith 
6477d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
64839a651e2SJacob Faibussowitsch       while (SAWs_Initialize()) {
649792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_Available_Port,(&port));
650792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port,(port));
6517d812c46SBarry Smith       }
6527d812c46SBarry Smith 
6539566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL));
654aa573868SBarry Smith       if (!silent) {
655792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
6569566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl));
657ffbd1cfbSBarry Smith       }
6587d812c46SBarry Smith     } else {
659792fecdfSBarry Smith       PetscCallSAWs(SAWs_Initialize,());
660aa573868SBarry Smith     }
6619566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister("@TechReport{ saws,\n"
6620af79b04SBarry Smith                                    "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6630af79b04SBarry Smith                                    "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6640af79b04SBarry Smith                                    "  Institution = {Argonne National Laboratory},\n"
6655f80ce2aSJacob Faibussowitsch                                    "  Year   = 2013\n}\n",NULL));
66611525c0dSBarry Smith   }
66711525c0dSBarry Smith   PetscFunctionReturn(0);
66811525c0dSBarry Smith }
66911525c0dSBarry Smith #endif
67011525c0dSBarry Smith 
6714dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6724dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6734dfee713SSatish Balay {
6744dfee713SSatish Balay   PetscFunctionBegin;
6754dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6764dfee713SSatish Balay   /* see MPI.py for details on this bug */
6774dfee713SSatish Balay   (void) setenv("HWLOC_COMPONENTS","-x86",1);
6784dfee713SSatish Balay #endif
6794dfee713SSatish Balay   PetscFunctionReturn(0);
6804dfee713SSatish Balay }
6814dfee713SSatish Balay 
682a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
683a56f64adSBarry Smith #include <adios.h>
68422580e64SBarry Smith #include <adios_read.h>
6857b56e58cSSatish Balay int64_t Petsc_adios_group;
686a56f64adSBarry Smith #endif
687a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
688cd1b99a6SBarry Smith #include <omp.h>
689f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
690cd1b99a6SBarry Smith #endif
691a56f64adSBarry Smith 
692a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_DEVICE)
693a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
694a4af0ceeSJacob Faibussowitsch #  if PetscDefined(HAVE_CUDA)
695a4af0ceeSJacob Faibussowitsch // REMOVE ME
696a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
697a4af0ceeSJacob Faibussowitsch #  endif
698a4af0ceeSJacob Faibussowitsch #  if PetscDefined(HAVE_HIP)
699a4af0ceeSJacob Faibussowitsch // REMOVE ME
700a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
701a4af0ceeSJacob Faibussowitsch #  endif
702a4af0ceeSJacob Faibussowitsch #endif
703a4af0ceeSJacob Faibussowitsch 
70427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
705bc8a24ecSBarry Smith #include <dlfcn.h>
706bc8a24ecSBarry Smith #endif
707a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
708a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
709a4af0ceeSJacob Faibussowitsch #endif
710a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
711a4af0ceeSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViennaCLInit();
712a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE;
713a4af0ceeSJacob Faibussowitsch #endif
714bc8a24ecSBarry Smith 
715660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE;
716660278c0SBarry Smith 
71785649d77SJunchao Zhang /*
71885649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
71985649d77SJunchao Zhang 
72085649d77SJunchao Zhang   prog:     program name
72102101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
72285649d77SJunchao Zhang   help:     program help message
72302101c96SSatish Balay   ftn:      is it called from Fortran initilization (petscinitializef_)?
72485649d77SJunchao Zhang   readarguments,len: used when fortran is true
72585649d77SJunchao Zhang */
72602101c96SSatish Balay PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char* prog,const char* file,const char *help,PetscBool ftn,PetscBool readarguments,PetscInt len)
72785649d77SJunchao Zhang {
72885649d77SJunchao Zhang   PetscMPIInt size;
72985649d77SJunchao Zhang   PetscBool   flg = PETSC_TRUE;
73085649d77SJunchao Zhang   char        hostname[256];
73185649d77SJunchao Zhang 
73227104ee2SJacob Faibussowitsch   PetscFunctionBegin;
73327104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
73439a651e2SJacob Faibussowitsch   /* these must be initialized in a routine, not as a constant declaration */
73539a651e2SJacob Faibussowitsch   PETSC_STDOUT = stdout;
73639a651e2SJacob Faibussowitsch   PETSC_STDERR = stderr;
73739a651e2SJacob Faibussowitsch 
7389566063dSJacob Faibussowitsch   /* PetscCall can be used from now */
73939a651e2SJacob Faibussowitsch   PetscErrorHandlingInitialized = PETSC_TRUE;
74039a651e2SJacob Faibussowitsch 
74185649d77SJunchao Zhang   /*
74285649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
74385649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
74485649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
74585649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
74685649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
74785649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
74885649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
74985649d77SJunchao Zhang       listed above and since that time are compatible.
75085649d77SJunchao Zhang 
75185649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
75285649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
75385649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
75485649d77SJunchao Zhang       to perform the checking.
75585649d77SJunchao Zhang 
75685649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
75785649d77SJunchao Zhang 
75885649d77SJunchao Zhang       Questions:
75985649d77SJunchao Zhang 
76085649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
76185649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
76285649d77SJunchao Zhang   */
76385649d77SJunchao Zhang 
76485649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
76585649d77SJunchao Zhang   {
76685649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
76785649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
76839a651e2SJacob Faibussowitsch 
7699566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength));
77085649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
77185649d77SJunchao Zhang #if defined(MPICH_VERSION)
77285649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000
77385649d77SJunchao Zhang     {
77485649d77SJunchao Zhang       char      *ver,*lf;
77585649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
77639a651e2SJacob Faibussowitsch 
7779566063dSJacob Faibussowitsch       PetscCall(PetscStrstr(mpilibraryversion,"MPICH Version:",&ver));
77839a651e2SJacob Faibussowitsch       if (ver) {
7799566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(ver,'\n',&lf));
78039a651e2SJacob Faibussowitsch         if (lf) {
78185649d77SJunchao Zhang           *lf = 0;
7829566063dSJacob Faibussowitsch           PetscCall(PetscStrendswith(ver,MPICH_VERSION,&flg));
78385649d77SJunchao Zhang         }
78485649d77SJunchao Zhang       }
78585649d77SJunchao Zhang       if (!flg) {
786d5b396e8SSatish Balay         PetscCall(PetscInfo(NULL,"PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n",mpilibraryversion,MPICH_VERSION));
78785649d77SJunchao Zhang         flg = PETSC_TRUE;
78885649d77SJunchao Zhang       }
78985649d77SJunchao Zhang     }
79085649d77SJunchao Zhang #endif
79185649d77SJunchao 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?) */
79285649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION)
79385649d77SJunchao Zhang     {
79485649d77SJunchao Zhang       char *ver,bs[MPI_MAX_LIBRARY_VERSION_STRING],*bsf;
79585649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
79685649d77SJunchao Zhang #define PSTRSZ 2
79785649d77SJunchao Zhang       char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI","FUJITSU MPI"};
79885649d77SJunchao Zhang       char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v","Library "};
79985649d77SJunchao Zhang       int i;
80085649d77SJunchao Zhang       for (i=0; i<PSTRSZ; i++) {
8019566063dSJacob Faibussowitsch         PetscCall(PetscStrstr(mpilibraryversion,ompistr1[i],&ver));
80239a651e2SJacob Faibussowitsch         if (ver) {
8039566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(bs,MPI_MAX_LIBRARY_VERSION_STRING,"%s%d.%d",ompistr2[i],OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION));
8049566063dSJacob Faibussowitsch           PetscCall(PetscStrstr(ver,bs,&bsf));
80539a651e2SJacob Faibussowitsch           if (bsf) flg = PETSC_TRUE;
80685649d77SJunchao Zhang           break;
80785649d77SJunchao Zhang         }
80885649d77SJunchao Zhang       }
80985649d77SJunchao Zhang       if (!flg) {
8107d3de750SJacob 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);
81185649d77SJunchao Zhang         flg = PETSC_TRUE;
81285649d77SJunchao Zhang       }
81385649d77SJunchao Zhang     }
81485649d77SJunchao Zhang #endif
81585649d77SJunchao Zhang   }
81685649d77SJunchao Zhang #endif
81785649d77SJunchao Zhang 
8186f5d4113SSatish Balay #if PetscDefined(HAVE_DLSYM) && defined(__USE_GNU)
81985649d77SJunchao 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 */
82039a651e2SJacob Faibussowitsch   PetscCheck(!dlsym(RTLD_DEFAULT,"ompi_mpi_init") || !dlsym(RTLD_DEFAULT,"MPID_Abort"),PETSC_COMM_SELF,PETSC_ERR_MPI_LIB_INCOMP,"Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly");
82185649d77SJunchao Zhang #endif
82285649d77SJunchao Zhang 
82385649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
82485649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
82585649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
82685649d77SJunchao Zhang #endif
82785649d77SJunchao Zhang 
8289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
82985649d77SJunchao Zhang 
83085649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
83185649d77SJunchao Zhang 
8329566063dSJacob Faibussowitsch   PetscCall(PetscSetProgramName(prog));
8339566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen));
8349566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout));
8359566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr));
8369566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscCommSpinLock));
83785649d77SJunchao Zhang 
83885649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
8399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN));
84085649d77SJunchao Zhang 
84185649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
8429566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS));
8439566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE));
84485649d77SJunchao Zhang   }
84585649d77SJunchao Zhang 
84685649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
8479566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
84885649d77SJunchao Zhang 
8499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank));
8509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize));
85185649d77SJunchao Zhang 
85285649d77SJunchao Zhang   MPIU_BOOL = MPI_INT;
85385649d77SJunchao Zhang   MPIU_ENUM = MPI_INT;
85485649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
85585649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
85685649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
85785649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
85885649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
85985649d77SJunchao Zhang #endif
86039a651e2SJacob Faibussowitsch   else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP_SYS,"Could not find MPI type for size_t");
86185649d77SJunchao Zhang 
86285649d77SJunchao Zhang   /*
86385649d77SJunchao Zhang      Initialized the global complex variable; this is because with
86485649d77SJunchao Zhang      shared libraries the constructors for global variables
86585649d77SJunchao Zhang      are not called; at least on IRIX.
86685649d77SJunchao Zhang   */
86785649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
86885649d77SJunchao Zhang   {
86985649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
87085649d77SJunchao Zhang     PetscComplex ic(0.0,1.0);
87185649d77SJunchao Zhang     PETSC_i = ic;
87285649d77SJunchao Zhang #else
87385649d77SJunchao Zhang     PETSC_i = _Complex_I;
87485649d77SJunchao Zhang #endif
87585649d77SJunchao Zhang   }
87685649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
87785649d77SJunchao Zhang 
87885649d77SJunchao Zhang   /*
87985649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
88085649d77SJunchao Zhang      half of the entries and maxes the second half.
88185649d77SJunchao Zhang   */
8829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP));
88385649d77SJunchao Zhang 
884*613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
8859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128));
8869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128));
8879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128));
8889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128));
889*613bf2b2SPierre Jolivet #if !defined(PETSC_USE_REAL___FLOAT128)
890*613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscSum_Local,1,&MPIU_SUM___FLOAT128));
89185649d77SJunchao Zhang #endif
892*613bf2b2SPierre Jolivet #endif
893*613bf2b2SPierre Jolivet #if defined(PETSC_USE_REAL___FP16)
8949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16));
8959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FP16));
89685649d77SJunchao Zhang #endif
89785649d77SJunchao Zhang 
89885649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
8999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(PetscSum_Local,1,&MPIU_SUM));
900*613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMax_Local,1,&MPIU_MAX));
901*613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMin_Local,1,&MPIU_MIN));
90285649d77SJunchao Zhang #endif
90385649d77SJunchao Zhang 
9049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR));
9059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR));
90685649d77SJunchao Zhang 
90785649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
90885649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
90985649d77SJunchao Zhang   {
91085649d77SJunchao Zhang     struct PetscRealInt { PetscReal v; PetscInt i; };
91185649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
91285649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscRealInt,v),offsetof(struct PetscRealInt,i)};
91385649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_REAL,MPIU_INT}, tmpStruct;
91485649d77SJunchao Zhang 
9159566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct));
9169566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscRealInt),&MPIU_REAL_INT));
9179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT));
91985649d77SJunchao Zhang   }
92085649d77SJunchao Zhang   {
92185649d77SJunchao Zhang     struct PetscScalarInt { PetscScalar v; PetscInt i; };
92285649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
92385649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscScalarInt,v),offsetof(struct PetscScalarInt,i)};
92485649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_SCALAR,MPIU_INT}, tmpStruct;
92585649d77SJunchao Zhang 
9269566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct));
9279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscScalarInt),&MPIU_SCALAR_INT));
9289566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9299566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
93085649d77SJunchao Zhang   }
93185649d77SJunchao Zhang #endif
93285649d77SJunchao Zhang 
93385649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT));
9359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2INT));
93685649d77SJunchao Zhang #endif
9379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4,MPI_INT,&MPI_4INT));
9389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPI_4INT));
9399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4,MPIU_INT,&MPIU_4INT));
9409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_4INT));
94185649d77SJunchao Zhang 
94285649d77SJunchao Zhang   /*
94385649d77SJunchao Zhang      Attributes to be set on PETSc communicators
94485649d77SJunchao Zhang   */
9459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0));
9469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0));
9479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0));
9489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0));
94985649d77SJunchao Zhang 
95085649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
9519566063dSJacob Faibussowitsch   if (ftn) PetscCall(PetscInitFortran_Private(readarguments,file,len));
95285649d77SJunchao Zhang   else
95385649d77SJunchao Zhang #endif
9549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInsert(NULL,&PetscGlobalArgc,&PetscGlobalArgs,file));
95585649d77SJunchao Zhang 
95685649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9579566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
95885649d77SJunchao Zhang 
95985649d77SJunchao Zhang   /*
96085649d77SJunchao Zhang      Check system options and print help
96185649d77SJunchao Zhang   */
9629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCheckInitial_Private(help));
96385649d77SJunchao Zhang 
964a4af0ceeSJacob Faibussowitsch   /*
965a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
966a4af0ceeSJacob Faibussowitsch 
967a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
968a4af0ceeSJacob Faibussowitsch    1. MPI initialized
969a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
970a4af0ceeSJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up its own SIGSEV handler via
971a4af0ceeSJacob Faibussowitsch       the push/pop interface.
972a4af0ceeSJacob Faibussowitsch   */
973a2158755SJunchao Zhang #if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP) || PetscDefined(HAVE_SYCL))
9749566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
975a4af0ceeSJacob Faibussowitsch #endif
976a4af0ceeSJacob Faibussowitsch 
977a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
978a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
9799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-log_summary",&flg));
9809566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsHasName(NULL,NULL,"-log_view",&flg));
9819566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsGetBool(NULL,NULL,"-viennacl_synchronize",&flg,NULL));
982a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
9839566063dSJacob Faibussowitsch   PetscCall(PetscViennaCLInit());
984a4af0ceeSJacob Faibussowitsch #endif
985a4af0ceeSJacob Faibussowitsch 
986a4af0ceeSJacob Faibussowitsch   /*
987a4af0ceeSJacob Faibussowitsch      Creates the logging data structures; this is enabled even if logging is not turned on
988a4af0ceeSJacob Faibussowitsch      This is the last thing we do before returning to the user code to prevent having the
989a4af0ceeSJacob Faibussowitsch      logging numbers contaminated by any startup time associated with MPI
990a4af0ceeSJacob Faibussowitsch   */
991a4af0ceeSJacob Faibussowitsch #if defined(PETSC_USE_LOG)
9929566063dSJacob Faibussowitsch   PetscCall(PetscLogInitialize());
993a4af0ceeSJacob Faibussowitsch #endif
994a4af0ceeSJacob Faibussowitsch 
9959566063dSJacob Faibussowitsch   PetscCall(PetscCitationsInitialize());
99685649d77SJunchao Zhang 
99785649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
9989566063dSJacob Faibussowitsch   PetscCall(PetscInitializeSAWs(ftn ? NULL : help));
99927104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
10009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-stack_view",&flg));
10019566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscStackViewSAWs());
100285649d77SJunchao Zhang #endif
100385649d77SJunchao Zhang 
100485649d77SJunchao Zhang   /*
100585649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
100685649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
100785649d77SJunchao Zhang   */
10089566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_DynamicLibraries());
100985649d77SJunchao Zhang 
10109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
10119566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL,"PETSc successfully started: number of processors = %d\n",size));
10129566063dSJacob Faibussowitsch   PetscCall(PetscGetHostName(hostname,256));
10139566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL,"Running on machine: %s\n",hostname));
101485649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
101585649d77SJunchao Zhang   {
101685649d77SJunchao Zhang     PetscBool       omp_view_flag;
101785649d77SJunchao Zhang     char           *threads = getenv("OMP_NUM_THREADS");
101885649d77SJunchao Zhang 
101985649d77SJunchao Zhang     if (threads) {
10209566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL,"Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n",threads));
102185649d77SJunchao Zhang       (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
102285649d77SJunchao Zhang     } else {
10232f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt) omp_get_max_threads();
10249566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL,"Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n",PetscNumOMPThreads));
102585649d77SJunchao Zhang     }
1026d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");
10279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-omp_num_threads","Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS","None",PetscNumOMPThreads,&PetscNumOMPThreads,&flg));
10289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag));
1029d0609cedSBarry Smith     PetscOptionsEnd();
103085649d77SJunchao Zhang     if (flg) {
10319566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL,"Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n",PetscNumOMPThreads));
103285649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
103385649d77SJunchao Zhang     }
103485649d77SJunchao Zhang     if (omp_view_flag) {
10359566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %" PetscInt_FMT "\n",PetscNumOMPThreads));
103685649d77SJunchao Zhang     }
103785649d77SJunchao Zhang   }
103885649d77SJunchao Zhang #endif
103985649d77SJunchao Zhang 
104085649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
104185649d77SJunchao Zhang   /*
104285649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
104385649d77SJunchao Zhang 
104485649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
104585649d77SJunchao Zhang   */
10469566063dSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL));
104785649d77SJunchao Zhang #endif
104885649d77SJunchao Zhang 
104985649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10509566063dSJacob Faibussowitsch   PetscCall(PetscFPTCreate(10000));
105185649d77SJunchao Zhang #endif
105285649d77SJunchao Zhang 
105385649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
105485649d77SJunchao Zhang   {
105585649d77SJunchao Zhang     PetscViewer viewer;
10569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg));
105785649d77SJunchao Zhang     if (flg) {
10589566063dSJacob Faibussowitsch       PetscCall(PetscProcessPlacementView(viewer));
10599566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
106085649d77SJunchao Zhang     }
106185649d77SJunchao Zhang   }
106285649d77SJunchao Zhang #endif
106385649d77SJunchao Zhang 
106485649d77SJunchao Zhang   flg  = PETSC_TRUE;
10659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL));
10669566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE));
106785649d77SJunchao Zhang 
106885649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10699566063dSJacob Faibussowitsch   PetscCall(adios_init_noxml(PETSC_COMM_WORLD));
10709566063dSJacob Faibussowitsch   PetscCall(adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default));
10719566063dSJacob Faibussowitsch   PetscCall(adios_select_method(Petsc_adios_group,"MPI","",""));
10729566063dSJacob Faibussowitsch   PetscCall(adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,""));
107385649d77SJunchao Zhang #endif
107485649d77SJunchao Zhang 
107585649d77SJunchao Zhang #if defined(__VALGRIND_H)
107685649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND? PETSC_TRUE: PETSC_FALSE;
107785649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
10789566063dSJacob Faibussowitsch   if (PETSC_RUNNING_ON_VALGRIND) PetscCall(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"));
107985649d77SJunchao Zhang #endif
108085649d77SJunchao Zhang #endif
108185649d77SJunchao Zhang   /*
108285649d77SJunchao Zhang       Set flag that we are completely initialized
108385649d77SJunchao Zhang   */
108485649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
108585649d77SJunchao Zhang 
10869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-python",&flg));
10879566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPythonInitialize(NULL,NULL));
1088f1f2ae84SBarry Smith 
1089f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL,NULL,"-mpi_linear_solver_server",&flg));
1090f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin());
1091f1f2ae84SBarry Smith   else PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_SUP,"PETSc configured using -with-single-library=0; -mpi_linear_solver_server not supported in that case");
109227104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
109385649d77SJunchao Zhang }
109485649d77SJunchao Zhang 
1095e5c89e4eSSatish Balay /*@C
1096e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1097e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
1098e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1099e5c89e4eSSatish Balay    your program -- usually the very first line!
1100e5c89e4eSSatish Balay 
1101e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
1102e5c89e4eSSatish Balay 
1103e5c89e4eSSatish Balay    Input Parameters:
1104e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1105e5c89e4eSSatish Balay .  args - the command line arguments
1106be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1107be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1108be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1109c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
11100298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1111e5c89e4eSSatish Balay 
111205827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
111305827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
111405827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
111505827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
111605827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1117e5c89e4eSSatish Balay 
1118e5c89e4eSSatish Balay    Options Database Keys:
11197ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
11207ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1121e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11228a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
11238a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
11248a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
11251af3601dSBarry Smith .  -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr
11268a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1127bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1128e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1129e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1130e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
113179dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
113279dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
113379dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
1134aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
113592f119d6SBarry 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
113692f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
113792f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1138608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
113992f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1140e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1141e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1142e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1143e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1144e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
11450841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1146e5c89e4eSSatish Balay 
1147c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1148c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1149c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1150c5b5d8d5SVaclav Hapla -  -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet()
1151c5b5d8d5SVaclav Hapla 
1152c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1153c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1154c5b5d8d5SVaclav Hapla    They can be used also in option files.
1155c5b5d8d5SVaclav Hapla 
1156c5b5d8d5SVaclav Hapla    See PetscOptionsMonitorSet() to do monitoring programmatically.
1157c5b5d8d5SVaclav Hapla 
1158e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1159a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1160fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
1161217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1162217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1163495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1164e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
11659a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
1166156b51fbSBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each event, see PetscLogView().
1167156b51fbSBarry Smith .  -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see PetscLogView().
11689a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1169495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1170f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1171495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
1172495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
1173c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
117487c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
1175c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1176495fc317SBarry Smith 
11776a77f485SPierre Jolivet     Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time
1178e5c89e4eSSatish Balay 
1179ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1180ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1181ffbd1cfbSBarry 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
1182ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1183ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1184ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1185ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1186ffbd1cfbSBarry Smith 
1187e5c89e4eSSatish Balay    Environmental Variables:
1188e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
1189e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
1190e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
11914a971ea4SToby Isaac .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
11924a971ea4SToby Isaac .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1193e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
1194e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
1195e5c89e4eSSatish Balay 
1196e5c89e4eSSatish Balay    Level: beginner
1197e5c89e4eSSatish Balay 
1198e5c89e4eSSatish Balay    Notes:
1199e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
1200e5c89e4eSSatish Balay    it before PetscInitialize().
1201e5c89e4eSSatish Balay 
1202e5c89e4eSSatish Balay    Fortran Version:
1203e5c89e4eSSatish Balay    In Fortran this routine has the format
1204e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
1205e5c89e4eSSatish Balay 
1206e5c89e4eSSatish Balay +  ierr - error return code
1207c5b5d8d5SVaclav Hapla -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
1208c5b5d8d5SVaclav Hapla           Use PETSC_NULL_CHARACTER to not check for code specific file.
1209c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
1210e5c89e4eSSatish Balay 
1211e5c89e4eSSatish Balay    Important Fortran Note:
12120eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
12130298fd71SBarry Smith    null character string; you CANNOT just use NULL as
1214a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
1215e5c89e4eSSatish Balay 
121601cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
121701cb0274SBarry Smith    calling PetscInitialize().
1218e5c89e4eSSatish Balay 
1219db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()`
1220e5c89e4eSSatish Balay 
1221e5c89e4eSSatish Balay @*/
12227087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
1223e5c89e4eSSatish Balay {
122485649d77SJunchao Zhang   PetscMPIInt  flag;
122585649d77SJunchao Zhang   const char  *prog = "Unknown Name";
1226e5c89e4eSSatish Balay 
122727104ee2SJacob Faibussowitsch   PetscFunctionBegin;
122827104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
12299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Initialized(&flag));
1230e5c89e4eSSatish Balay   if (!flag) {
123139a651e2SJacob Faibussowitsch     PetscCheck(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");
12329566063dSJacob Faibussowitsch     PetscCall(PetscPreMPIInit_Private());
12335e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12345e765c61SJed Brown     {
123539a651e2SJacob Faibussowitsch       PetscMPIInt PETSC_UNUSED provided;
12369566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided));
12375e765c61SJed Brown     }
12385e765c61SJed Brown #else
12399566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Init(argc,args));
12405e765c61SJed Brown #endif
1241e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1242e5c89e4eSSatish Balay   }
12434dfee713SSatish Balay 
124485649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1245e5c89e4eSSatish Balay   if (argc && args) {
1246e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1247e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1248e5c89e4eSSatish Balay   }
12499566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_Common(prog,file,help,PETSC_FALSE/*C*/,PETSC_FALSE,0));
125027104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
1251e5c89e4eSSatish Balay }
1252e5c89e4eSSatish Balay 
1253a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
125495c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1255ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1256ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
125795c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
12584097062eSBarry Smith #endif
1259e5c89e4eSSatish Balay 
1260008a6e76SBarry Smith /*
1261008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1262008a6e76SBarry Smith */
1263008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1264008a6e76SBarry Smith {
1265008a6e76SBarry Smith   PetscFunctionBegin;
1266*613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
12679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128));
12689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128));
1269*613bf2b2SPierre Jolivet #if !defined(PETSC_USE_REAL___FLOAT128)
1270*613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_SUM___FLOAT128));
1271008a6e76SBarry Smith #endif
1272*613bf2b2SPierre Jolivet #endif
1273*613bf2b2SPierre Jolivet #if defined(PETSC_USE_REAL___FP16)
12749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FP16));
1275008a6e76SBarry Smith #endif
1276008a6e76SBarry Smith 
1277de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
12789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_SUM));
1279*613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MAX));
1280*613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MIN));
1281008a6e76SBarry Smith #endif
1282008a6e76SBarry Smith 
12839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR));
12849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT));
12859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT));
128640df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
12879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2INT));
1288008a6e76SBarry Smith #endif
12899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPI_4INT));
12909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_4INT));
12919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
1292008a6e76SBarry Smith   PetscFunctionReturn(0);
1293008a6e76SBarry Smith }
1294008a6e76SBarry Smith 
1295a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
1296a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1297a4af0ceeSJacob Faibussowitsch #endif
1298a4af0ceeSJacob Faibussowitsch 
1299e5c89e4eSSatish Balay /*@C
1300e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1301e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1302e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1303e5c89e4eSSatish Balay 
1304e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1305e5c89e4eSSatish Balay 
1306e5c89e4eSSatish Balay    Options Database Keys:
130726a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1308e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
13097eb1d149SBarry 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
1310e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
131192f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1312e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
131392f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1314e5c89e4eSSatish Balay 
1315e5c89e4eSSatish Balay    Level: beginner
1316e5c89e4eSSatish Balay 
1317e5c89e4eSSatish Balay    Note:
1318e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1319e5c89e4eSSatish Balay 
1320db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()`
1321e5c89e4eSSatish Balay @*/
13227087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1323e5c89e4eSSatish Balay {
13244bb5149bSJed Brown   PetscMPIInt rank;
1325a8d2bbe5SBarry Smith   PetscInt    nopt;
13262bf49c77SBarry Smith   PetscBool   flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1327dff31646SBarry Smith   PetscBool   flg;
132810463e74SBarry Smith #if defined(PETSC_USE_LOG)
132910463e74SBarry Smith   char        mname[PETSC_MAX_PATH_LEN];
133010463e74SBarry Smith #endif
1331e5c89e4eSSatish Balay 
13323cbbc5ffSBarry Smith   PetscFunctionBegin;
133339a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscInitialize() must be called before PetscFinalize()");
13349566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL,"PetscFinalize() called\n"));
1335b022a5c1SBarry Smith 
1336f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL,NULL,"-mpi_linear_solver_server",&flg));
1337f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd());
1338f1f2ae84SBarry Smith 
13399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1340a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13419566063dSJacob Faibussowitsch   PetscCall(adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE));
13429566063dSJacob Faibussowitsch   PetscCall(adios_finalize(rank));
1343a56f64adSBarry Smith #endif
13449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-citations",&flg));
1345dff31646SBarry Smith   if (flg) {
13461f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
13471f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
13481f817a21SBarry Smith 
13499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL));
13501f817a21SBarry Smith     if (filename[0]) {
13519566063dSJacob Faibussowitsch       PetscCall(PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd));
13521f817a21SBarry Smith     }
13539566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferGet(PetscCitationsList,1,&cits));
1354dff31646SBarry Smith     cits[0] = 0;
13559566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList,&cits));
13569566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n"));
13579566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n"));
13589566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits));
13599566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n"));
13609566063dSJacob Faibussowitsch     PetscCall(PetscFClose(PETSC_COMM_WORLD,fd));
13619566063dSJacob Faibussowitsch     PetscCall(PetscFree(cits));
1362dff31646SBarry Smith   }
13639566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&PetscCitationsList));
1364dff31646SBarry Smith 
1365c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
136604102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
136704102261SBarry Smith   {
136804102261SBarry Smith     PetscInt nmax = 2;
136904102261SBarry Smith     char     **buffs;
13709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2,&buffs));
13719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1));
137204102261SBarry Smith     if (flg1) {
137328b400f6SJacob Faibussowitsch       PetscCheck(nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
137404102261SBarry Smith       if (nmax == 1) {
13759566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(128,&buffs[1]));
13769566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1],32));
13779566063dSJacob Faibussowitsch         PetscCall(PetscStrcat(buffs[1]," has completed"));
137804102261SBarry Smith       }
13799566063dSJacob Faibussowitsch       PetscCall(PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL));
13809566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
13819566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
138204102261SBarry Smith     }
13839566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
138404102261SBarry Smith   }
1385763ec7b1SBarry Smith   {
1386763ec7b1SBarry Smith     PetscInt nmax = 2;
1387763ec7b1SBarry Smith     char     **buffs;
13889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2,&buffs));
13899566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1));
1390763ec7b1SBarry Smith     if (flg1) {
139128b400f6SJacob Faibussowitsch       PetscCheck(nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1392763ec7b1SBarry Smith       if (nmax == 1) {
13939566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(128,&buffs[1]));
13949566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1],32));
13959566063dSJacob Faibussowitsch         PetscCall(PetscStrcat(buffs[1]," has completed"));
1396763ec7b1SBarry Smith       }
13979566063dSJacob Faibussowitsch       PetscCall(PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL));
13989566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
13999566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
1400763ec7b1SBarry Smith     }
14019566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
1402763ec7b1SBarry Smith   }
140304102261SBarry Smith #endif
140404102261SBarry Smith 
14052d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14069566063dSJacob Faibussowitsch   PetscCall(PetscFPTDestroy());
14072d53ad75SBarry Smith #endif
14082d53ad75SBarry Smith 
1409e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1410dff31646SBarry Smith   flg = PETSC_FALSE;
14119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL));
14121baa6e33SBarry Smith   if (flg) PetscCall(PetscOptionsSAWsDestroy());
1413d5649816SBarry Smith #endif
1414d5649816SBarry Smith 
1415681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1416681455b2SBarry Smith   flg1 = PETSC_FALSE;
14179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL));
1418681455b2SBarry Smith   if (flg1) {
1419681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
14209566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL));
1421681455b2SBarry Smith   }
1422681455b2SBarry Smith #endif
1423681455b2SBarry Smith 
142467584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL));
1426e5c89e4eSSatish Balay   if (!flg2) {
142790d69ab7SBarry Smith     flg2 = PETSC_FALSE;
14289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL));
1429e5c89e4eSSatish Balay   }
1430e5c89e4eSSatish Balay   if (flg2) {
14319566063dSJacob Faibussowitsch     PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n"));
1432e5c89e4eSSatish Balay   }
143367584ceeSBarry Smith #endif
1434e5c89e4eSSatish Balay 
1435e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
143690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL));
1438e5c89e4eSSatish Balay   if (flg1) {
1439e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
14409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD));
14419566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops));
1442e5c89e4eSSatish Balay   }
1443e5c89e4eSSatish Balay #endif
1444e5c89e4eSSatish Balay 
1445e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1446e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1447e5c89e4eSSatish Balay   mname[0] = 0;
14489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1));
1449e5c89e4eSSatish Balay   if (flg1) {
14509566063dSJacob Faibussowitsch     if (mname[0]) PetscCall(PetscLogMPEDump(mname));
14519566063dSJacob Faibussowitsch     else          PetscCall(PetscLogMPEDump(0));
1452e5c89e4eSSatish Balay   }
1453e5c89e4eSSatish Balay #endif
1454356e5837SBarry Smith #endif
1455a297a907SKarl Rupp 
1456dd710f27SBarry Smith   /*
1457dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1458dd710f27SBarry Smith   */
14599566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1460dd710f27SBarry Smith 
1461356e5837SBarry Smith #if defined(PETSC_USE_LOG)
14629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE));
14639566063dSJacob Faibussowitsch   PetscCall(PetscLogViewFromOptions());
14649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPopGetViewerOff());
146587c3beb6SLisandro Dalcin 
1466356e5837SBarry Smith   mname[0] = 0;
14679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1));
1468e5c89e4eSSatish Balay   if (flg1) {
146991eabc43SBarry Smith     PetscViewer viewer;
14709566063dSJacob Faibussowitsch     PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n"));
147191eabc43SBarry Smith     if (mname[0]) {
14729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer));
14739566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
14749566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
147533f85c2fSBarry Smith     } else {
147633f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
14779566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT));
14789566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
14799566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
148033f85c2fSBarry Smith     }
1481e5c89e4eSSatish Balay   }
1482a297a907SKarl Rupp 
1483dd710f27SBarry Smith   /*
1484dd710f27SBarry Smith      Free any objects created by the last block of code.
1485dd710f27SBarry Smith   */
14869566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1487dd710f27SBarry Smith 
1488dd710f27SBarry Smith   mname[0] = 0;
14899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1));
14909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2));
14919566063dSJacob Faibussowitsch   if (flg1 || flg2) PetscCall(PetscLogDump(mname));
1492e5c89e4eSSatish Balay #endif
149310463e74SBarry Smith 
149490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL));
14969566063dSJacob Faibussowitsch   if (!flg1) PetscCall(PetscPopSignalHandler());
149790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL));
14991baa6e33SBarry Smith   if (flg1) PetscCall(PetscMPIDump(stdout));
150090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
150190d69ab7SBarry Smith   flg2 = PETSC_FALSE;
15028bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
15039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1));
15049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1));
15059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL));
1506e4c476e2SSatish Balay 
1507e5c89e4eSSatish Balay   if (flg2) {
1508be56827dSJed Brown     PetscViewer viewer;
15099566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
15109566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer,PETSCVIEWERASCII));
15119566063dSJacob Faibussowitsch     PetscCall(PetscOptionsView(NULL,viewer));
15129566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1513e5c89e4eSSatish Balay   }
1514e5c89e4eSSatish Balay 
1515e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
15169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-nox",&flg1));
15179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1));
1518e5c89e4eSSatish Balay 
151933fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
15209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1));
15219245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1522e5c89e4eSSatish Balay   if (flg3) {
15233de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1524be56827dSJed Brown       PetscViewer viewer;
15259566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
15269566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer,PETSCVIEWERASCII));
15279566063dSJacob Faibussowitsch       PetscCall(PetscOptionsView(NULL,viewer));
15289566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
1529e5c89e4eSSatish Balay     }
15309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsAllUsed(NULL,&nopt));
15313de2bfdfSBarry Smith     if (nopt) {
15329566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n"));
15339566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n"));
15343de2bfdfSBarry Smith       if (nopt == 1) {
15359566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n"));
1536e5c89e4eSSatish Balay       } else {
15379566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"There are %" PetscInt_FMT " unused database options. They are:\n",nopt));
1538e5c89e4eSSatish Balay       }
15393de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15409566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n"));
1541df12ba86SBarry Smith     }
15429566063dSJacob Faibussowitsch     PetscCall(PetscOptionsLeft(NULL));
1543e5c89e4eSSatish Balay   }
1544e5c89e4eSSatish Balay 
1545e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1546d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15479566063dSJacob Faibussowitsch     PetscCall(PetscStackSAWsViewOff());
1548792fecdfSBarry Smith     PetscCallSAWs(SAWs_Finalize,());
1549d45a07a7SBarry Smith   }
1550ec957eceSBarry Smith #endif
1551ec957eceSBarry Smith 
15524097062eSBarry Smith #if defined(PETSC_USE_LOG)
155310463e74SBarry Smith   /*
1554dbc8283eSBarry Smith        List all objects the user may have forgot to free
15552eff7a51SBarry Smith   */
155605df10baSBarry Smith   if (PetscObjectsLog) {
15579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1));
1558a64a8e02SBarry Smith     if (flg1) {
1559a64a8e02SBarry Smith       MPI_Comm local_comm;
15607eb1d149SBarry Smith       char     string[64];
1561a64a8e02SBarry Smith 
15629566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL));
1563afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD,&local_comm));
15649566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm,1));
15659566063dSJacob Faibussowitsch       PetscCall(PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
15669566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm,1));
15679566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
15680a1571b3SBarry Smith     }
156905df10baSBarry Smith   }
15704097062eSBarry Smith #endif
15714097062eSBarry Smith 
15724097062eSBarry Smith #if defined(PETSC_USE_LOG)
1573dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1574dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
15759566063dSJacob Faibussowitsch   PetscCall(PetscFree(PetscObjects));
15764097062eSBarry Smith #endif
15772eff7a51SBarry Smith 
157833f85c2fSBarry Smith   /*
157933f85c2fSBarry Smith      Destroy any packages that registered a finalize
158033f85c2fSBarry Smith   */
15819566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalizeAll());
158233f85c2fSBarry Smith 
1583101409b8SToby Isaac #if defined(PETSC_USE_LOG)
15849566063dSJacob Faibussowitsch   PetscCall(PetscLogFinalize());
1585101409b8SToby Isaac #endif
1586101409b8SToby Isaac 
158733f85c2fSBarry Smith   /*
158848dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
158948dd1dffSBarry Smith   */
15902e956fe4SStefano Zampini   if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll());
159137e93019SBarry Smith 
15924028d114SSatish Balay   if (petsc_history) {
15939566063dSJacob Faibussowitsch     PetscCall(PetscCloseHistoryFile(&petsc_history));
159402c9f0b5SLisandro Dalcin     petsc_history = NULL;
1595e5c89e4eSSatish Balay   }
15969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
15979566063dSJacob Faibussowitsch   PetscCall(PetscInfoDestroy());
1598e5c89e4eSSatish Balay 
159967584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
160092f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1601e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
160292f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1603e5c89e4eSSatish Balay     FILE *fd;
1604ed9cf6e9SBarry Smith     int  err;
1605e5c89e4eSSatish Balay 
1606dc92acbaSJed Brown     flg2 = PETSC_FALSE;
160792f119d6SBarry Smith     flg3 = PETSC_FALSE;
16089566063dSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL));
16099566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL));
161092f119d6SBarry Smith     fname[0] = 0;
16119566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1));
1612e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1613e5c89e4eSSatish Balay 
1614589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
161528b400f6SJacob Faibussowitsch       fd   = fopen(sname,"w"); PetscCheck(fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
16169566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(fd));
1617ed9cf6e9SBarry Smith       err  = fclose(fd);
161828b400f6SJacob Faibussowitsch       PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
161992f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1620e5c89e4eSSatish Balay       MPI_Comm local_comm;
1621e5c89e4eSSatish Balay 
1622afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD,&local_comm));
16239566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm,1));
16249566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(stdout));
16259566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm,1));
16269566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1627e5c89e4eSSatish Balay     }
1628e5c89e4eSSatish Balay     fname[0] = 0;
16299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1));
1630e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1631e5c89e4eSSatish Balay 
1632589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
163328b400f6SJacob Faibussowitsch       fd   = fopen(sname,"w"); PetscCheck(fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
16349566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(fd));
1635ed9cf6e9SBarry Smith       err  = fclose(fd);
163628b400f6SJacob Faibussowitsch       PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
163792f119d6SBarry Smith     } else if (flg1) {
163892f119d6SBarry Smith       MPI_Comm local_comm;
163992f119d6SBarry Smith 
1640afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD,&local_comm));
16419566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm,1));
16429566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(stdout));
16439566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm,1));
16449566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1645e5c89e4eSSatish Balay     }
1646e5c89e4eSSatish Balay   }
164767584ceeSBarry Smith #endif
164820e2c332SMatthew G. Knepley 
16495486ca60SMatthew G. Knepley   /*
16505486ca60SMatthew G. Knepley      Close any open dynamic libraries
16515486ca60SMatthew G. Knepley   */
16529566063dSJacob Faibussowitsch   PetscCall(PetscFinalize_DynamicLibraries());
16535486ca60SMatthew G. Knepley 
1654e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroyDefault());
1656e5c89e4eSSatish Balay 
1657e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
165802c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1659e5c89e4eSSatish Balay 
1660c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1661c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
16629566063dSJacob Faibussowitsch     PetscCall(PetscKokkosFinalize_Private());
1663c2b86a48SJunchao Zhang     PetscBeganKokkos = PETSC_FALSE;
166445639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1665c2b86a48SJunchao Zhang   }
1666c2b86a48SJunchao Zhang #endif
1667c2b86a48SJunchao Zhang 
166871438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
166971438e86SJunchao Zhang   if (PetscBeganNvshmem) {
16709566063dSJacob Faibussowitsch     PetscCall(PetscNvshmemFinalize());
167171438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
167271438e86SJunchao Zhang   }
167371438e86SJunchao Zhang #endif
1674a0e72f99SJunchao Zhang 
16759566063dSJacob Faibussowitsch   PetscCall(PetscFreeMPIResources());
1676e5c89e4eSSatish Balay 
1677dbc8283eSBarry Smith   /*
1678efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1679efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1680efb80d3cSBarry Smith 
1681efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1682efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1683dbc8283eSBarry Smith  */
1684b770b1f6SSatish Balay   {
1685dbc8283eSBarry Smith     PetscCommCounter *counter;
1686dbc8283eSBarry Smith     PetscMPIInt      flg;
1687dbc8283eSBarry Smith     MPI_Comm         icomm;
1688265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
16899566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg));
1690dbc8283eSBarry Smith     if (flg) {
1691265f3f35SJed Brown       icomm = ucomm.comm;
16929566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg));
169328b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1694dbc8283eSBarry Smith 
16959566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval));
16969566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval));
16979566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1698dbc8283eSBarry Smith     }
16999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg));
1700dbc8283eSBarry Smith     if (flg) {
1701265f3f35SJed Brown       icomm = ucomm.comm;
17029566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg));
170328b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1704dbc8283eSBarry Smith 
17059566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval));
17069566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval));
17079566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1708dbc8283eSBarry Smith     }
1709b770b1f6SSatish Balay   }
1710dbc8283eSBarry Smith 
17119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
17129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
17139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
17149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
1715480cf27aSJed Brown 
17169566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
17179566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
17189566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
17199566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock));
1720ef19f930SBarry Smith 
1721e5c89e4eSSatish Balay   if (PetscBeganMPI) {
172299b1327fSBarry Smith     PetscMPIInt flag;
17239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalized(&flag));
172439a651e2SJacob Faibussowitsch     PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
172539a651e2SJacob Faibussowitsch     /* wait until the very last moment to disable error handling */
172639a651e2SJacob Faibussowitsch     PetscErrorHandlingInitialized = PETSC_FALSE;
17279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalize());
172839a651e2SJacob Faibussowitsch   } else PetscErrorHandlingInitialized = PETSC_FALSE;
172939a651e2SJacob Faibussowitsch 
1730e5c89e4eSSatish Balay /*
1731e5c89e4eSSatish Balay 
1732e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1733e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1734e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1735e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1736e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17370e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1738e5c89e4eSSatish Balay    memory was not freed.
1739e5c89e4eSSatish Balay 
1740e5c89e4eSSatish Balay */
17419566063dSJacob Faibussowitsch   PetscCall(PetscMallocClear());
17429566063dSJacob Faibussowitsch   PetscCall(PetscStackReset());
1743a297a907SKarl Rupp 
1744e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1745e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
174656883f60SBarry Smith #if defined(PETSC_USE_GCOV)
174756883f60SBarry Smith   /*
174856883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
174956883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
175056883f60SBarry Smith    */
175156883f60SBarry Smith   __gcov_flush();
175256883f60SBarry Smith #endif
17531724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
17541724198aSStefano Zampini   PetscStackClearTop;
175527104ee2SJacob Faibussowitsch   return 0;
1756e5c89e4eSSatish Balay }
1757e5c89e4eSSatish Balay 
175843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
17598cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
176043db4dbbSBarry Smith {
176143db4dbbSBarry Smith   if (*a == *b) return 1;
176243db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
176343db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
176443db4dbbSBarry Smith   return 0;
176543db4dbbSBarry Smith }
1766a70650f6SBarry Smith #endif
176743db4dbbSBarry Smith 
176843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
17698cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
177043db4dbbSBarry Smith {
177143db4dbbSBarry Smith   if (*a == *b) return 1;
177243db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
177343db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
177443db4dbbSBarry Smith   return 0;
177543db4dbbSBarry Smith }
177643db4dbbSBarry Smith #endif
1777