xref: /petsc/src/sys/objects/pinit.c (revision ad2e3d557eac138370cb75454ee4614fd15ef442)
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>
762e5d2d2SJDBetteridge #include <petsc/private/garbagecollector.h>
8a0e72f99SJunchao Zhang 
96edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
106edef35eSSatish Balay   #include <petsc/private/valgrind/valgrind.h>
116edef35eSSatish Balay #endif
126edef35eSSatish Balay 
13fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS)
1485649d77SJunchao Zhang   #include <petsc/private/fortranimpl.h>
1585649d77SJunchao Zhang #endif
1685649d77SJunchao Zhang 
177ce81a4bSJacob Faibussowitsch #if PetscDefined(USE_COVERAGE)
1856883f60SBarry Smith EXTERN_C_BEGIN
19aaf3846cSSatish Balay   #if defined(PETSC_HAVE___GCOV_DUMP)
20aaf3846cSSatish Balay     #define __gcov_flush(x) __gcov_dump(x)
21aaf3846cSSatish Balay   #endif
2256883f60SBarry Smith void __gcov_flush(void);
2356883f60SBarry Smith EXTERN_C_END
2456883f60SBarry Smith #endif
258101f56cSMatthew Knepley 
262d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
2795c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
282d53ad75SBarry Smith PetscFPT              PetscFPTData = 0;
292d53ad75SBarry Smith #endif
302d53ad75SBarry Smith 
3127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
3216ad0300SBarry Smith   #include <petscviewersaws.h>
33a6790183SBarry Smith #endif
34eb02082bSJunchao Zhang 
3595c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
36e5c89e4eSSatish Balay 
3795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm, int);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm, int);
4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE **);
420069ddf5SShri Abhyankar 
436de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
44e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4527104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
466de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
476de5d289SStefano Zampini #else
486de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
496de5d289SStefano Zampini #endif
50e5c89e4eSSatish Balay 
51480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval      = MPI_KEYVAL_INVALID;
52480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval    = MPI_KEYVAL_INVALID;
53480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval    = MPI_KEYVAL_INVALID;
545f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval      = MPI_KEYVAL_INVALID;
5562e5d2d2SJDBetteridge PetscMPIInt Petsc_CreationIdx_keyval  = MPI_KEYVAL_INVALID;
5662e5d2d2SJDBetteridge PetscMPIInt Petsc_Garbage_HMap_keyval = MPI_KEYVAL_INVALID;
57480cf27aSJed Brown 
585ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedWD_keyval  = MPI_KEYVAL_INVALID;
595ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedTmp_keyval = MPI_KEYVAL_INVALID;
605ea2b939SDuncan Campbell 
61e5c89e4eSSatish Balay /*
62e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
63e5c89e4eSSatish Balay */
6402c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL};
6502c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL};
66e5c89e4eSSatish Balay 
67ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
68ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
690f8e0872SSatish Balay 
70a2f94806SJed Brown PetscInt PetscHotRegionDepth;
71a2f94806SJed Brown 
726edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
736edef35eSSatish Balay 
74b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
75b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
76b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
77b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
78b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
79b22622e2STadeu Manoel #endif
80b22622e2STadeu Manoel 
81aec76313SJacob Faibussowitsch /*@C
82945d1669SBarry Smith   PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8372a42c3cSBarry Smith 
8472a42c3cSBarry Smith   Collective
8572a42c3cSBarry Smith 
8610450e9eSJacob Faibussowitsch   Input Parameters:
8710450e9eSJacob Faibussowitsch + argc     - number of args
8810450e9eSJacob Faibussowitsch . args     - array of command line arguments
8910450e9eSJacob Faibussowitsch . filename - optional name of the program file, pass `NULL` to ignore
9010450e9eSJacob Faibussowitsch - help     - optional help, pass `NULL` to ignore
9110450e9eSJacob Faibussowitsch 
9272a42c3cSBarry Smith   Level: advanced
9372a42c3cSBarry Smith 
9495452b02SPatrick Sanan   Notes:
95a56f64adSBarry Smith   this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
960f11a792SBarry Smith   indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
97a56f64adSBarry Smith   be called multiple times from Julia without the problem of trying to initialize MPI more than once.
980f11a792SBarry Smith 
9910450e9eSJacob Faibussowitsch   Developer Notes:
10010450e9eSJacob Faibussowitsch   Turns off PETSc signal handling to allow Julia to manage signals
1011ea5a559SBarry Smith 
102db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()`
1030f11a792SBarry Smith */
104d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoPointers(int argc, char **args, const char *filename, const char *help)
105d71ae5a4SJacob Faibussowitsch {
10672a42c3cSBarry Smith   int    myargc = argc;
10772a42c3cSBarry Smith   char **myargs = args;
10872a42c3cSBarry Smith 
10972a42c3cSBarry Smith   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&myargc, &myargs, filename, help));
1119566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
112df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11472a42c3cSBarry Smith }
11572a42c3cSBarry Smith 
116f0865b08SBarry Smith /*
117a56f64adSBarry Smith       Used by Julia interface to get communicator
118f0865b08SBarry Smith */
119d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
120d71ae5a4SJacob Faibussowitsch {
121f0865b08SBarry Smith   PetscFunctionBegin;
1224f572ea9SToby Isaac   if (PetscInitializeCalled) PetscAssertPointer(comm, 1);
123f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125f0865b08SBarry Smith }
126f0865b08SBarry Smith 
127e5c89e4eSSatish Balay /*@C
128811af0c4SBarry Smith   PetscInitializeNoArguments - Calls `PetscInitialize()` from C/C++ without
129e5c89e4eSSatish Balay   the command line arguments.
130e5c89e4eSSatish Balay 
131e5c89e4eSSatish Balay   Collective
132e5c89e4eSSatish Balay 
133e5c89e4eSSatish Balay   Level: advanced
134e5c89e4eSSatish Balay 
135db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`
136e5c89e4eSSatish Balay @*/
137d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoArguments(void)
138d71ae5a4SJacob Faibussowitsch {
139e5c89e4eSSatish Balay   int    argc = 0;
14002c9f0b5SLisandro Dalcin   char **args = NULL;
141e5c89e4eSSatish Balay 
142e5c89e4eSSatish Balay   PetscFunctionBegin;
1439566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, NULL, NULL));
1443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145e5c89e4eSSatish Balay }
146e5c89e4eSSatish Balay 
147e5c89e4eSSatish Balay /*@
148e5c89e4eSSatish Balay   PetscInitialized - Determine whether PETSc is initialized.
149e5c89e4eSSatish Balay 
15010450e9eSJacob Faibussowitsch   Output Parameter:
15110450e9eSJacob Faibussowitsch . isInitialized - `PETSC_TRUE` if PETSc is initialized, `PETSC_FALSE` otherwise
15210450e9eSJacob Faibussowitsch 
15393b6d2d1SJed Brown   Level: beginner
154e5c89e4eSSatish Balay 
155db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
156e5c89e4eSSatish Balay @*/
157d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialized(PetscBool *isInitialized)
158d71ae5a4SJacob Faibussowitsch {
15927104ee2SJacob Faibussowitsch   PetscFunctionBegin;
1604f572ea9SToby Isaac   if (PetscInitializeCalled) PetscAssertPointer(isInitialized, 1);
161e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
163e5c89e4eSSatish Balay }
164e5c89e4eSSatish Balay 
165e5c89e4eSSatish Balay /*@
166811af0c4SBarry Smith   PetscFinalized - Determine whether `PetscFinalize()` has been called yet
167e5c89e4eSSatish Balay 
16810450e9eSJacob Faibussowitsch   Output Parameter:
16910450e9eSJacob Faibussowitsch . isFinalized - `PETSC_TRUE` if PETSc is finalized, `PETSC_FALSE` otherwise
17010450e9eSJacob Faibussowitsch 
171e5c89e4eSSatish Balay   Level: developer
172e5c89e4eSSatish Balay 
173db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
174e5c89e4eSSatish Balay @*/
175d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalized(PetscBool *isFinalized)
176d71ae5a4SJacob Faibussowitsch {
17727104ee2SJacob Faibussowitsch   PetscFunctionBegin;
1784f572ea9SToby Isaac   if (!PetscFinalizeCalled) PetscAssertPointer(isFinalized, 1);
179e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
181e5c89e4eSSatish Balay }
182e5c89e4eSSatish Balay 
18357171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]);
184e5c89e4eSSatish Balay 
185e5c89e4eSSatish Balay /*
186e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
187e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
188e5c89e4eSSatish Balay */
189367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP               = 0;
19062e5d2d2SJDBetteridge MPI_Op Petsc_Garbage_SetIntersectOp = 0;
191e5c89e4eSSatish Balay 
192d71ae5a4SJacob Faibussowitsch PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, int *cnt, MPI_Datatype *datatype)
193d71ae5a4SJacob Faibussowitsch {
194e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out, i, count = *cnt;
195e5c89e4eSSatish Balay 
196e5c89e4eSSatish Balay   PetscFunctionBegin;
197e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
1983ba16761SJacob Faibussowitsch     PetscErrorCode ierr = (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
1993ba16761SJacob Faibussowitsch     (void)ierr;
20041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
201e5c89e4eSSatish Balay   }
202e5c89e4eSSatish Balay 
203e5c89e4eSSatish Balay   for (i = 0; i < count; i++) {
204e5c89e4eSSatish Balay     xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]);
205e5c89e4eSSatish Balay     xout[2 * i + 1] += xin[2 * i + 1];
206e5c89e4eSSatish Balay   }
207812af9f3SBarry Smith   PetscFunctionReturnVoid();
208e5c89e4eSSatish Balay }
209e5c89e4eSSatish Balay 
210e5c89e4eSSatish Balay /*
211e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
212e5c89e4eSSatish Balay sum of the second entry.
213b693b147SBarry Smith 
21476f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
215367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
216b693b147SBarry Smith there would be no place to store the both needed results.
217e5c89e4eSSatish Balay */
218d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt sizes[], PetscInt *max, PetscInt *sum)
219d71ae5a4SJacob Faibussowitsch {
220e5c89e4eSSatish Balay   PetscFunctionBegin;
221d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
222d6e4c47cSJed Brown   {
2239371c9d4SSatish Balay     struct {
2249371c9d4SSatish Balay       PetscInt max, sum;
2259371c9d4SSatish Balay     } work;
2269566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce_scatter_block((void *)sizes, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm));
227d6e4c47cSJed Brown     *max = work.max;
228d6e4c47cSJed Brown     *sum = work.sum;
229d6e4c47cSJed Brown   }
230d6e4c47cSJed Brown #else
231d6e4c47cSJed Brown   {
232d6e4c47cSJed Brown     PetscMPIInt size, rank;
2339371c9d4SSatish Balay     struct {
2349371c9d4SSatish Balay       PetscInt max, sum;
2359371c9d4SSatish Balay     } *work;
2369566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
2379566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
2389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &work));
2391c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce((void *)sizes, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm));
2406ac3741eSJed Brown     *max = work[rank].max;
2416ac3741eSJed Brown     *sum = work[rank].sum;
2429566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
243d6e4c47cSJed Brown   }
244d6e4c47cSJed Brown #endif
2453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
246e5c89e4eSSatish Balay }
247e5c89e4eSSatish Balay 
2489e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
249613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
250613bf2b2SPierre Jolivet     #include <quadmath.h>
251613bf2b2SPierre Jolivet   #endif
2529e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0;
253de272c7aSSatish Balay   #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
25406a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
255613bf2b2SPierre Jolivet   #endif
256e5c89e4eSSatish Balay 
257d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
258d71ae5a4SJacob Faibussowitsch {
259e5c89e4eSSatish Balay   PetscInt i, count = *cnt;
260e5c89e4eSSatish Balay 
261e5c89e4eSSatish Balay   PetscFunctionBegin;
2627c2de775SJed Brown   if (*datatype == MPIU_REAL) {
263e2e03761SBarry Smith     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
264a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2657c2de775SJed Brown   }
2667c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
267c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2687c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
269a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2707c2de775SJed Brown   }
2717c2de775SJed Brown   #endif
272613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
273613bf2b2SPierre Jolivet   else if (*datatype == MPIU___FLOAT128) {
274613bf2b2SPierre Jolivet     __float128 *xin = (__float128 *)in, *xout = (__float128 *)out;
275613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
27674c01ffaSSatish Balay     #if defined(PETSC_HAVE_COMPLEX)
2779371c9d4SSatish Balay   } else if (*datatype == MPIU___COMPLEX128) {
278613bf2b2SPierre Jolivet     __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out;
279613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
28074c01ffaSSatish Balay     #endif
281613bf2b2SPierre Jolivet   }
282613bf2b2SPierre Jolivet   #endif
2839e517322SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FP16)
2849e517322SPierre Jolivet   else if (*datatype == MPIU___FP16) {
2859e517322SPierre Jolivet     __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out;
2869e517322SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2879e517322SPierre Jolivet   }
2889e517322SPierre Jolivet   #endif
2897c2de775SJed Brown   else {
2909e517322SPierre Jolivet   #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16)
2913ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SElF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
2929e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FP16)
2933ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types"));
2949e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FLOAT128)
2953ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types"));
2969e517322SPierre Jolivet   #else
2973ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types"));
298613bf2b2SPierre Jolivet   #endif
29941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
300e2e03761SBarry Smith   }
301812af9f3SBarry Smith   PetscFunctionReturnVoid();
302e5c89e4eSSatish Balay }
303e5c89e4eSSatish Balay #endif
304e5c89e4eSSatish Balay 
305570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
306d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
307d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
308d9822059SBarry Smith 
309d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
310d71ae5a4SJacob Faibussowitsch {
311d9822059SBarry Smith   PetscInt i, count = *cnt;
312d9822059SBarry Smith 
313d9822059SBarry Smith   PetscFunctionBegin;
3147c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3158c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
316a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]);
3177c2de775SJed Brown   }
3187c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3197c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3207c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
321ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3227c2de775SJed Brown   }
3237c2de775SJed Brown   #endif
3247c2de775SJed Brown   else {
3253ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
32641e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3278c764dc5SJose Roman   }
328d9822059SBarry Smith   PetscFunctionReturnVoid();
329d9822059SBarry Smith }
330d9822059SBarry Smith 
331d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
332d71ae5a4SJacob Faibussowitsch {
333d9822059SBarry Smith   PetscInt i, count = *cnt;
334d9822059SBarry Smith 
335d9822059SBarry Smith   PetscFunctionBegin;
3367c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3378c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
338a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]);
3397c2de775SJed Brown   }
3407c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3417c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3427c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
343ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3447c2de775SJed Brown   }
3457c2de775SJed Brown   #endif
3467c2de775SJed Brown   else {
3473ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"));
34841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3498c764dc5SJose Roman   }
350d9822059SBarry Smith   PetscFunctionReturnVoid();
351d9822059SBarry Smith }
352d9822059SBarry Smith #endif
353d9822059SBarry Smith 
354480cf27aSJed Brown /*
355480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
356480cf27aSJed Brown 
357ff0e51ddSBarry 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.
358480cf27aSJed Brown 
35912801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
360480cf27aSJed Brown 
361480cf27aSJed Brown */
362d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state)
363d71ae5a4SJacob Faibussowitsch {
36405342407SJunchao Zhang   PetscCommCounter      *counter = (PetscCommCounter *)count_val;
36557f21012SBarry Smith   struct PetscCommStash *comms   = counter->comms, *pcomm;
366480cf27aSJed Brown 
367480cf27aSJed Brown   PetscFunctionBegin;
3689566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm));
3699566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter->iflags));
37057f21012SBarry Smith   while (comms) {
3719566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&comms->comm));
37257f21012SBarry Smith     pcomm = comms;
37357f21012SBarry Smith     comms = comms->next;
3749566063dSJacob Faibussowitsch     PetscCall(PetscFree(pcomm));
37557f21012SBarry Smith   }
3769566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter));
377480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
378480cf27aSJed Brown }
379480cf27aSJed Brown 
380480cf27aSJed Brown /*
38147435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
382da3039f7SJed Brown   calls MPI_Comm_free().
383da3039f7SJed Brown 
384da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
385480cf27aSJed Brown 
386ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
387480cf27aSJed Brown 
38812801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
389480cf27aSJed Brown 
390480cf27aSJed Brown */
391d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
392d71ae5a4SJacob Faibussowitsch {
3939371c9d4SSatish Balay   union
394480cf27aSJed Brown   {
3959371c9d4SSatish Balay     MPI_Comm comm;
3969371c9d4SSatish Balay     void    *ptr;
3979371c9d4SSatish Balay   } icomm;
398480cf27aSJed Brown 
399480cf27aSJed Brown   PetscFunctionBegin;
40012801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval");
401ec4fadc2SJed Brown   icomm.ptr = attr_val;
40276bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
40333779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
40433779a13SJunchao Zhang     PetscMPIInt flg;
4059371c9d4SSatish Balay     union
4069371c9d4SSatish Balay     {
4079371c9d4SSatish Balay       MPI_Comm comm;
4089371c9d4SSatish Balay       void    *ptr;
4099371c9d4SSatish Balay     } ocomm;
4109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg));
41133779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute");
41233779a13SJunchao 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");
41333779a13SJunchao Zhang   }
4149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval));
4159566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm));
416da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
417b89831e5SBarry Smith }
418da3039f7SJed Brown 
419da3039f7SJed Brown /*
42033779a13SJunchao 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.
421da3039f7SJed Brown  */
422d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
423d71ae5a4SJacob Faibussowitsch {
424da3039f7SJed Brown   PetscFunctionBegin;
4259566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm));
426480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
427480cf27aSJed Brown }
428480cf27aSJed Brown 
42933779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *);
43042218b76SBarry Smith 
431951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
4328cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *);
4338cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
4348cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
435e39fd77fSBarry Smith #endif
436e39fd77fSBarry Smith 
4370f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE;
43812801b39SBarry Smith 
439eb27c7c8SSatish Balay PETSC_INTERN int    PetscGlobalArgc;
440eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4416ae9a8a6SBarry Smith int                 PetscGlobalArgc = 0;
44202c9f0b5SLisandro Dalcin char              **PetscGlobalArgs = NULL;
443dff31646SBarry Smith PetscSegBuffer      PetscCitationsList;
444e5c89e4eSSatish Balay 
445d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void)
446d71ae5a4SJacob Faibussowitsch {
447051e4cf2SJed Brown   PetscFunctionBegin;
4489566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList));
44930c35bf2SSatish Balay 
45030c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\
45130c35bf2SSatish Balay   Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\
45230c35bf2SSatish Balay     and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\
4533f81df79SBarry Smith     and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\
45430c35bf2SSatish Balay     and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\
45530c35bf2SSatish Balay     and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\
45630c35bf2SSatish Balay     and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\
45730c35bf2SSatish Balay     and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\
45830c35bf2SSatish Balay   Title = {{PETSc/TAO} Users Manual},\n\
459c4dc7257SSatish Balay   Number = {ANL-21/39 - Revision 3.19},\n\
4609a2cd68eSMatthew Knepley   Doi = {10.2172/1968587},\n\
46130c35bf2SSatish Balay   Institution = {Argonne National Laboratory},\n\
462c4dc7257SSatish Balay   Year = {2023}\n}\n",
4639371c9d4SSatish Balay                                    NULL));
46430c35bf2SSatish Balay 
46530c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\
46630c35bf2SSatish Balay   Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\
46730c35bf2SSatish Balay   Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\
46830c35bf2SSatish Balay   Booktitle = {Modern Software Tools in Scientific Computing},\n\
46930c35bf2SSatish Balay   Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\
47030c35bf2SSatish Balay   Pages = {163--202},\n\
47130c35bf2SSatish Balay   Publisher = {Birkh{\\\"{a}}user Press},\n\
4729371c9d4SSatish Balay   Year = {1997}\n}\n",
4739371c9d4SSatish Balay                                    NULL));
47430c35bf2SSatish Balay 
4753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
476051e4cf2SJed Brown }
477e5c89e4eSSatish Balay 
4782d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4792d747510SLisandro Dalcin 
480d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[])
481d71ae5a4SJacob Faibussowitsch {
4822d747510SLisandro Dalcin   PetscFunctionBegin;
4839566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(programname, name, sizeof(programname)));
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4852d747510SLisandro Dalcin }
4862d747510SLisandro Dalcin 
4872d747510SLisandro Dalcin /*@C
4882d747510SLisandro Dalcin   PetscGetProgramName - Gets the name of the running program.
4892d747510SLisandro Dalcin 
4902d747510SLisandro Dalcin   Not Collective
4912d747510SLisandro Dalcin 
4922d747510SLisandro Dalcin   Input Parameter:
4932d747510SLisandro Dalcin . len - length of the string name
4942d747510SLisandro Dalcin 
4952d747510SLisandro Dalcin   Output Parameter:
496811af0c4SBarry Smith . name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN`
4972d747510SLisandro Dalcin 
4982d747510SLisandro Dalcin   Level: advanced
4992d747510SLisandro Dalcin 
50021532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()`
5012d747510SLisandro Dalcin @*/
502d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProgramName(char name[], size_t len)
503d71ae5a4SJacob Faibussowitsch {
5042d747510SLisandro Dalcin   PetscFunctionBegin;
5059566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(name, programname, len));
5063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5072d747510SLisandro Dalcin }
5082d747510SLisandro Dalcin 
509e5c89e4eSSatish Balay /*@C
510e5c89e4eSSatish Balay   PetscGetArgs - Allows you to access the raw command line arguments anywhere
511811af0c4SBarry Smith   after PetscInitialize() is called but before `PetscFinalize()`.
512e5c89e4eSSatish Balay 
513e5c89e4eSSatish Balay   Not Collective
514e5c89e4eSSatish Balay 
515e5c89e4eSSatish Balay   Output Parameters:
516e5c89e4eSSatish Balay + argc - count of number of command line arguments
517e5c89e4eSSatish Balay - args - the command line arguments
518e5c89e4eSSatish Balay 
519e5c89e4eSSatish Balay   Level: intermediate
520e5c89e4eSSatish Balay 
521e5c89e4eSSatish Balay   Notes:
522e5c89e4eSSatish Balay   This is usually used to pass the command line arguments into other libraries
523e5c89e4eSSatish Balay   that are called internally deep in PETSc or the application.
524e5c89e4eSSatish Balay 
525f177e3b1SBarry Smith   The first argument contains the program name as is normal for C arguments.
526f177e3b1SBarry Smith 
52721532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()`
528e5c89e4eSSatish Balay @*/
529d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args)
530d71ae5a4SJacob Faibussowitsch {
531e5c89e4eSSatish Balay   PetscFunctionBegin;
53239a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
533e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
534e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
5353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
536e5c89e4eSSatish Balay }
537e5c89e4eSSatish Balay 
538793721a6SBarry Smith /*@C
539793721a6SBarry Smith   PetscGetArguments - Allows you to access the  command line arguments anywhere
540811af0c4SBarry Smith   after `PetscInitialize()` is called but before `PetscFinalize()`.
541793721a6SBarry Smith 
542793721a6SBarry Smith   Not Collective
543793721a6SBarry Smith 
5442fe279fdSBarry Smith   Output Parameter:
545793721a6SBarry Smith . args - the command line arguments
546793721a6SBarry Smith 
547793721a6SBarry Smith   Level: intermediate
548793721a6SBarry Smith 
54921532e8aSBarry Smith   Note:
55021532e8aSBarry Smith   This does NOT start with the program name and IS `NULL` terminated (final arg is void)
551793721a6SBarry Smith 
55221532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`, `PetscInitialize()`
553793721a6SBarry Smith @*/
554d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args)
555d71ae5a4SJacob Faibussowitsch {
556793721a6SBarry Smith   PetscInt i, argc = PetscGlobalArgc;
557793721a6SBarry Smith 
558793721a6SBarry Smith   PetscFunctionBegin;
55939a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
5609371c9d4SSatish Balay   if (!argc) {
5619371c9d4SSatish Balay     *args = NULL;
5623ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5639371c9d4SSatish Balay   }
5649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(argc, args));
5659566063dSJacob Faibussowitsch   for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i]));
5662d747510SLisandro Dalcin   (*args)[argc - 1] = NULL;
5673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
568793721a6SBarry Smith }
569793721a6SBarry Smith 
570793721a6SBarry Smith /*@C
571811af0c4SBarry Smith   PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()`
572793721a6SBarry Smith 
573793721a6SBarry Smith   Not Collective
574793721a6SBarry Smith 
5752fe279fdSBarry Smith   Output Parameter:
576793721a6SBarry Smith . args - the command line arguments
577793721a6SBarry Smith 
578793721a6SBarry Smith   Level: intermediate
579793721a6SBarry Smith 
580db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()`
581793721a6SBarry Smith @*/
582d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args)
583d71ae5a4SJacob Faibussowitsch {
58439a651e2SJacob Faibussowitsch   PetscFunctionBegin;
58539a651e2SJacob Faibussowitsch   if (args) {
586793721a6SBarry Smith     PetscInt i = 0;
587793721a6SBarry Smith 
5889566063dSJacob Faibussowitsch     while (args[i]) PetscCall(PetscFree(args[i++]));
5899566063dSJacob Faibussowitsch     PetscCall(PetscFree(args));
59039a651e2SJacob Faibussowitsch   }
5913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
592793721a6SBarry Smith }
593793721a6SBarry Smith 
59427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
59530befbd2SBarry Smith   #include <petscconfiginfo.h>
59630befbd2SBarry Smith 
597d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
598d71ae5a4SJacob Faibussowitsch {
59927104ee2SJacob Faibussowitsch   PetscFunctionBegin;
60011525c0dSBarry Smith   if (!PetscGlobalRank) {
60130befbd2SBarry Smith     char      cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64];
60211525c0dSBarry Smith     int       port;
603ffbd1cfbSBarry Smith     PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE;
60411525c0dSBarry Smith     size_t    applinelen, introlen;
605ffbd1cfbSBarry Smith     char      sawsurl[256];
60611525c0dSBarry Smith 
6079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg));
60811525c0dSBarry Smith     if (flg) {
60911525c0dSBarry Smith       char sawslog[PETSC_MAX_PATH_LEN];
61011525c0dSBarry Smith 
6119566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL));
61211525c0dSBarry Smith       if (sawslog[0]) {
613792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog));
61411525c0dSBarry Smith       } else {
615792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL));
61611525c0dSBarry Smith       }
61711525c0dSBarry Smith     }
6189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg));
61948a46eb9SPierre Jolivet     if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert));
6209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL));
621ffbd1cfbSBarry Smith     if (selectport) {
622792fecdfSBarry Smith       PetscCallSAWs(SAWs_Get_Available_Port, (&port));
623792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Port, (port));
624ffbd1cfbSBarry Smith     } else {
6259566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg));
62648a46eb9SPierre Jolivet       if (flg) PetscCallSAWs(SAWs_Set_Port, (port));
627ffbd1cfbSBarry Smith     }
6289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg));
62911525c0dSBarry Smith     if (flg) {
630792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Document_Root, (root));
6319566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(root, ".", &rootlocal));
6329c1e0ce8SBarry Smith     } else {
6339566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg));
6349c1e0ce8SBarry Smith       if (flg) {
6359566063dSJacob Faibussowitsch         PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root)));
636792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Document_Root, (root));
6379c1e0ce8SBarry Smith       }
63811525c0dSBarry Smith     }
6399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2));
64011525c0dSBarry Smith     if (flg2) {
64111525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
64228b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option");
6439566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root));
6449566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(jsdir, 'r', &flg));
64528b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory");
646792fecdfSBarry Smith       PetscCallSAWs(SAWs_Push_Local_Header, ());
64711525c0dSBarry Smith     }
6489566063dSJacob Faibussowitsch     PetscCall(PetscGetProgramName(programname, sizeof(programname)));
6499566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(help, &applinelen));
65011525c0dSBarry Smith     introlen = 4096 + applinelen;
65130a8c9c0SSurtai Han     applinelen += 1024;
6529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(applinelen, &appline));
6539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(introlen, &intro));
65411525c0dSBarry Smith 
65511525c0dSBarry Smith     if (rootlocal) {
6569566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname));
6579566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(appline, 'r', &rootlocal));
65811525c0dSBarry Smith     }
6599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetAll(NULL, &options));
66011525c0dSBarry Smith     if (rootlocal && help) {
6619566063dSJacob 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));
66211525c0dSBarry Smith     } else if (help) {
6639566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help));
66411525c0dSBarry Smith     } else {
6659566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options));
66611525c0dSBarry Smith     }
6679566063dSJacob Faibussowitsch     PetscCall(PetscFree(options));
6689566063dSJacob Faibussowitsch     PetscCall(PetscGetVersion(version, sizeof(version)));
6699371c9d4SSatish Balay     PetscCall(PetscSNPrintf(intro, introlen,
6709371c9d4SSatish Balay                             "<body>\n"
671a17b96a8SKyle 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"
672df62c222SBarry 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"
6739371c9d4SSatish Balay                             "%s",
6749371c9d4SSatish Balay                             version, petscconfigureoptions, appline));
675792fecdfSBarry Smith     PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro));
6769566063dSJacob Faibussowitsch     PetscCall(PetscFree(intro));
6779566063dSJacob Faibussowitsch     PetscCall(PetscFree(appline));
678ffbd1cfbSBarry Smith     if (selectport) {
679aa573868SBarry Smith       PetscBool silent;
6807d812c46SBarry Smith 
6817d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
68239a651e2SJacob Faibussowitsch       while (SAWs_Initialize()) {
683792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_Available_Port, (&port));
684792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port, (port));
6857d812c46SBarry Smith       }
6867d812c46SBarry Smith 
6879566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL));
688aa573868SBarry Smith       if (!silent) {
689792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl));
6909566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl));
691ffbd1cfbSBarry Smith       }
6927d812c46SBarry Smith     } else {
693792fecdfSBarry Smith       PetscCallSAWs(SAWs_Initialize, ());
694aa573868SBarry Smith     }
6959566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister("@TechReport{ saws,\n"
6960af79b04SBarry Smith                                      "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6970af79b04SBarry Smith                                      "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6980af79b04SBarry Smith                                      "  Institution = {Argonne National Laboratory},\n"
6999371c9d4SSatish Balay                                      "  Year   = 2013\n}\n",
7009371c9d4SSatish Balay                                      NULL));
70111525c0dSBarry Smith   }
7023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70311525c0dSBarry Smith }
70411525c0dSBarry Smith #endif
70511525c0dSBarry Smith 
7064dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
707d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
708d71ae5a4SJacob Faibussowitsch {
7094dfee713SSatish Balay   PetscFunctionBegin;
7104dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
7114dfee713SSatish Balay   /* see MPI.py for details on this bug */
7124dfee713SSatish Balay   (void)setenv("HWLOC_COMPONENTS", "-x86", 1);
7134dfee713SSatish Balay #endif
7143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7154dfee713SSatish Balay }
7164dfee713SSatish Balay 
717a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
718a56f64adSBarry Smith   #include <adios.h>
71922580e64SBarry Smith   #include <adios_read.h>
7207b56e58cSSatish Balay int64_t Petsc_adios_group;
721a56f64adSBarry Smith #endif
722a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
723cd1b99a6SBarry Smith   #include <omp.h>
724f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
725cd1b99a6SBarry Smith #endif
726a56f64adSBarry Smith 
727a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
728a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA)
7290e6b6b59SJacob Faibussowitsch   #include <petscdevice_cuda.h>
730a4af0ceeSJacob Faibussowitsch // REMOVE ME
731a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
732a4af0ceeSJacob Faibussowitsch #endif
733a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP)
7340e6b6b59SJacob Faibussowitsch   #include <petscdevice_hip.h>
735a4af0ceeSJacob Faibussowitsch // REMOVE ME
736a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
737a4af0ceeSJacob Faibussowitsch #endif
738a4af0ceeSJacob Faibussowitsch 
73927104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
740bc8a24ecSBarry Smith   #include <dlfcn.h>
741bc8a24ecSBarry Smith #endif
742a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
743a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
7443274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscViennaCLInit(void);
745a4af0ceeSJacob Faibussowitsch PetscBool                   PetscViennaCLSynchronize = PETSC_FALSE;
746a4af0ceeSJacob Faibussowitsch #endif
747bc8a24ecSBarry Smith 
748660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE;
749660278c0SBarry Smith 
75085649d77SJunchao Zhang /*
75185649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
75285649d77SJunchao Zhang 
75385649d77SJunchao Zhang   prog:     program name
75402101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
75585649d77SJunchao Zhang   help:     program help message
756da81f932SPierre Jolivet   ftn:      is it called from Fortran initialization (petscinitializef_)?
75785649d77SJunchao Zhang   readarguments,len: used when fortran is true
75885649d77SJunchao Zhang */
759d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len)
760d71ae5a4SJacob Faibussowitsch {
76185649d77SJunchao Zhang   PetscMPIInt size;
76285649d77SJunchao Zhang   PetscBool   flg = PETSC_TRUE;
76385649d77SJunchao Zhang   char        hostname[256];
76485649d77SJunchao Zhang 
76527104ee2SJacob Faibussowitsch   PetscFunctionBegin;
7663ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
76739a651e2SJacob Faibussowitsch   /* these must be initialized in a routine, not as a constant declaration */
76839a651e2SJacob Faibussowitsch   PETSC_STDOUT = stdout;
76939a651e2SJacob Faibussowitsch   PETSC_STDERR = stderr;
77039a651e2SJacob Faibussowitsch 
7719566063dSJacob Faibussowitsch   /* PetscCall can be used from now */
77239a651e2SJacob Faibussowitsch   PetscErrorHandlingInitialized = PETSC_TRUE;
77339a651e2SJacob Faibussowitsch 
77485649d77SJunchao Zhang   /*
77585649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
77685649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
77785649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
77885649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
77985649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
78085649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
78185649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
78285649d77SJunchao Zhang       listed above and since that time are compatible.
78385649d77SJunchao Zhang 
78485649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
78585649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
78685649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
78785649d77SJunchao Zhang       to perform the checking.
78885649d77SJunchao Zhang 
78985649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
79085649d77SJunchao Zhang 
79185649d77SJunchao Zhang       Questions:
79285649d77SJunchao Zhang 
79385649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
79485649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
79585649d77SJunchao Zhang   */
79685649d77SJunchao Zhang 
79785649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
79885649d77SJunchao Zhang   {
79985649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
80085649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
80139a651e2SJacob Faibussowitsch 
8029566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength));
80385649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
80485649d77SJunchao Zhang   #if defined(MPICH_VERSION)
80585649d77SJunchao Zhang     #if MPICH_NUMVERSION < 30100000
80685649d77SJunchao Zhang     {
80785649d77SJunchao Zhang       char     *ver, *lf;
80885649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
80939a651e2SJacob Faibussowitsch 
8109566063dSJacob Faibussowitsch       PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver));
81139a651e2SJacob Faibussowitsch       if (ver) {
8129566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(ver, '\n', &lf));
81339a651e2SJacob Faibussowitsch         if (lf) {
81485649d77SJunchao Zhang           *lf = 0;
8159566063dSJacob Faibussowitsch           PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg));
81685649d77SJunchao Zhang         }
81785649d77SJunchao Zhang       }
81885649d77SJunchao Zhang       if (!flg) {
819d5b396e8SSatish Balay         PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION));
82085649d77SJunchao Zhang         flg = PETSC_TRUE;
82185649d77SJunchao Zhang       }
82285649d77SJunchao Zhang     }
82385649d77SJunchao Zhang     #endif
82485649d77SJunchao 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?) */
82585649d77SJunchao Zhang   #elif defined(OMPI_MAJOR_VERSION)
82685649d77SJunchao Zhang     {
82785649d77SJunchao Zhang       char     *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf;
82885649d77SJunchao Zhang       PetscBool flg                                              = PETSC_FALSE;
82985649d77SJunchao Zhang     #define PSTRSZ 2
83085649d77SJunchao Zhang       char      ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"};
83185649d77SJunchao Zhang       char      ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "};
83285649d77SJunchao Zhang       int       i;
83385649d77SJunchao Zhang       for (i = 0; i < PSTRSZ; i++) {
8349566063dSJacob Faibussowitsch         PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver));
83539a651e2SJacob Faibussowitsch         if (ver) {
8369566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION));
8379566063dSJacob Faibussowitsch           PetscCall(PetscStrstr(ver, bs, &bsf));
83839a651e2SJacob Faibussowitsch           if (bsf) flg = PETSC_TRUE;
83985649d77SJunchao Zhang           break;
84085649d77SJunchao Zhang         }
84185649d77SJunchao Zhang       }
84285649d77SJunchao Zhang       if (!flg) {
8433ba16761SJacob Faibussowitsch         PetscCall(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));
84485649d77SJunchao Zhang         flg = PETSC_TRUE;
84585649d77SJunchao Zhang       }
84685649d77SJunchao Zhang     }
84785649d77SJunchao Zhang   #endif
84885649d77SJunchao Zhang   }
84985649d77SJunchao Zhang #endif
85085649d77SJunchao Zhang 
851e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__))
85285649d77SJunchao 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 */
85339a651e2SJacob 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");
85485649d77SJunchao Zhang #endif
85585649d77SJunchao Zhang 
85685649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
85785649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
85885649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
85985649d77SJunchao Zhang #endif
86085649d77SJunchao Zhang 
8619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
86285649d77SJunchao Zhang 
86385649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
86485649d77SJunchao Zhang 
8659566063dSJacob Faibussowitsch   PetscCall(PetscSetProgramName(prog));
8669566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen));
8679566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout));
8689566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr));
8699566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscCommSpinLock));
87085649d77SJunchao Zhang 
87185649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
8729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN));
87385649d77SJunchao Zhang 
87485649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
8759566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS));
8769566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE));
87785649d77SJunchao Zhang   }
87885649d77SJunchao Zhang 
87985649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
8809566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
88185649d77SJunchao Zhang 
8829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank));
8839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize));
88485649d77SJunchao Zhang 
88585649d77SJunchao Zhang   MPIU_BOOL        = MPI_INT;
88685649d77SJunchao Zhang   MPIU_ENUM        = MPI_INT;
88785649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64;
88885649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
88985649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
89085649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
89185649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
89285649d77SJunchao Zhang #endif
89339a651e2SJacob Faibussowitsch   else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t");
89485649d77SJunchao Zhang 
89585649d77SJunchao Zhang     /*
89685649d77SJunchao Zhang      Initialized the global complex variable; this is because with
89785649d77SJunchao Zhang      shared libraries the constructors for global variables
89885649d77SJunchao Zhang      are not called; at least on IRIX.
89985649d77SJunchao Zhang   */
90085649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
90185649d77SJunchao Zhang   {
90285649d77SJunchao Zhang   #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
90385649d77SJunchao Zhang     PetscComplex ic(0.0, 1.0);
90485649d77SJunchao Zhang     PETSC_i = ic;
90585649d77SJunchao Zhang   #else
90685649d77SJunchao Zhang     PETSC_i = _Complex_I;
90785649d77SJunchao Zhang   #endif
90885649d77SJunchao Zhang   }
90985649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
91085649d77SJunchao Zhang 
91185649d77SJunchao Zhang   /*
91285649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
91385649d77SJunchao Zhang      half of the entries and maxes the second half.
91485649d77SJunchao Zhang   */
9159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP));
91685649d77SJunchao Zhang 
917613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
9189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128));
9199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128));
92074c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
9219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128));
9229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128));
92385649d77SJunchao Zhang   #endif
92474c01ffaSSatish Balay #endif
9259e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
9269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16));
9279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FP16));
92885649d77SJunchao Zhang #endif
92985649d77SJunchao Zhang 
93085649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
9319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM));
932613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX));
933613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN));
9349e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
9359e517322SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128));
93685649d77SJunchao Zhang #endif
93785649d77SJunchao Zhang 
9389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR));
93962e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp));
9409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR));
94185649d77SJunchao Zhang 
94285649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
94385649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
94485649d77SJunchao Zhang   {
94585649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
94693d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)};
94785649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_REAL, MPIU_INT}, tmpStruct;
94885649d77SJunchao Zhang 
9499566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
95093d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT));
9519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9529566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT));
95385649d77SJunchao Zhang   }
95485649d77SJunchao Zhang   {
95585649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
95693d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)};
95785649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_SCALAR, MPIU_INT}, tmpStruct;
95885649d77SJunchao Zhang 
9599566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
96093d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT));
9619566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9629566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
96385649d77SJunchao Zhang   }
96485649d77SJunchao Zhang #endif
96585649d77SJunchao Zhang 
96685649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT));
9689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2INT));
96985649d77SJunchao Zhang #endif
9709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT));
9719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPI_4INT));
9729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT));
9739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_4INT));
97485649d77SJunchao Zhang 
97585649d77SJunchao Zhang   /*
97685649d77SJunchao Zhang      Attributes to be set on PETSc communicators
97785649d77SJunchao Zhang   */
9789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0));
9799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0));
9809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0));
9819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0));
98262e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, (void *)0));
98362e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, (void *)0));
98485649d77SJunchao Zhang 
985fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS)
9869566063dSJacob Faibussowitsch   if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len));
98785649d77SJunchao Zhang   else
98885649d77SJunchao Zhang #endif
9899566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file));
99085649d77SJunchao Zhang 
99185649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9929566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
99385649d77SJunchao Zhang 
99485649d77SJunchao Zhang   /*
99585649d77SJunchao Zhang      Check system options and print help
99685649d77SJunchao Zhang   */
9979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCheckInitial_Private(help));
99885649d77SJunchao Zhang 
999a4af0ceeSJacob Faibussowitsch   /*
10000e6b6b59SJacob Faibussowitsch     Creates the logging data structures; this is enabled even if logging is not turned on
10010e6b6b59SJacob Faibussowitsch     This is the last thing we do before returning to the user code to prevent having the
10020e6b6b59SJacob Faibussowitsch     logging numbers contaminated by any startup time associated with MPI
10030e6b6b59SJacob Faibussowitsch   */
10040e6b6b59SJacob Faibussowitsch   PetscCall(PetscLogInitialize());
10050e6b6b59SJacob Faibussowitsch 
10060e6b6b59SJacob Faibussowitsch   /*
1007a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
1008a4af0ceeSJacob Faibussowitsch 
1009a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
1010a4af0ceeSJacob Faibussowitsch    1. MPI initialized
1011a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
10120e6b6b59SJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up
10130e6b6b59SJacob Faibussowitsch       its own SIGSEV handler via the push/pop interface.
10140e6b6b59SJacob Faibussowitsch    4. Logging initialized
1015a4af0ceeSJacob Faibussowitsch   */
10169566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
1017a4af0ceeSJacob Faibussowitsch 
1018a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
1019a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
1020d75802c7SJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg));
10219566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL));
1022a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
10239566063dSJacob Faibussowitsch   PetscCall(PetscViennaCLInit());
1024a4af0ceeSJacob Faibussowitsch #endif
1025a4af0ceeSJacob Faibussowitsch 
10269566063dSJacob Faibussowitsch   PetscCall(PetscCitationsInitialize());
102785649d77SJunchao Zhang 
102885649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
10299566063dSJacob Faibussowitsch   PetscCall(PetscInitializeSAWs(ftn ? NULL : help));
103027104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
10319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg));
10329566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscStackViewSAWs());
103385649d77SJunchao Zhang #endif
103485649d77SJunchao Zhang 
103585649d77SJunchao Zhang   /*
103685649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
103785649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
103885649d77SJunchao Zhang   */
10399566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_DynamicLibraries());
104085649d77SJunchao Zhang 
10419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10429566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size));
104396dcfd6cSBarry Smith   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
10449566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname));
104585649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
104685649d77SJunchao Zhang   {
104785649d77SJunchao Zhang     PetscBool omp_view_flag;
104885649d77SJunchao Zhang     char     *threads = getenv("OMP_NUM_THREADS");
104985649d77SJunchao Zhang 
105085649d77SJunchao Zhang     if (threads) {
10519566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads));
105285649d77SJunchao Zhang       (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads);
105385649d77SJunchao Zhang     } else {
10542f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt)omp_get_max_threads();
10559566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads));
105685649d77SJunchao Zhang     }
1057d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys");
10589566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg));
10599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag));
1060d0609cedSBarry Smith     PetscOptionsEnd();
106185649d77SJunchao Zhang     if (flg) {
10629566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads));
106385649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
106485649d77SJunchao Zhang     }
106548a46eb9SPierre Jolivet     if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads));
106685649d77SJunchao Zhang   }
106785649d77SJunchao Zhang #endif
106885649d77SJunchao Zhang 
106985649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
107085649d77SJunchao Zhang   /*
107185649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
107285649d77SJunchao Zhang 
107385649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
107485649d77SJunchao Zhang   */
10759566063dSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL));
107685649d77SJunchao Zhang #endif
107785649d77SJunchao Zhang 
107885649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10799566063dSJacob Faibussowitsch   PetscCall(PetscFPTCreate(10000));
108085649d77SJunchao Zhang #endif
108185649d77SJunchao Zhang 
108285649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
108385649d77SJunchao Zhang   {
108485649d77SJunchao Zhang     PetscViewer viewer;
10859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg));
108685649d77SJunchao Zhang     if (flg) {
10879566063dSJacob Faibussowitsch       PetscCall(PetscProcessPlacementView(viewer));
10889566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
108985649d77SJunchao Zhang     }
109085649d77SJunchao Zhang   }
109185649d77SJunchao Zhang #endif
109285649d77SJunchao Zhang 
109385649d77SJunchao Zhang   flg = PETSC_TRUE;
10949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL));
10959566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE));
109685649d77SJunchao Zhang 
109785649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10983ba16761SJacob Faibussowitsch   PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD);
10993ba16761SJacob Faibussowitsch   PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default);
11003ba16761SJacob Faibussowitsch   PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", "");
11013ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, "");
110285649d77SJunchao Zhang #endif
110385649d77SJunchao Zhang 
110485649d77SJunchao Zhang #if defined(__VALGRIND_H)
110585649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE;
110685649d77SJunchao Zhang   #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
11079566063dSJacob 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"));
110885649d77SJunchao Zhang   #endif
110985649d77SJunchao Zhang #endif
111085649d77SJunchao Zhang   /*
111185649d77SJunchao Zhang       Set flag that we are completely initialized
111285649d77SJunchao Zhang   */
111385649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
111485649d77SJunchao Zhang 
11159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg));
11169566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPythonInitialize(NULL, NULL));
1117f1f2ae84SBarry Smith 
1118f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1119f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin());
1120f1f2ae84SBarry 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");
11213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112285649d77SJunchao Zhang }
112385649d77SJunchao Zhang 
112410450e9eSJacob Faibussowitsch // "Unknown section 'Environmental Variables'"
112510450e9eSJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-section-header-unknown
1126e5c89e4eSSatish Balay /*@C
1127e5c89e4eSSatish Balay   PetscInitialize - Initializes the PETSc database and MPI.
1128811af0c4SBarry Smith   `PetscInitialize()` calls MPI_Init() if that has yet to be called,
1129e5c89e4eSSatish Balay   so this routine should always be called near the beginning of
1130e5c89e4eSSatish Balay   your program -- usually the very first line!
1131e5c89e4eSSatish Balay 
1132811af0c4SBarry Smith   Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set
1133e5c89e4eSSatish Balay 
1134e5c89e4eSSatish Balay   Input Parameters:
1135e5c89e4eSSatish Balay + argc - count of number of command line arguments
1136e5c89e4eSSatish Balay . args - the command line arguments
1137be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1138be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1139be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1140c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
11410298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message
1142e5c89e4eSSatish Balay 
1143811af0c4SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that
1144811af0c4SBarry Smith    communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a
1145811af0c4SBarry Smith    four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not,
1146811af0c4SBarry Smith    then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even
114705827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1148e5c89e4eSSatish Balay 
1149e5c89e4eSSatish Balay   Options Database Keys:
11507ca660e7SBarry Smith + -help [intro]                                       - prints help method for each option; if intro is given the program stops after printing the introductory help message
11517ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...]        - Starts program in debugger
1152e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11538a690491SBarry Smith . -on_error_emacs <machinename>                       - causes emacsclient to jump to error file
1154811af0c4SBarry Smith . -on_error_abort                                     - calls `abort()` when error detected (no traceback)
1155811af0c4SBarry Smith . -on_error_mpiabort                                  - calls `MPI_abort()` when error detected
11561af3601dSBarry Smith . -error_output_stdout                                - prints PETSc error messages to stdout instead of the default stderr
11578a690491SBarry Smith . -error_output_none                                  - does not print the error messages (but handles errors in the same way as if this was not called)
1158bf4d2887SBarry Smith . -debugger_ranks [rank1,rank2,...]                   - Indicates ranks to start in debugger
1159e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds)            - Pauses debugger
1160e5c89e4eSSatish Balay . -stop_for_debugger                                  - Print message on how to attach debugger manually to
1161e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
1162aee23540SBarry Smith . -malloc_dump                                        - prints a list of all unfreed memory at the end of the run
116392f119d6SBarry 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
1164811af0c4SBarry Smith . -malloc_view                                        - show a list of all allocated memory during `PetscFinalize()`
116592f119d6SBarry Smith . -malloc_view_threshold <t>                          - only list memory allocations of size greater than t with -malloc_view
1166608c71bfSMatthew G. Knepley . -malloc_requested_size                              - malloc logging will record the requested size rather than size after alignment
116792f119d6SBarry Smith . -fp_trap                                            - Stops on floating point exceptions
1168e5c89e4eSSatish Balay . -no_signal_handler                                  - Indicates not to trap error signals
1169e5c89e4eSSatish Balay . -shared_tmp                                         - indicates /tmp directory is shared by all processors
1170e5c89e4eSSatish Balay . -not_shared_tmp                                     - each processor has own /tmp
1171e5c89e4eSSatish Balay . -tmp                                                - alternative name of /tmp directory
1172e5c89e4eSSatish Balay . -get_total_flops                                    - returns total flops done by all processors
11730841954dSBarry Smith - -memory_view                                        - Print memory usage at end of run
1174e5c89e4eSSatish Balay 
1175c5b5d8d5SVaclav Hapla   Options Database Keys for Option Database:
1176c5b5d8d5SVaclav Hapla + -skip_petscrc           - skip the default option files ~/.petscrc, .petscrc, petscrc
1177c5b5d8d5SVaclav Hapla . -options_monitor        - monitor all set options to standard output for the whole program run
1178811af0c4SBarry Smith - -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()`
1179c5b5d8d5SVaclav Hapla 
1180c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1181c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1182c5b5d8d5SVaclav Hapla    They can be used also in option files.
1183c5b5d8d5SVaclav Hapla 
1184811af0c4SBarry Smith    See `PetscOptionsMonitorSet()` to do monitoring programmatically.
1185c5b5d8d5SVaclav Hapla 
1186e5c89e4eSSatish Balay   Options Database Keys for Profiling:
1187a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1188811af0c4SBarry Smith + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`.
1189217044c2SLisandro Dalcin . -log_sync                                            - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1190217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1191495fc317SBarry Smith . -log_trace [filename]                                - Print traces of all PETSc calls to the screen (useful to determine where a program
1192811af0c4SBarry Smith         hangs without running in the debugger).  See `PetscLogTraceBegin()`.
1193*ad2e3d55SToby Isaac . -log_view [:filename:format][,[:filename:format]...] - Prints summary of flop and timing information to screen or file, see `PetscLogView()` (up to 4 viewers)
1194811af0c4SBarry Smith . -log_view_memory                                     - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`.
1195811af0c4SBarry Smith . -log_view_gpu_time                                   - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView().
1196f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes>                  - excludes subset of object classes from logging
1197b665b14eSToby Isaac . -log [filename]                                      - Logs profiling information in a dump file, see `PetscLogDump()`.
1198b665b14eSToby Isaac . -log_all [filename]                                  - Same as `-log`.
1199c2f74817SBarry Smith . -log_mpe [filename]                                  - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
12002611ad71SToby Isaac . -log_perfstubs                                       - Starts a log handler with the perfstubs interface (which is used by TAU)
1201811af0c4SBarry Smith . -viewfromoptions on,off                              - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off
1202c2f74817SBarry Smith - -check_pointer_intensity 0,1,2                       - if pointers are checked for validity (debug version only), using 0 will result in faster code
1203495fc317SBarry Smith 
1204ffbd1cfbSBarry Smith   Options Database Keys for SAWs:
1205ffbd1cfbSBarry Smith + -saws_port <portnumber>        - port number to publish SAWs data, default is 8080
1206ffbd1cfbSBarry 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
1207ffbd1cfbSBarry Smith                                    this is useful when you are running many jobs that utilize SAWs at the same time
1208ffbd1cfbSBarry Smith . -saws_log <filename>           - save a log of all SAWs communication
1209ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1210ffbd1cfbSBarry Smith - -saws_root <directory>         - allow SAWs to have access to the given directory to search for requested resources and files
1211ffbd1cfbSBarry Smith 
1212e5c89e4eSSatish Balay   Environmental Variables:
1213811af0c4SBarry Smith +   `PETSC_TMP` - alternative tmp directory
1214811af0c4SBarry Smith .   `PETSC_SHARED_TMP` - tmp is shared by all processes
1215811af0c4SBarry Smith .   `PETSC_NOT_SHARED_TMP` - each process has its own private tmp
1216811af0c4SBarry Smith .   `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs
1217811af0c4SBarry Smith .   `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1218811af0c4SBarry Smith .   `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer
1219811af0c4SBarry Smith -   `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to
1220e5c89e4eSSatish Balay 
1221e5c89e4eSSatish Balay   Level: beginner
1222e5c89e4eSSatish Balay 
1223811af0c4SBarry Smith   Note:
1224811af0c4SBarry Smith   If for some reason you must call `MPI_Init()` separately, call
1225811af0c4SBarry Smith   it before `PetscInitialize()`.
1226e5c89e4eSSatish Balay 
1227811af0c4SBarry Smith   Fortran Notes:
1228811af0c4SBarry Smith   In Fortran this routine can be called with
1229811af0c4SBarry Smith .vb
1230811af0c4SBarry Smith        call PetscInitialize(ierr)
1231811af0c4SBarry Smith        call PetscInitialize(file,ierr) or
1232811af0c4SBarry Smith        call PetscInitialize(file,help,ierr)
1233811af0c4SBarry Smith .ve
1234e5c89e4eSSatish Balay 
1235811af0c4SBarry Smith   If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after
1236811af0c4SBarry Smith   calling `PetscInitialize()`.
1237e5c89e4eSSatish Balay 
12385dedd997SBarry Smith   Options Database Key for Developers:
12395dedd997SBarry Smith . -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form:
12405dedd997SBarry Smith     "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py)
12415dedd997SBarry Smith 
1242db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()`
1243e5c89e4eSSatish Balay @*/
1244d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[])
1245d71ae5a4SJacob Faibussowitsch {
124685649d77SJunchao Zhang   PetscMPIInt flag;
124721ef0414SBarry Smith   const char *prog = "Unknown Name", *mpienv;
1248e5c89e4eSSatish Balay 
124927104ee2SJacob Faibussowitsch   PetscFunctionBegin;
12503ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
12519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Initialized(&flag));
1252e5c89e4eSSatish Balay   if (!flag) {
125339a651e2SJacob 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");
12549566063dSJacob Faibussowitsch     PetscCall(PetscPreMPIInit_Private());
12555e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12565e765c61SJed Brown     {
125739a651e2SJacob Faibussowitsch       PetscMPIInt PETSC_UNUSED provided;
12589566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED, &provided));
12595e765c61SJed Brown     }
12605e765c61SJed Brown #else
12619566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Init(argc, args));
12625e765c61SJed Brown #endif
126321ef0414SBarry Smith     if (PetscDefined(HAVE_MPIUNI)) {
126421ef0414SBarry Smith       mpienv = getenv("PMI_SIZE");
126521ef0414SBarry Smith       if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE");
126621ef0414SBarry Smith       if (mpienv) {
126721ef0414SBarry Smith         PetscInt isize;
126821ef0414SBarry Smith         PetscCall(PetscOptionsStringToInt(mpienv, &isize));
126921ef0414SBarry Smith         if (isize != 1) printf("You are using an MPI-uni (sequential) install of PETSc but trying to launch parallel jobs; you need full MPI version of PETSc\n");
127021ef0414SBarry Smith         PetscCheck(isize == 1, MPI_COMM_SELF, PETSC_ERR_MPI, "You are using an MPI-uni (sequential) install of PETSc but trying to launch parallel jobs; you need full MPI version of PETSc");
127121ef0414SBarry Smith       }
127221ef0414SBarry Smith     }
1273e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1274e5c89e4eSSatish Balay   }
12754dfee713SSatish Balay 
127685649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1277e5c89e4eSSatish Balay   if (argc && args) {
1278e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1279e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1280e5c89e4eSSatish Balay   }
1281811af0c4SBarry Smith   PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0));
12823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1283e5c89e4eSSatish Balay }
1284e5c89e4eSSatish Balay 
128595c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1286ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsCounts;
1287ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsMaxCounts;
128895c0884eSLisandro Dalcin PETSC_INTERN PetscBool    PetscObjectsLog;
1289e5c89e4eSSatish Balay 
1290008a6e76SBarry Smith /*
1291008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1292008a6e76SBarry Smith */
1293d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void)
1294d71ae5a4SJacob Faibussowitsch {
1295008a6e76SBarry Smith   PetscFunctionBegin;
1296613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
12979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128));
129874c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
12999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128));
1300008a6e76SBarry Smith   #endif
130174c01ffaSSatish Balay #endif
13029e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
13039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FP16));
1304008a6e76SBarry Smith #endif
1305008a6e76SBarry Smith 
1306de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
13079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_SUM));
1308613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MAX));
1309613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MIN));
13109e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
13119e517322SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128));
1312008a6e76SBarry Smith #endif
1313008a6e76SBarry Smith 
13149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR));
13159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT));
13169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT));
131740df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
13189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2INT));
1319008a6e76SBarry Smith #endif
13209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPI_4INT));
13219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_4INT));
13229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
132362e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp));
13243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1325008a6e76SBarry Smith }
1326008a6e76SBarry Smith 
1327a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1328a4af0ceeSJacob Faibussowitsch 
1329e5c89e4eSSatish Balay /*@C
1330e5c89e4eSSatish Balay   PetscFinalize - Checks for options to be called at the conclusion
1331811af0c4SBarry Smith   of the program. `MPI_Finalize()` is called only if the user had not
1332811af0c4SBarry Smith   called `MPI_Init()` before calling `PetscInitialize()`.
1333e5c89e4eSSatish Balay 
1334811af0c4SBarry Smith   Collective on `PETSC_COMM_WORLD`
1335e5c89e4eSSatish Balay 
1336e5c89e4eSSatish Balay   Options Database Keys:
1337811af0c4SBarry Smith + -options_view                    - Calls `PetscOptionsView()`
1338e5c89e4eSSatish Balay . -options_left                    - Prints unused options that remain in the database
13397eb1d149SBarry 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
1340e5c89e4eSSatish Balay . -mpidump                         - Calls PetscMPIDump()
1341811af0c4SBarry Smith . -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed
13424bd3d7f8SBarry Smith . -memory_view                     - Prints total memory usage
13434bd3d7f8SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and in what functions
1344e5c89e4eSSatish Balay 
1345e5c89e4eSSatish Balay   Level: beginner
1346e5c89e4eSSatish Balay 
1347e5c89e4eSSatish Balay   Note:
1348811af0c4SBarry Smith   See `PetscInitialize()` for other runtime options.
1349e5c89e4eSSatish Balay 
1350db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()`
1351e5c89e4eSSatish Balay @*/
1352d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void)
1353d71ae5a4SJacob Faibussowitsch {
13544bb5149bSJed Brown   PetscMPIInt rank;
1355a8d2bbe5SBarry Smith   PetscInt    nopt;
13562bf49c77SBarry Smith   PetscBool   flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE;
1357dff31646SBarry Smith   PetscBool   flg;
135810463e74SBarry Smith   char        mname[PETSC_MAX_PATH_LEN];
1359e5c89e4eSSatish Balay 
13603cbbc5ffSBarry Smith   PetscFunctionBegin;
136139a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()");
13629566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PetscFinalize() called\n"));
1363b022a5c1SBarry Smith 
1364f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1365f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd());
1366f1f2ae84SBarry Smith 
136762e5d2d2SJDBetteridge   /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */
136862e5d2d2SJDBetteridge   {
136962e5d2d2SJDBetteridge     union
137062e5d2d2SJDBetteridge     {
137162e5d2d2SJDBetteridge       MPI_Comm comm;
137262e5d2d2SJDBetteridge       void    *ptr;
137362e5d2d2SJDBetteridge     } ucomm;
137462e5d2d2SJDBetteridge     PetscMPIInt flg;
137562e5d2d2SJDBetteridge     void       *tmp;
137662e5d2d2SJDBetteridge 
137762e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
137862e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
137962e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF));
138062e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
138162e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
138262e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD));
138362e5d2d2SJDBetteridge   }
138462e5d2d2SJDBetteridge 
13859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1386a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13873ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE);
13883ba16761SJacob Faibussowitsch   PetscCallExternal(adios_finalize, rank);
1389a56f64adSBarry Smith #endif
13909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg));
1391dff31646SBarry Smith   if (flg) {
13921f817a21SBarry Smith     char *cits, filename[PETSC_MAX_PATH_LEN];
13931f817a21SBarry Smith     FILE *fd = PETSC_STDOUT;
13941f817a21SBarry Smith 
13959566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL));
139648a46eb9SPierre Jolivet     if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd));
13979566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits));
1398dff31646SBarry Smith     cits[0] = 0;
13999566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits));
14009566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n"));
14019566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14029566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits));
14039566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14049566063dSJacob Faibussowitsch     PetscCall(PetscFClose(PETSC_COMM_WORLD, fd));
14059566063dSJacob Faibussowitsch     PetscCall(PetscFree(cits));
1406dff31646SBarry Smith   }
14079566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&PetscCitationsList));
1408dff31646SBarry Smith 
14092d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14109566063dSJacob Faibussowitsch   PetscCall(PetscFPTDestroy());
14112d53ad75SBarry Smith #endif
14122d53ad75SBarry Smith 
1413e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1414dff31646SBarry Smith   flg = PETSC_FALSE;
14159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL));
14161baa6e33SBarry Smith   if (flg) PetscCall(PetscOptionsSAWsDestroy());
1417d5649816SBarry Smith #endif
1418d5649816SBarry Smith 
1419681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1420681455b2SBarry Smith   flg1 = PETSC_FALSE;
14219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL));
1422681455b2SBarry Smith   if (flg1) {
1423681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
14249566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL));
1425681455b2SBarry Smith   }
1426681455b2SBarry Smith #endif
1427681455b2SBarry Smith 
142867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL));
143048a46eb9SPierre Jolivet   if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n"));
143167584ceeSBarry Smith #endif
1432e5c89e4eSSatish Balay 
14332611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
143490d69ab7SBarry Smith     flg1 = PETSC_FALSE;
14359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL));
1436e5c89e4eSSatish Balay     if (flg1) {
1437e5c89e4eSSatish Balay       PetscLogDouble flops = 0;
14389566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD));
14399566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops));
1440e5c89e4eSSatish Balay     }
14412611ad71SToby Isaac   }
1442e5c89e4eSSatish Balay 
14432611ad71SToby Isaac   if (PetscDefined(USE_LOG) && PetscDefined(HAVE_MPE)) {
1444e5c89e4eSSatish Balay     mname[0] = 0;
14459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1));
14462611ad71SToby Isaac     if (flg1) PetscCall(PetscLogMPEDump(mname[0] ? mname : NULL));
1447e5c89e4eSSatish Balay   }
1448a297a907SKarl Rupp 
1449dd710f27SBarry Smith   /*
1450dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1451dd710f27SBarry Smith   */
14529566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1453dd710f27SBarry Smith 
14542611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
14559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE));
14569566063dSJacob Faibussowitsch     PetscCall(PetscLogViewFromOptions());
14579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsPopGetViewerOff());
145887c3beb6SLisandro Dalcin 
1459dd710f27SBarry Smith     /*
1460dd710f27SBarry Smith        Free any objects created by the last block of code.
1461dd710f27SBarry Smith        */
14629566063dSJacob Faibussowitsch     PetscCall(PetscObjectRegisterDestroyAll());
1463dd710f27SBarry Smith 
1464dd710f27SBarry Smith     mname[0] = 0;
14659566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1));
14669566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2));
14679566063dSJacob Faibussowitsch     if (flg1 || flg2) PetscCall(PetscLogDump(mname));
14682611ad71SToby Isaac   }
146910463e74SBarry Smith 
147090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL));
14729566063dSJacob Faibussowitsch   if (!flg1) PetscCall(PetscPopSignalHandler());
147390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL));
14751baa6e33SBarry Smith   if (flg1) PetscCall(PetscMPIDump(stdout));
147690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
147790d69ab7SBarry Smith   flg2 = PETSC_FALSE;
14788bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
14799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1));
14809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
14819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL));
1482e4c476e2SSatish Balay 
1483e5c89e4eSSatish Balay   if (flg2) {
1484be56827dSJed Brown     PetscViewer viewer;
14859566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
14869566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
14879566063dSJacob Faibussowitsch     PetscCall(PetscOptionsView(NULL, viewer));
14889566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1489e5c89e4eSSatish Balay   }
1490e5c89e4eSSatish Balay 
1491e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
14929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1));
14939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1));
1494e5c89e4eSSatish Balay 
149533fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
14969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1));
149759f199a7SJed Brown   if (!flg1) flg3 = PETSC_TRUE;
1498e5c89e4eSSatish Balay   if (flg3) {
14993de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1500be56827dSJed Brown       PetscViewer viewer;
15019566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15029566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15039566063dSJacob Faibussowitsch       PetscCall(PetscOptionsView(NULL, viewer));
15049566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
1505e5c89e4eSSatish Balay     }
15069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsAllUsed(NULL, &nopt));
15073de2bfdfSBarry Smith     if (nopt) {
15089566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n"));
15099566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n"));
15103de2bfdfSBarry Smith       if (nopt == 1) {
15119566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n"));
1512e5c89e4eSSatish Balay       } else {
15139566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt));
1514e5c89e4eSSatish Balay       }
15153de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15169566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n"));
1517df12ba86SBarry Smith     }
15189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsLeft(NULL));
1519e5c89e4eSSatish Balay   }
1520e5c89e4eSSatish Balay 
1521e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1522d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15239566063dSJacob Faibussowitsch     PetscCall(PetscStackSAWsViewOff());
1524792fecdfSBarry Smith     PetscCallSAWs(SAWs_Finalize, ());
1525d45a07a7SBarry Smith   }
1526ec957eceSBarry Smith #endif
1527ec957eceSBarry Smith 
152810463e74SBarry Smith   /*
1529dbc8283eSBarry Smith        List all objects the user may have forgot to free
15302eff7a51SBarry Smith   */
15312611ad71SToby Isaac   if (PetscDefined(USE_LOG) && PetscObjectsLog) {
15329566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
1533a64a8e02SBarry Smith     if (flg1) {
1534a64a8e02SBarry Smith       MPI_Comm local_comm;
15357eb1d149SBarry Smith       char     string[64];
1536a64a8e02SBarry Smith 
15379566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL));
1538afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
15399566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
15409566063dSJacob Faibussowitsch       PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
15419566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
15429566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
15430a1571b3SBarry Smith     }
154405df10baSBarry Smith   }
15454097062eSBarry Smith 
1546dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1547dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
15489566063dSJacob Faibussowitsch   PetscCall(PetscFree(PetscObjects));
15492eff7a51SBarry Smith 
155033f85c2fSBarry Smith   /*
155133f85c2fSBarry Smith      Destroy any packages that registered a finalize
155233f85c2fSBarry Smith   */
15539566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalizeAll());
155433f85c2fSBarry Smith 
15559566063dSJacob Faibussowitsch   PetscCall(PetscLogFinalize());
1556101409b8SToby Isaac 
155733f85c2fSBarry Smith   /*
155848dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
155948dd1dffSBarry Smith   */
15602e956fe4SStefano Zampini   if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll());
156137e93019SBarry Smith 
15624028d114SSatish Balay   if (petsc_history) {
15639566063dSJacob Faibussowitsch     PetscCall(PetscCloseHistoryFile(&petsc_history));
156402c9f0b5SLisandro Dalcin     petsc_history = NULL;
1565e5c89e4eSSatish Balay   }
15669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
15679566063dSJacob Faibussowitsch   PetscCall(PetscInfoDestroy());
1568e5c89e4eSSatish Balay 
156967584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
157092f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1571e5c89e4eSSatish Balay     char  fname[PETSC_MAX_PATH_LEN];
157292f119d6SBarry Smith     char  sname[PETSC_MAX_PATH_LEN];
1573e5c89e4eSSatish Balay     FILE *fd;
1574ed9cf6e9SBarry Smith     int   err;
1575e5c89e4eSSatish Balay 
1576dc92acbaSJed Brown     flg2 = PETSC_FALSE;
157792f119d6SBarry Smith     flg3 = PETSC_FALSE;
15789566063dSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL));
15799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL));
158092f119d6SBarry Smith     fname[0] = 0;
15819566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1));
1582e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
15833ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
15849371c9d4SSatish Balay       fd = fopen(sname, "w");
15859371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
15869566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(fd));
1587ed9cf6e9SBarry Smith       err = fclose(fd);
158828b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
158992f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1590e5c89e4eSSatish Balay       MPI_Comm local_comm;
1591e5c89e4eSSatish Balay 
1592afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
15939566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
15949566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(stdout));
15959566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
15969566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1597e5c89e4eSSatish Balay     }
1598e5c89e4eSSatish Balay     fname[0] = 0;
15999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1));
1600e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
16013ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
16029371c9d4SSatish Balay       fd = fopen(sname, "w");
16039371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16049566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(fd));
1605ed9cf6e9SBarry Smith       err = fclose(fd);
160628b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
160792f119d6SBarry Smith     } else if (flg1) {
160892f119d6SBarry Smith       MPI_Comm local_comm;
160992f119d6SBarry Smith 
1610afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16119566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16129566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(stdout));
16139566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16149566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1615e5c89e4eSSatish Balay     }
1616e5c89e4eSSatish Balay   }
161767584ceeSBarry Smith #endif
161820e2c332SMatthew G. Knepley 
16195486ca60SMatthew G. Knepley   /*
16205486ca60SMatthew G. Knepley      Close any open dynamic libraries
16215486ca60SMatthew G. Knepley   */
16229566063dSJacob Faibussowitsch   PetscCall(PetscFinalize_DynamicLibraries());
16235486ca60SMatthew G. Knepley 
1624e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroyDefault());
1626e5c89e4eSSatish Balay 
1627e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
162802c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1629e5c89e4eSSatish Balay 
1630c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1631c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
16329566063dSJacob Faibussowitsch     PetscCall(PetscKokkosFinalize_Private());
1633c2b86a48SJunchao Zhang     PetscBeganKokkos       = PETSC_FALSE;
163445639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1635c2b86a48SJunchao Zhang   }
1636c2b86a48SJunchao Zhang #endif
1637c2b86a48SJunchao Zhang 
163871438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
163971438e86SJunchao Zhang   if (PetscBeganNvshmem) {
16409566063dSJacob Faibussowitsch     PetscCall(PetscNvshmemFinalize());
164171438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
164271438e86SJunchao Zhang   }
164371438e86SJunchao Zhang #endif
1644a0e72f99SJunchao Zhang 
16459566063dSJacob Faibussowitsch   PetscCall(PetscFreeMPIResources());
1646e5c89e4eSSatish Balay 
1647dbc8283eSBarry Smith   /*
1648efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1649efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1650efb80d3cSBarry Smith 
1651efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1652efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1653dbc8283eSBarry Smith  */
1654b770b1f6SSatish Balay   {
1655dbc8283eSBarry Smith     PetscCommCounter *counter;
1656dbc8283eSBarry Smith     PetscMPIInt       flg;
1657dbc8283eSBarry Smith     MPI_Comm          icomm;
16589371c9d4SSatish Balay     union
16599371c9d4SSatish Balay     {
16609371c9d4SSatish Balay       MPI_Comm comm;
16619371c9d4SSatish Balay       void    *ptr;
16629371c9d4SSatish Balay     } ucomm;
16639566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
1664dbc8283eSBarry Smith     if (flg) {
1665265f3f35SJed Brown       icomm = ucomm.comm;
16669566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
166728b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1668dbc8283eSBarry Smith 
16699566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval));
16709566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
16719566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1672dbc8283eSBarry Smith     }
16739566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
1674dbc8283eSBarry Smith     if (flg) {
1675265f3f35SJed Brown       icomm = ucomm.comm;
16769566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
167728b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1678dbc8283eSBarry Smith 
16799566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval));
16809566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
16819566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1682dbc8283eSBarry Smith     }
1683b770b1f6SSatish Balay   }
1684dbc8283eSBarry Smith 
16859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
16869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
16879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
16889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
168962e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval));
169062e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval));
1691480cf27aSJed Brown 
16925ea2b939SDuncan Campbell   // Free keyvals which may be silently created by some routines
16935ea2b939SDuncan Campbell   if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval));
16945ea2b939SDuncan Campbell   if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval));
16955ea2b939SDuncan Campbell 
16969566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
16979566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
16989566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
16999566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock));
1700ef19f930SBarry Smith 
1701e5c89e4eSSatish Balay   if (PetscBeganMPI) {
170299b1327fSBarry Smith     PetscMPIInt flag;
17039566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalized(&flag));
170439a651e2SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
170539a651e2SJacob Faibussowitsch     /* wait until the very last moment to disable error handling */
170639a651e2SJacob Faibussowitsch     PetscErrorHandlingInitialized = PETSC_FALSE;
17079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalize());
170839a651e2SJacob Faibussowitsch   } else PetscErrorHandlingInitialized = PETSC_FALSE;
170939a651e2SJacob Faibussowitsch 
1710e5c89e4eSSatish Balay   /*
1711e5c89e4eSSatish Balay 
1712e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1713e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1714e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1715e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1716e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17170e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1718e5c89e4eSSatish Balay    memory was not freed.
1719e5c89e4eSSatish Balay 
1720e5c89e4eSSatish Balay */
17219566063dSJacob Faibussowitsch   PetscCall(PetscMallocClear());
17229566063dSJacob Faibussowitsch   PetscCall(PetscStackReset());
1723a297a907SKarl Rupp 
1724e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1725e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
17267ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE)
172756883f60SBarry Smith   /*
172856883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
172956883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
173056883f60SBarry Smith    */
173156883f60SBarry Smith   __gcov_flush();
173256883f60SBarry Smith #endif
17331724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
17341724198aSStefano Zampini   PetscStackClearTop;
17353ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1736e5c89e4eSSatish Balay }
1737e5c89e4eSSatish Balay 
173843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
1739d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b)
1740d71ae5a4SJacob Faibussowitsch {
174143db4dbbSBarry Smith   if (*a == *b) return 1;
174243db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
174343db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
174443db4dbbSBarry Smith   return 0;
174543db4dbbSBarry Smith }
1746a70650f6SBarry Smith #endif
174743db4dbbSBarry Smith 
174843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
1749d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b)
1750d71ae5a4SJacob Faibussowitsch {
175143db4dbbSBarry Smith   if (*a == *b) return 1;
175243db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
175343db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
175443db4dbbSBarry Smith   return 0;
175543db4dbbSBarry Smith }
175643db4dbbSBarry Smith #endif
1757