xref: /petsc/src/sys/objects/pinit.c (revision 835d5ba202af4ef59a256bb814952841bbb5ae5d)
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*/
6473903fcSJunchao Zhang #include <petsc/private/logimpl.h>
7665c2dedSJed Brown #include <petscviewer.h>
862e5d2d2SJDBetteridge #include <petsc/private/garbagecollector.h>
9a0e72f99SJunchao Zhang 
106edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
116edef35eSSatish Balay   #include <petsc/private/valgrind/valgrind.h>
126edef35eSSatish Balay #endif
136edef35eSSatish Balay 
14fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS)
1585649d77SJunchao Zhang   #include <petsc/private/fortranimpl.h>
1685649d77SJunchao Zhang #endif
1785649d77SJunchao Zhang 
187ce81a4bSJacob Faibussowitsch #if PetscDefined(USE_COVERAGE)
1956883f60SBarry Smith EXTERN_C_BEGIN
20aaf3846cSSatish Balay   #if defined(PETSC_HAVE___GCOV_DUMP)
21aaf3846cSSatish Balay     #define __gcov_flush(x) __gcov_dump(x)
22aaf3846cSSatish Balay   #endif
2356883f60SBarry Smith void __gcov_flush(void);
2456883f60SBarry Smith EXTERN_C_END
2556883f60SBarry Smith #endif
268101f56cSMatthew Knepley 
272d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
2895c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
292d53ad75SBarry Smith PetscFPT              PetscFPTData = 0;
302d53ad75SBarry Smith #endif
312d53ad75SBarry Smith 
3227104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
3316ad0300SBarry Smith   #include <petscviewersaws.h>
34a6790183SBarry Smith #endif
35eb02082bSJunchao Zhang 
3695c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
37e5c89e4eSSatish Balay 
3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm, int);
4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm, int);
4295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE **);
430069ddf5SShri Abhyankar 
446de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
45e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
47*835d5ba2SJunchao Zhang PetscMPIInt PETSC_MPI_THREAD_REQUIRED = PETSC_DECIDE;
486de5d289SStefano Zampini #else
49*835d5ba2SJunchao Zhang PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_SINGLE;
506de5d289SStefano Zampini #endif
51e5c89e4eSSatish Balay 
52480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval      = MPI_KEYVAL_INVALID;
53480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval    = MPI_KEYVAL_INVALID;
54480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval    = MPI_KEYVAL_INVALID;
555f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval      = MPI_KEYVAL_INVALID;
5662e5d2d2SJDBetteridge PetscMPIInt Petsc_CreationIdx_keyval  = MPI_KEYVAL_INVALID;
5762e5d2d2SJDBetteridge PetscMPIInt Petsc_Garbage_HMap_keyval = MPI_KEYVAL_INVALID;
58480cf27aSJed Brown 
595ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedWD_keyval  = MPI_KEYVAL_INVALID;
605ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedTmp_keyval = MPI_KEYVAL_INVALID;
615ea2b939SDuncan Campbell 
62e5c89e4eSSatish Balay /*
63e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
64e5c89e4eSSatish Balay */
6502c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL};
6602c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL};
67e5c89e4eSSatish Balay 
68ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
69ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
700f8e0872SSatish Balay 
71a2f94806SJed Brown PetscInt PetscHotRegionDepth;
72a2f94806SJed Brown 
736edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
746edef35eSSatish Balay 
75b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
76b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
77b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
78b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
79b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
80b22622e2STadeu Manoel #endif
81b22622e2STadeu Manoel 
82aec76313SJacob Faibussowitsch /*@C
83945d1669SBarry Smith   PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8472a42c3cSBarry Smith 
8572a42c3cSBarry Smith   Collective
8672a42c3cSBarry Smith 
8710450e9eSJacob Faibussowitsch   Input Parameters:
8810450e9eSJacob Faibussowitsch + argc     - number of args
8910450e9eSJacob Faibussowitsch . args     - array of command line arguments
9010450e9eSJacob Faibussowitsch . filename - optional name of the program file, pass `NULL` to ignore
9110450e9eSJacob Faibussowitsch - help     - optional help, pass `NULL` to ignore
9210450e9eSJacob Faibussowitsch 
9372a42c3cSBarry Smith   Level: advanced
9472a42c3cSBarry Smith 
9595452b02SPatrick Sanan   Notes:
96a56f64adSBarry Smith   this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
970f11a792SBarry Smith   indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
98a56f64adSBarry Smith   be called multiple times from Julia without the problem of trying to initialize MPI more than once.
990f11a792SBarry Smith 
10010450e9eSJacob Faibussowitsch   Developer Notes:
10110450e9eSJacob Faibussowitsch   Turns off PETSc signal handling to allow Julia to manage signals
1021ea5a559SBarry Smith 
103db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()`
1040f11a792SBarry Smith */
105d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoPointers(int argc, char **args, const char *filename, const char *help)
106d71ae5a4SJacob Faibussowitsch {
10772a42c3cSBarry Smith   int    myargc = argc;
10872a42c3cSBarry Smith   char **myargs = args;
10972a42c3cSBarry Smith 
11072a42c3cSBarry Smith   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&myargc, &myargs, filename, help));
1129566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
113df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11572a42c3cSBarry Smith }
11672a42c3cSBarry Smith 
117e5c89e4eSSatish Balay /*@C
118811af0c4SBarry Smith   PetscInitializeNoArguments - Calls `PetscInitialize()` from C/C++ without
119e5c89e4eSSatish Balay   the command line arguments.
120e5c89e4eSSatish Balay 
121e5c89e4eSSatish Balay   Collective
122e5c89e4eSSatish Balay 
123e5c89e4eSSatish Balay   Level: advanced
124e5c89e4eSSatish Balay 
125db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`
126e5c89e4eSSatish Balay @*/
127d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoArguments(void)
128d71ae5a4SJacob Faibussowitsch {
129e5c89e4eSSatish Balay   int    argc = 0;
13002c9f0b5SLisandro Dalcin   char **args = NULL;
131e5c89e4eSSatish Balay 
132e5c89e4eSSatish Balay   PetscFunctionBegin;
1339566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, NULL, NULL));
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135e5c89e4eSSatish Balay }
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay /*@
138e5c89e4eSSatish Balay   PetscInitialized - Determine whether PETSc is initialized.
139e5c89e4eSSatish Balay 
14010450e9eSJacob Faibussowitsch   Output Parameter:
14110450e9eSJacob Faibussowitsch . isInitialized - `PETSC_TRUE` if PETSc is initialized, `PETSC_FALSE` otherwise
14210450e9eSJacob Faibussowitsch 
14393b6d2d1SJed Brown   Level: beginner
144e5c89e4eSSatish Balay 
145db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
146e5c89e4eSSatish Balay @*/
147d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialized(PetscBool *isInitialized)
148d71ae5a4SJacob Faibussowitsch {
14927104ee2SJacob Faibussowitsch   PetscFunctionBegin;
1504f572ea9SToby Isaac   if (PetscInitializeCalled) PetscAssertPointer(isInitialized, 1);
151e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153e5c89e4eSSatish Balay }
154e5c89e4eSSatish Balay 
155e5c89e4eSSatish Balay /*@
156811af0c4SBarry Smith   PetscFinalized - Determine whether `PetscFinalize()` has been called yet
157e5c89e4eSSatish Balay 
15810450e9eSJacob Faibussowitsch   Output Parameter:
15910450e9eSJacob Faibussowitsch . isFinalized - `PETSC_TRUE` if PETSc is finalized, `PETSC_FALSE` otherwise
16010450e9eSJacob Faibussowitsch 
161e5c89e4eSSatish Balay   Level: developer
162e5c89e4eSSatish Balay 
163db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
164e5c89e4eSSatish Balay @*/
165d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalized(PetscBool *isFinalized)
166d71ae5a4SJacob Faibussowitsch {
16727104ee2SJacob Faibussowitsch   PetscFunctionBegin;
1684f572ea9SToby Isaac   if (!PetscFinalizeCalled) PetscAssertPointer(isFinalized, 1);
169e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171e5c89e4eSSatish Balay }
172e5c89e4eSSatish Balay 
17357171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]);
174e5c89e4eSSatish Balay 
175e5c89e4eSSatish Balay /*
176e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
177e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
178e5c89e4eSSatish Balay */
179367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP               = 0;
18062e5d2d2SJDBetteridge MPI_Op Petsc_Garbage_SetIntersectOp = 0;
181e5c89e4eSSatish Balay 
182d71ae5a4SJacob Faibussowitsch PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, int *cnt, MPI_Datatype *datatype)
183d71ae5a4SJacob Faibussowitsch {
184e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out, i, count = *cnt;
185e5c89e4eSSatish Balay 
186e5c89e4eSSatish Balay   PetscFunctionBegin;
187e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
1883ba16761SJacob Faibussowitsch     PetscErrorCode ierr = (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
1893ba16761SJacob Faibussowitsch     (void)ierr;
19041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
191e5c89e4eSSatish Balay   }
192e5c89e4eSSatish Balay 
193e5c89e4eSSatish Balay   for (i = 0; i < count; i++) {
194e5c89e4eSSatish Balay     xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]);
195e5c89e4eSSatish Balay     xout[2 * i + 1] += xin[2 * i + 1];
196e5c89e4eSSatish Balay   }
197812af9f3SBarry Smith   PetscFunctionReturnVoid();
198e5c89e4eSSatish Balay }
199e5c89e4eSSatish Balay 
200e5c89e4eSSatish Balay /*
201e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
202e5c89e4eSSatish Balay sum of the second entry.
203b693b147SBarry Smith 
20476f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
205367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
206b693b147SBarry Smith there would be no place to store the both needed results.
207e5c89e4eSSatish Balay */
208d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt sizes[], PetscInt *max, PetscInt *sum)
209d71ae5a4SJacob Faibussowitsch {
210e5c89e4eSSatish Balay   PetscFunctionBegin;
211d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
212d6e4c47cSJed Brown   {
2139371c9d4SSatish Balay     struct {
2149371c9d4SSatish Balay       PetscInt max, sum;
2159371c9d4SSatish Balay     } work;
2169566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce_scatter_block((void *)sizes, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm));
217d6e4c47cSJed Brown     *max = work.max;
218d6e4c47cSJed Brown     *sum = work.sum;
219d6e4c47cSJed Brown   }
220d6e4c47cSJed Brown #else
221d6e4c47cSJed Brown   {
222d6e4c47cSJed Brown     PetscMPIInt size, rank;
2239371c9d4SSatish Balay     struct {
2249371c9d4SSatish Balay       PetscInt max, sum;
2259371c9d4SSatish Balay     } *work;
2269566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
2279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
2289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &work));
2291c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce((void *)sizes, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm));
2306ac3741eSJed Brown     *max = work[rank].max;
2316ac3741eSJed Brown     *sum = work[rank].sum;
2329566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
233d6e4c47cSJed Brown   }
234d6e4c47cSJed Brown #endif
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
236e5c89e4eSSatish Balay }
237e5c89e4eSSatish Balay 
2389e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
239613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
240613bf2b2SPierre Jolivet     #include <quadmath.h>
241613bf2b2SPierre Jolivet   #endif
2429e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0;
243de272c7aSSatish Balay   #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
24406a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
245613bf2b2SPierre Jolivet   #endif
246e5c89e4eSSatish Balay 
247d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
248d71ae5a4SJacob Faibussowitsch {
249e5c89e4eSSatish Balay   PetscInt i, count = *cnt;
250e5c89e4eSSatish Balay 
251e5c89e4eSSatish Balay   PetscFunctionBegin;
2527c2de775SJed Brown   if (*datatype == MPIU_REAL) {
253e2e03761SBarry Smith     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
254a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2557c2de775SJed Brown   }
2567c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
257c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2587c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
259a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2607c2de775SJed Brown   }
2617c2de775SJed Brown   #endif
262613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
263613bf2b2SPierre Jolivet   else if (*datatype == MPIU___FLOAT128) {
264613bf2b2SPierre Jolivet     __float128 *xin = (__float128 *)in, *xout = (__float128 *)out;
265613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
26674c01ffaSSatish Balay     #if defined(PETSC_HAVE_COMPLEX)
2679371c9d4SSatish Balay   } else if (*datatype == MPIU___COMPLEX128) {
268613bf2b2SPierre Jolivet     __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out;
269613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
27074c01ffaSSatish Balay     #endif
271613bf2b2SPierre Jolivet   }
272613bf2b2SPierre Jolivet   #endif
2739e517322SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FP16)
2749e517322SPierre Jolivet   else if (*datatype == MPIU___FP16) {
2759e517322SPierre Jolivet     __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out;
2769e517322SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2779e517322SPierre Jolivet   }
2789e517322SPierre Jolivet   #endif
2797c2de775SJed Brown   else {
2809e517322SPierre Jolivet   #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16)
2813ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SElF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
2829e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FP16)
2833ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types"));
2849e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FLOAT128)
2853ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types"));
2869e517322SPierre Jolivet   #else
2873ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types"));
288613bf2b2SPierre Jolivet   #endif
28941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
290e2e03761SBarry Smith   }
291812af9f3SBarry Smith   PetscFunctionReturnVoid();
292e5c89e4eSSatish Balay }
293e5c89e4eSSatish Balay #endif
294e5c89e4eSSatish Balay 
295570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
296d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
297d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
298d9822059SBarry Smith 
299d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
300d71ae5a4SJacob Faibussowitsch {
301d9822059SBarry Smith   PetscInt i, count = *cnt;
302d9822059SBarry Smith 
303d9822059SBarry Smith   PetscFunctionBegin;
3047c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3058c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
306a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]);
3077c2de775SJed Brown   }
3087c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3097c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3107c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
311ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3127c2de775SJed Brown   }
3137c2de775SJed Brown   #endif
3147c2de775SJed Brown   else {
3153ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
31641e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3178c764dc5SJose Roman   }
318d9822059SBarry Smith   PetscFunctionReturnVoid();
319d9822059SBarry Smith }
320d9822059SBarry Smith 
321d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
322d71ae5a4SJacob Faibussowitsch {
323d9822059SBarry Smith   PetscInt i, count = *cnt;
324d9822059SBarry Smith 
325d9822059SBarry Smith   PetscFunctionBegin;
3267c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3278c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
328a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]);
3297c2de775SJed Brown   }
3307c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3317c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3327c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
333ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3347c2de775SJed Brown   }
3357c2de775SJed Brown   #endif
3367c2de775SJed Brown   else {
3373ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"));
33841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3398c764dc5SJose Roman   }
340d9822059SBarry Smith   PetscFunctionReturnVoid();
341d9822059SBarry Smith }
342d9822059SBarry Smith #endif
343d9822059SBarry Smith 
344480cf27aSJed Brown /*
345480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
346480cf27aSJed Brown 
347ff0e51ddSBarry 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.
348480cf27aSJed Brown 
34912801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
350480cf27aSJed Brown 
351480cf27aSJed Brown */
352d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state)
353d71ae5a4SJacob Faibussowitsch {
35405342407SJunchao Zhang   PetscCommCounter      *counter = (PetscCommCounter *)count_val;
35557f21012SBarry Smith   struct PetscCommStash *comms   = counter->comms, *pcomm;
356480cf27aSJed Brown 
357480cf27aSJed Brown   PetscFunctionBegin;
3589566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm));
3599566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter->iflags));
36057f21012SBarry Smith   while (comms) {
3619566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&comms->comm));
36257f21012SBarry Smith     pcomm = comms;
36357f21012SBarry Smith     comms = comms->next;
3649566063dSJacob Faibussowitsch     PetscCall(PetscFree(pcomm));
36557f21012SBarry Smith   }
3669566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter));
367480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
368480cf27aSJed Brown }
369480cf27aSJed Brown 
370480cf27aSJed Brown /*
37147435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
372da3039f7SJed Brown   calls MPI_Comm_free().
373da3039f7SJed Brown 
374da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
375480cf27aSJed Brown 
376ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
377480cf27aSJed Brown 
37812801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
379480cf27aSJed Brown 
380480cf27aSJed Brown */
381d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
382d71ae5a4SJacob Faibussowitsch {
3839371c9d4SSatish Balay   union
384480cf27aSJed Brown   {
3859371c9d4SSatish Balay     MPI_Comm comm;
3869371c9d4SSatish Balay     void    *ptr;
3879371c9d4SSatish Balay   } icomm;
388480cf27aSJed Brown 
389480cf27aSJed Brown   PetscFunctionBegin;
39012801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval");
391ec4fadc2SJed Brown   icomm.ptr = attr_val;
39276bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
39333779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
39433779a13SJunchao Zhang     PetscMPIInt flg;
3959371c9d4SSatish Balay     union
3969371c9d4SSatish Balay     {
3979371c9d4SSatish Balay       MPI_Comm comm;
3989371c9d4SSatish Balay       void    *ptr;
3999371c9d4SSatish Balay     } ocomm;
4009566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg));
40133779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute");
40233779a13SJunchao 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");
40333779a13SJunchao Zhang   }
4049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval));
4059566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm));
406da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
407b89831e5SBarry Smith }
408da3039f7SJed Brown 
409da3039f7SJed Brown /*
41033779a13SJunchao 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.
411da3039f7SJed Brown  */
412d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
413d71ae5a4SJacob Faibussowitsch {
414da3039f7SJed Brown   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm));
416480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
417480cf27aSJed Brown }
418480cf27aSJed Brown 
41933779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *);
42042218b76SBarry Smith 
421951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
4228cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *);
4238cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
4248cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
425e39fd77fSBarry Smith #endif
426e39fd77fSBarry Smith 
4270f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE;
42812801b39SBarry Smith 
429eb27c7c8SSatish Balay PETSC_INTERN int    PetscGlobalArgc;
430eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4316ae9a8a6SBarry Smith int                 PetscGlobalArgc = 0;
43202c9f0b5SLisandro Dalcin char              **PetscGlobalArgs = NULL;
433dff31646SBarry Smith PetscSegBuffer      PetscCitationsList;
434e5c89e4eSSatish Balay 
435d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void)
436d71ae5a4SJacob Faibussowitsch {
437051e4cf2SJed Brown   PetscFunctionBegin;
4389566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList));
43930c35bf2SSatish Balay 
44030c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\
44130c35bf2SSatish Balay   Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\
44230c35bf2SSatish Balay     and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\
4433f81df79SBarry Smith     and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\
44430c35bf2SSatish Balay     and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\
44530c35bf2SSatish Balay     and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\
44630c35bf2SSatish Balay     and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\
44730c35bf2SSatish Balay     and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\
44830c35bf2SSatish Balay   Title = {{PETSc/TAO} Users Manual},\n\
44938aca504SSatish Balay   Number = {ANL-21/39 - Revision 3.20},\n\
4509a2cd68eSMatthew Knepley   Doi = {10.2172/1968587},\n\
45130c35bf2SSatish Balay   Institution = {Argonne National Laboratory},\n\
452c4dc7257SSatish Balay   Year = {2023}\n}\n",
4539371c9d4SSatish Balay                                    NULL));
45430c35bf2SSatish Balay 
45530c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\
45630c35bf2SSatish Balay   Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\
45730c35bf2SSatish Balay   Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\
45830c35bf2SSatish Balay   Booktitle = {Modern Software Tools in Scientific Computing},\n\
45930c35bf2SSatish Balay   Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\
46030c35bf2SSatish Balay   Pages = {163--202},\n\
46130c35bf2SSatish Balay   Publisher = {Birkh{\\\"{a}}user Press},\n\
4629371c9d4SSatish Balay   Year = {1997}\n}\n",
4639371c9d4SSatish Balay                                    NULL));
46430c35bf2SSatish Balay 
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
466051e4cf2SJed Brown }
467e5c89e4eSSatish Balay 
4682d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4692d747510SLisandro Dalcin 
470d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[])
471d71ae5a4SJacob Faibussowitsch {
4722d747510SLisandro Dalcin   PetscFunctionBegin;
4739566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(programname, name, sizeof(programname)));
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4752d747510SLisandro Dalcin }
4762d747510SLisandro Dalcin 
4772d747510SLisandro Dalcin /*@C
4782d747510SLisandro Dalcin   PetscGetProgramName - Gets the name of the running program.
4792d747510SLisandro Dalcin 
4802d747510SLisandro Dalcin   Not Collective
4812d747510SLisandro Dalcin 
4822d747510SLisandro Dalcin   Input Parameter:
4832d747510SLisandro Dalcin . len - length of the string name
4842d747510SLisandro Dalcin 
4852d747510SLisandro Dalcin   Output Parameter:
486811af0c4SBarry Smith . name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN`
4872d747510SLisandro Dalcin 
4882d747510SLisandro Dalcin   Level: advanced
4892d747510SLisandro Dalcin 
49021532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()`
4912d747510SLisandro Dalcin @*/
492d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProgramName(char name[], size_t len)
493d71ae5a4SJacob Faibussowitsch {
4942d747510SLisandro Dalcin   PetscFunctionBegin;
4959566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(name, programname, len));
4963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4972d747510SLisandro Dalcin }
4982d747510SLisandro Dalcin 
499e5c89e4eSSatish Balay /*@C
500e5c89e4eSSatish Balay   PetscGetArgs - Allows you to access the raw command line arguments anywhere
501811af0c4SBarry Smith   after PetscInitialize() is called but before `PetscFinalize()`.
502e5c89e4eSSatish Balay 
503e5c89e4eSSatish Balay   Not Collective
504e5c89e4eSSatish Balay 
505e5c89e4eSSatish Balay   Output Parameters:
506e5c89e4eSSatish Balay + argc - count of number of command line arguments
507e5c89e4eSSatish Balay - args - the command line arguments
508e5c89e4eSSatish Balay 
509e5c89e4eSSatish Balay   Level: intermediate
510e5c89e4eSSatish Balay 
511e5c89e4eSSatish Balay   Notes:
512e5c89e4eSSatish Balay   This is usually used to pass the command line arguments into other libraries
513e5c89e4eSSatish Balay   that are called internally deep in PETSc or the application.
514e5c89e4eSSatish Balay 
515f177e3b1SBarry Smith   The first argument contains the program name as is normal for C arguments.
516f177e3b1SBarry Smith 
51721532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()`
518e5c89e4eSSatish Balay @*/
519d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args)
520d71ae5a4SJacob Faibussowitsch {
521e5c89e4eSSatish Balay   PetscFunctionBegin;
52239a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
523e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
524e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
5253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
526e5c89e4eSSatish Balay }
527e5c89e4eSSatish Balay 
528793721a6SBarry Smith /*@C
529793721a6SBarry Smith   PetscGetArguments - Allows you to access the  command line arguments anywhere
530811af0c4SBarry Smith   after `PetscInitialize()` is called but before `PetscFinalize()`.
531793721a6SBarry Smith 
532793721a6SBarry Smith   Not Collective
533793721a6SBarry Smith 
5342fe279fdSBarry Smith   Output Parameter:
535793721a6SBarry Smith . args - the command line arguments
536793721a6SBarry Smith 
537793721a6SBarry Smith   Level: intermediate
538793721a6SBarry Smith 
53921532e8aSBarry Smith   Note:
54021532e8aSBarry Smith   This does NOT start with the program name and IS `NULL` terminated (final arg is void)
541793721a6SBarry Smith 
54221532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`, `PetscInitialize()`
543793721a6SBarry Smith @*/
544d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args)
545d71ae5a4SJacob Faibussowitsch {
546793721a6SBarry Smith   PetscInt i, argc = PetscGlobalArgc;
547793721a6SBarry Smith 
548793721a6SBarry Smith   PetscFunctionBegin;
54939a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
5509371c9d4SSatish Balay   if (!argc) {
5519371c9d4SSatish Balay     *args = NULL;
5523ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5539371c9d4SSatish Balay   }
5549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(argc, args));
5559566063dSJacob Faibussowitsch   for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i]));
5562d747510SLisandro Dalcin   (*args)[argc - 1] = NULL;
5573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
558793721a6SBarry Smith }
559793721a6SBarry Smith 
560793721a6SBarry Smith /*@C
561811af0c4SBarry Smith   PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()`
562793721a6SBarry Smith 
563793721a6SBarry Smith   Not Collective
564793721a6SBarry Smith 
5652fe279fdSBarry Smith   Output Parameter:
566793721a6SBarry Smith . args - the command line arguments
567793721a6SBarry Smith 
568793721a6SBarry Smith   Level: intermediate
569793721a6SBarry Smith 
570db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()`
571793721a6SBarry Smith @*/
572d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args)
573d71ae5a4SJacob Faibussowitsch {
57439a651e2SJacob Faibussowitsch   PetscFunctionBegin;
57539a651e2SJacob Faibussowitsch   if (args) {
576793721a6SBarry Smith     PetscInt i = 0;
577793721a6SBarry Smith 
5789566063dSJacob Faibussowitsch     while (args[i]) PetscCall(PetscFree(args[i++]));
5799566063dSJacob Faibussowitsch     PetscCall(PetscFree(args));
58039a651e2SJacob Faibussowitsch   }
5813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
582793721a6SBarry Smith }
583793721a6SBarry Smith 
58427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
58530befbd2SBarry Smith   #include <petscconfiginfo.h>
58630befbd2SBarry Smith 
587d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
588d71ae5a4SJacob Faibussowitsch {
58927104ee2SJacob Faibussowitsch   PetscFunctionBegin;
59011525c0dSBarry Smith   if (!PetscGlobalRank) {
59130befbd2SBarry Smith     char      cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64];
59211525c0dSBarry Smith     int       port;
593ffbd1cfbSBarry Smith     PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE;
59411525c0dSBarry Smith     size_t    applinelen, introlen;
595ffbd1cfbSBarry Smith     char      sawsurl[256];
59611525c0dSBarry Smith 
5979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg));
59811525c0dSBarry Smith     if (flg) {
59911525c0dSBarry Smith       char sawslog[PETSC_MAX_PATH_LEN];
60011525c0dSBarry Smith 
6019566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL));
60211525c0dSBarry Smith       if (sawslog[0]) {
603792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog));
60411525c0dSBarry Smith       } else {
605792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL));
60611525c0dSBarry Smith       }
60711525c0dSBarry Smith     }
6089566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg));
60948a46eb9SPierre Jolivet     if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert));
6109566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL));
611ffbd1cfbSBarry Smith     if (selectport) {
612792fecdfSBarry Smith       PetscCallSAWs(SAWs_Get_Available_Port, (&port));
613792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Port, (port));
614ffbd1cfbSBarry Smith     } else {
6159566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg));
61648a46eb9SPierre Jolivet       if (flg) PetscCallSAWs(SAWs_Set_Port, (port));
617ffbd1cfbSBarry Smith     }
6189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg));
61911525c0dSBarry Smith     if (flg) {
620792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Document_Root, (root));
6219566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(root, ".", &rootlocal));
6229c1e0ce8SBarry Smith     } else {
6239566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg));
6249c1e0ce8SBarry Smith       if (flg) {
6259566063dSJacob Faibussowitsch         PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root)));
626792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Document_Root, (root));
6279c1e0ce8SBarry Smith       }
62811525c0dSBarry Smith     }
6299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2));
63011525c0dSBarry Smith     if (flg2) {
63111525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
63228b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option");
6339566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root));
6349566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(jsdir, 'r', &flg));
63528b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory");
636792fecdfSBarry Smith       PetscCallSAWs(SAWs_Push_Local_Header, ());
63711525c0dSBarry Smith     }
6389566063dSJacob Faibussowitsch     PetscCall(PetscGetProgramName(programname, sizeof(programname)));
6399566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(help, &applinelen));
64011525c0dSBarry Smith     introlen = 4096 + applinelen;
64130a8c9c0SSurtai Han     applinelen += 1024;
6429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(applinelen, &appline));
6439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(introlen, &intro));
64411525c0dSBarry Smith 
64511525c0dSBarry Smith     if (rootlocal) {
6469566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname));
6479566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(appline, 'r', &rootlocal));
64811525c0dSBarry Smith     }
6499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetAll(NULL, &options));
65011525c0dSBarry Smith     if (rootlocal && help) {
6519566063dSJacob 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));
65211525c0dSBarry Smith     } else if (help) {
6539566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help));
65411525c0dSBarry Smith     } else {
6559566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options));
65611525c0dSBarry Smith     }
6579566063dSJacob Faibussowitsch     PetscCall(PetscFree(options));
6589566063dSJacob Faibussowitsch     PetscCall(PetscGetVersion(version, sizeof(version)));
6599371c9d4SSatish Balay     PetscCall(PetscSNPrintf(intro, introlen,
6609371c9d4SSatish Balay                             "<body>\n"
661a17b96a8SKyle 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"
662df62c222SBarry 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"
6639371c9d4SSatish Balay                             "%s",
6649371c9d4SSatish Balay                             version, petscconfigureoptions, appline));
665792fecdfSBarry Smith     PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro));
6669566063dSJacob Faibussowitsch     PetscCall(PetscFree(intro));
6679566063dSJacob Faibussowitsch     PetscCall(PetscFree(appline));
668ffbd1cfbSBarry Smith     if (selectport) {
669aa573868SBarry Smith       PetscBool silent;
6707d812c46SBarry Smith 
6717d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
67239a651e2SJacob Faibussowitsch       while (SAWs_Initialize()) {
673792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_Available_Port, (&port));
674792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port, (port));
6757d812c46SBarry Smith       }
6767d812c46SBarry Smith 
6779566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL));
678aa573868SBarry Smith       if (!silent) {
679792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl));
6809566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl));
681ffbd1cfbSBarry Smith       }
6827d812c46SBarry Smith     } else {
683792fecdfSBarry Smith       PetscCallSAWs(SAWs_Initialize, ());
684aa573868SBarry Smith     }
6859566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister("@TechReport{ saws,\n"
6860af79b04SBarry Smith                                      "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6870af79b04SBarry Smith                                      "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6880af79b04SBarry Smith                                      "  Institution = {Argonne National Laboratory},\n"
6899371c9d4SSatish Balay                                      "  Year   = 2013\n}\n",
6909371c9d4SSatish Balay                                      NULL));
69111525c0dSBarry Smith   }
6923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69311525c0dSBarry Smith }
69411525c0dSBarry Smith #endif
69511525c0dSBarry Smith 
6964dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
697d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
698d71ae5a4SJacob Faibussowitsch {
6994dfee713SSatish Balay   PetscFunctionBegin;
7004dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
7014dfee713SSatish Balay   /* see MPI.py for details on this bug */
7024dfee713SSatish Balay   (void)setenv("HWLOC_COMPONENTS", "-x86", 1);
7034dfee713SSatish Balay #endif
7043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7054dfee713SSatish Balay }
7064dfee713SSatish Balay 
707a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
708a56f64adSBarry Smith   #include <adios.h>
70922580e64SBarry Smith   #include <adios_read.h>
7107b56e58cSSatish Balay int64_t Petsc_adios_group;
711a56f64adSBarry Smith #endif
712a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
713cd1b99a6SBarry Smith   #include <omp.h>
714f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
715cd1b99a6SBarry Smith #endif
716a56f64adSBarry Smith 
717a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
718a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA)
7190e6b6b59SJacob Faibussowitsch   #include <petscdevice_cuda.h>
720a4af0ceeSJacob Faibussowitsch // REMOVE ME
721a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
722a4af0ceeSJacob Faibussowitsch #endif
723a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP)
7240e6b6b59SJacob Faibussowitsch   #include <petscdevice_hip.h>
725a4af0ceeSJacob Faibussowitsch // REMOVE ME
726a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
727a4af0ceeSJacob Faibussowitsch #endif
728a4af0ceeSJacob Faibussowitsch 
72927104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
730bc8a24ecSBarry Smith   #include <dlfcn.h>
731bc8a24ecSBarry Smith #endif
732a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
733a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
7343274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscViennaCLInit(void);
735a4af0ceeSJacob Faibussowitsch PetscBool                   PetscViennaCLSynchronize = PETSC_FALSE;
736a4af0ceeSJacob Faibussowitsch #endif
737bc8a24ecSBarry Smith 
738660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE;
739660278c0SBarry Smith 
74085649d77SJunchao Zhang /*
74185649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
74285649d77SJunchao Zhang 
74385649d77SJunchao Zhang   prog:     program name
74402101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
74585649d77SJunchao Zhang   help:     program help message
746da81f932SPierre Jolivet   ftn:      is it called from Fortran initialization (petscinitializef_)?
74785649d77SJunchao Zhang   readarguments,len: used when fortran is true
74885649d77SJunchao Zhang */
749d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len)
750d71ae5a4SJacob Faibussowitsch {
75185649d77SJunchao Zhang   PetscMPIInt size;
75285649d77SJunchao Zhang   PetscBool   flg = PETSC_TRUE;
75385649d77SJunchao Zhang   char        hostname[256];
75485649d77SJunchao Zhang 
75527104ee2SJacob Faibussowitsch   PetscFunctionBegin;
7563ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
75739a651e2SJacob Faibussowitsch   /* these must be initialized in a routine, not as a constant declaration */
75839a651e2SJacob Faibussowitsch   PETSC_STDOUT = stdout;
75939a651e2SJacob Faibussowitsch   PETSC_STDERR = stderr;
76039a651e2SJacob Faibussowitsch 
7619566063dSJacob Faibussowitsch   /* PetscCall can be used from now */
76239a651e2SJacob Faibussowitsch   PetscErrorHandlingInitialized = PETSC_TRUE;
76339a651e2SJacob Faibussowitsch 
76485649d77SJunchao Zhang   /*
76585649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
76685649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
76785649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
76885649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
76985649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
77085649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
77185649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
77285649d77SJunchao Zhang       listed above and since that time are compatible.
77385649d77SJunchao Zhang 
77485649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
77585649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
77685649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
77785649d77SJunchao Zhang       to perform the checking.
77885649d77SJunchao Zhang 
77985649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
78085649d77SJunchao Zhang 
78185649d77SJunchao Zhang       Questions:
78285649d77SJunchao Zhang 
78385649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
78485649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
78585649d77SJunchao Zhang   */
78685649d77SJunchao Zhang 
78785649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
78885649d77SJunchao Zhang   {
78985649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
79085649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
79139a651e2SJacob Faibussowitsch 
7929566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength));
79385649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
79485649d77SJunchao Zhang   #if defined(MPICH_VERSION)
79585649d77SJunchao Zhang     #if MPICH_NUMVERSION < 30100000
79685649d77SJunchao Zhang     {
79785649d77SJunchao Zhang       char     *ver, *lf;
79885649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
79939a651e2SJacob Faibussowitsch 
8009566063dSJacob Faibussowitsch       PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver));
80139a651e2SJacob Faibussowitsch       if (ver) {
8029566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(ver, '\n', &lf));
80339a651e2SJacob Faibussowitsch         if (lf) {
80485649d77SJunchao Zhang           *lf = 0;
8059566063dSJacob Faibussowitsch           PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg));
80685649d77SJunchao Zhang         }
80785649d77SJunchao Zhang       }
80885649d77SJunchao Zhang       if (!flg) {
809d5b396e8SSatish Balay         PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION));
81085649d77SJunchao Zhang         flg = PETSC_TRUE;
81185649d77SJunchao Zhang       }
81285649d77SJunchao Zhang     }
81385649d77SJunchao Zhang     #endif
81485649d77SJunchao 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?) */
81585649d77SJunchao Zhang   #elif defined(OMPI_MAJOR_VERSION)
81685649d77SJunchao Zhang     {
81785649d77SJunchao Zhang       char     *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf;
81885649d77SJunchao Zhang       PetscBool flg                                              = PETSC_FALSE;
81985649d77SJunchao Zhang     #define PSTRSZ 2
82085649d77SJunchao Zhang       char      ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"};
82185649d77SJunchao Zhang       char      ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "};
82285649d77SJunchao Zhang       int       i;
82385649d77SJunchao Zhang       for (i = 0; i < PSTRSZ; i++) {
8249566063dSJacob Faibussowitsch         PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver));
82539a651e2SJacob Faibussowitsch         if (ver) {
8269566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION));
8279566063dSJacob Faibussowitsch           PetscCall(PetscStrstr(ver, bs, &bsf));
82839a651e2SJacob Faibussowitsch           if (bsf) flg = PETSC_TRUE;
82985649d77SJunchao Zhang           break;
83085649d77SJunchao Zhang         }
83185649d77SJunchao Zhang       }
83285649d77SJunchao Zhang       if (!flg) {
8333ba16761SJacob 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));
83485649d77SJunchao Zhang         flg = PETSC_TRUE;
83585649d77SJunchao Zhang       }
83685649d77SJunchao Zhang     }
83785649d77SJunchao Zhang   #endif
83885649d77SJunchao Zhang   }
83985649d77SJunchao Zhang #endif
84085649d77SJunchao Zhang 
841e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__))
84285649d77SJunchao 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 */
84339a651e2SJacob 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");
84485649d77SJunchao Zhang #endif
84585649d77SJunchao Zhang 
84685649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
84785649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
84885649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
84985649d77SJunchao Zhang #endif
85085649d77SJunchao Zhang 
8519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
85285649d77SJunchao Zhang 
85385649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
85485649d77SJunchao Zhang 
8559566063dSJacob Faibussowitsch   PetscCall(PetscSetProgramName(prog));
8569566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen));
8579566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout));
8589566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr));
8599566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscCommSpinLock));
86085649d77SJunchao Zhang 
86185649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
8629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN));
86385649d77SJunchao Zhang 
86485649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
8659566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS));
8669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE));
86785649d77SJunchao Zhang   }
86885649d77SJunchao Zhang 
86985649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
8709566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
87185649d77SJunchao Zhang 
8729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank));
8739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize));
87485649d77SJunchao Zhang 
87585649d77SJunchao Zhang   MPIU_BOOL        = MPI_INT;
87685649d77SJunchao Zhang   MPIU_ENUM        = MPI_INT;
87785649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64;
87885649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
87985649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
88085649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
88185649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
88285649d77SJunchao Zhang #endif
88339a651e2SJacob Faibussowitsch   else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t");
88485649d77SJunchao Zhang 
88585649d77SJunchao Zhang     /*
88685649d77SJunchao Zhang      Initialized the global complex variable; this is because with
88785649d77SJunchao Zhang      shared libraries the constructors for global variables
88885649d77SJunchao Zhang      are not called; at least on IRIX.
88985649d77SJunchao Zhang   */
89085649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
89185649d77SJunchao Zhang   {
89285649d77SJunchao Zhang   #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
89385649d77SJunchao Zhang     PetscComplex ic(0.0, 1.0);
89485649d77SJunchao Zhang     PETSC_i = ic;
89585649d77SJunchao Zhang   #else
89685649d77SJunchao Zhang     PETSC_i = _Complex_I;
89785649d77SJunchao Zhang   #endif
89885649d77SJunchao Zhang   }
89985649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
90085649d77SJunchao Zhang 
90185649d77SJunchao Zhang   /*
90285649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
90385649d77SJunchao Zhang      half of the entries and maxes the second half.
90485649d77SJunchao Zhang   */
9059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP));
90685649d77SJunchao Zhang 
907613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
9089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128));
9099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128));
91074c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
9119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128));
9129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128));
91385649d77SJunchao Zhang   #endif
91474c01ffaSSatish Balay #endif
9159e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
9169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16));
9179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FP16));
91885649d77SJunchao Zhang #endif
91985649d77SJunchao Zhang 
92085649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
9219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM));
922613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX));
923613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN));
9249e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
9259e517322SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128));
92685649d77SJunchao Zhang #endif
92785649d77SJunchao Zhang 
9289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR));
92962e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp));
9309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR));
93185649d77SJunchao Zhang 
93285649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
93385649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
93485649d77SJunchao Zhang   {
93585649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
93693d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)};
93785649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_REAL, MPIU_INT}, tmpStruct;
93885649d77SJunchao Zhang 
9399566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
94093d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT));
9419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9429566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT));
94385649d77SJunchao Zhang   }
94485649d77SJunchao Zhang   {
94585649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
94693d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)};
94785649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_SCALAR, 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_scalar_int), &MPIU_SCALAR_INT));
9519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9529566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
95385649d77SJunchao Zhang   }
95485649d77SJunchao Zhang #endif
95585649d77SJunchao Zhang 
95685649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT));
9589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2INT));
95985649d77SJunchao Zhang #endif
9609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT));
9619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPI_4INT));
9629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT));
9639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_4INT));
96485649d77SJunchao Zhang 
96585649d77SJunchao Zhang   /*
96685649d77SJunchao Zhang      Attributes to be set on PETSc communicators
96785649d77SJunchao Zhang   */
9689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0));
9699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0));
9709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0));
9719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0));
97262e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, (void *)0));
97362e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, (void *)0));
97485649d77SJunchao Zhang 
975fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS)
9769566063dSJacob Faibussowitsch   if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len));
97785649d77SJunchao Zhang   else
97885649d77SJunchao Zhang #endif
9799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file));
98085649d77SJunchao Zhang 
98185649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9829566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
98385649d77SJunchao Zhang 
98485649d77SJunchao Zhang   /*
98585649d77SJunchao Zhang      Check system options and print help
98685649d77SJunchao Zhang   */
9879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCheckInitial_Private(help));
98885649d77SJunchao Zhang 
989a4af0ceeSJacob Faibussowitsch   /*
9900e6b6b59SJacob Faibussowitsch     Creates the logging data structures; this is enabled even if logging is not turned on
9910e6b6b59SJacob Faibussowitsch     This is the last thing we do before returning to the user code to prevent having the
9920e6b6b59SJacob Faibussowitsch     logging numbers contaminated by any startup time associated with MPI
9930e6b6b59SJacob Faibussowitsch   */
9940e6b6b59SJacob Faibussowitsch   PetscCall(PetscLogInitialize());
9950e6b6b59SJacob Faibussowitsch 
9960e6b6b59SJacob Faibussowitsch   /*
997a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
998a4af0ceeSJacob Faibussowitsch 
999a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
1000a4af0ceeSJacob Faibussowitsch    1. MPI initialized
1001a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
10020e6b6b59SJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up
10030e6b6b59SJacob Faibussowitsch       its own SIGSEV handler via the push/pop interface.
10040e6b6b59SJacob Faibussowitsch    4. Logging initialized
1005a4af0ceeSJacob Faibussowitsch   */
10069566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
1007a4af0ceeSJacob Faibussowitsch 
1008a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
1009a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
1010d75802c7SJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg));
10119566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL));
1012a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
10139566063dSJacob Faibussowitsch   PetscCall(PetscViennaCLInit());
1014a4af0ceeSJacob Faibussowitsch #endif
1015a4af0ceeSJacob Faibussowitsch 
10169566063dSJacob Faibussowitsch   PetscCall(PetscCitationsInitialize());
101785649d77SJunchao Zhang 
101885649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
10199566063dSJacob Faibussowitsch   PetscCall(PetscInitializeSAWs(ftn ? NULL : help));
102027104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
10219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg));
10229566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscStackViewSAWs());
102385649d77SJunchao Zhang #endif
102485649d77SJunchao Zhang 
102585649d77SJunchao Zhang   /*
102685649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
102785649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
102885649d77SJunchao Zhang   */
10299566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_DynamicLibraries());
103085649d77SJunchao Zhang 
10319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10329566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size));
103396dcfd6cSBarry Smith   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
10349566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname));
103585649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
103685649d77SJunchao Zhang   {
103785649d77SJunchao Zhang     PetscBool omp_view_flag;
103885649d77SJunchao Zhang     char     *threads = getenv("OMP_NUM_THREADS");
103985649d77SJunchao Zhang 
104085649d77SJunchao Zhang     if (threads) {
10419566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads));
104285649d77SJunchao Zhang       (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads);
104385649d77SJunchao Zhang     } else {
10442f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt)omp_get_max_threads();
10459566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads));
104685649d77SJunchao Zhang     }
1047d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys");
10489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg));
10499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag));
1050d0609cedSBarry Smith     PetscOptionsEnd();
105185649d77SJunchao Zhang     if (flg) {
10529566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads));
105385649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
105485649d77SJunchao Zhang     }
105548a46eb9SPierre Jolivet     if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads));
105685649d77SJunchao Zhang   }
105785649d77SJunchao Zhang #endif
105885649d77SJunchao Zhang 
105985649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
106085649d77SJunchao Zhang   /*
106185649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
106285649d77SJunchao Zhang 
106385649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
106485649d77SJunchao Zhang   */
10659566063dSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL));
106685649d77SJunchao Zhang #endif
106785649d77SJunchao Zhang 
106885649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10699566063dSJacob Faibussowitsch   PetscCall(PetscFPTCreate(10000));
107085649d77SJunchao Zhang #endif
107185649d77SJunchao Zhang 
107285649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
107385649d77SJunchao Zhang   {
107485649d77SJunchao Zhang     PetscViewer viewer;
10759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg));
107685649d77SJunchao Zhang     if (flg) {
10779566063dSJacob Faibussowitsch       PetscCall(PetscProcessPlacementView(viewer));
10789566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
107985649d77SJunchao Zhang     }
108085649d77SJunchao Zhang   }
108185649d77SJunchao Zhang #endif
108285649d77SJunchao Zhang 
108385649d77SJunchao Zhang   flg = PETSC_TRUE;
10849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL));
10859566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE));
108685649d77SJunchao Zhang 
108785649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10883ba16761SJacob Faibussowitsch   PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD);
10893ba16761SJacob Faibussowitsch   PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default);
10903ba16761SJacob Faibussowitsch   PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", "");
10913ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, "");
109285649d77SJunchao Zhang #endif
109385649d77SJunchao Zhang 
109485649d77SJunchao Zhang #if defined(__VALGRIND_H)
109585649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE;
109685649d77SJunchao Zhang   #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
10979566063dSJacob 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"));
109885649d77SJunchao Zhang   #endif
109985649d77SJunchao Zhang #endif
110085649d77SJunchao Zhang   /*
110185649d77SJunchao Zhang       Set flag that we are completely initialized
110285649d77SJunchao Zhang   */
110385649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
110485649d77SJunchao Zhang 
11059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg));
11069566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPythonInitialize(NULL, NULL));
1107f1f2ae84SBarry Smith 
1108f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1109f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin());
1110f1f2ae84SBarry 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");
11113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111285649d77SJunchao Zhang }
111385649d77SJunchao Zhang 
111410450e9eSJacob Faibussowitsch // "Unknown section 'Environmental Variables'"
111510450e9eSJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-section-header-unknown
1116e5c89e4eSSatish Balay /*@C
1117e5c89e4eSSatish Balay   PetscInitialize - Initializes the PETSc database and MPI.
1118811af0c4SBarry Smith   `PetscInitialize()` calls MPI_Init() if that has yet to be called,
1119e5c89e4eSSatish Balay   so this routine should always be called near the beginning of
1120e5c89e4eSSatish Balay   your program -- usually the very first line!
1121e5c89e4eSSatish Balay 
1122811af0c4SBarry Smith   Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set
1123e5c89e4eSSatish Balay 
1124e5c89e4eSSatish Balay   Input Parameters:
1125e5c89e4eSSatish Balay + argc - count of number of command line arguments
1126e5c89e4eSSatish Balay . args - the command line arguments
1127be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1128be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1129be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1130c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
11310298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message
1132e5c89e4eSSatish Balay 
1133811af0c4SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that
1134811af0c4SBarry Smith    communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a
1135811af0c4SBarry Smith    four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not,
1136811af0c4SBarry Smith    then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even
113705827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1138e5c89e4eSSatish Balay 
1139e5c89e4eSSatish Balay   Options Database Keys:
11407ca660e7SBarry Smith + -help [intro]                                       - prints help method for each option; if intro is given the program stops after printing the introductory help message
11417ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...]        - Starts program in debugger
1142e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11438a690491SBarry Smith . -on_error_emacs <machinename>                       - causes emacsclient to jump to error file
1144811af0c4SBarry Smith . -on_error_abort                                     - calls `abort()` when error detected (no traceback)
1145811af0c4SBarry Smith . -on_error_mpiabort                                  - calls `MPI_abort()` when error detected
11461af3601dSBarry Smith . -error_output_stdout                                - prints PETSc error messages to stdout instead of the default stderr
11478a690491SBarry Smith . -error_output_none                                  - does not print the error messages (but handles errors in the same way as if this was not called)
1148bf4d2887SBarry Smith . -debugger_ranks [rank1,rank2,...]                   - Indicates ranks to start in debugger
1149e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds)            - Pauses debugger
1150e5c89e4eSSatish Balay . -stop_for_debugger                                  - Print message on how to attach debugger manually to
1151e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
1152aee23540SBarry Smith . -malloc_dump                                        - prints a list of all unfreed memory at the end of the run
115392f119d6SBarry 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
1154811af0c4SBarry Smith . -malloc_view                                        - show a list of all allocated memory during `PetscFinalize()`
115592f119d6SBarry Smith . -malloc_view_threshold <t>                          - only list memory allocations of size greater than t with -malloc_view
1156608c71bfSMatthew G. Knepley . -malloc_requested_size                              - malloc logging will record the requested size rather than size after alignment
115792f119d6SBarry Smith . -fp_trap                                            - Stops on floating point exceptions
1158e5c89e4eSSatish Balay . -no_signal_handler                                  - Indicates not to trap error signals
1159e5c89e4eSSatish Balay . -shared_tmp                                         - indicates /tmp directory is shared by all processors
1160e5c89e4eSSatish Balay . -not_shared_tmp                                     - each processor has own /tmp
1161e5c89e4eSSatish Balay . -tmp                                                - alternative name of /tmp directory
1162e5c89e4eSSatish Balay . -get_total_flops                                    - returns total flops done by all processors
11630841954dSBarry Smith - -memory_view                                        - Print memory usage at end of run
1164e5c89e4eSSatish Balay 
1165c5b5d8d5SVaclav Hapla   Options Database Keys for Option Database:
1166c5b5d8d5SVaclav Hapla + -skip_petscrc           - skip the default option files ~/.petscrc, .petscrc, petscrc
1167c5b5d8d5SVaclav Hapla . -options_monitor        - monitor all set options to standard output for the whole program run
1168811af0c4SBarry Smith - -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()`
1169c5b5d8d5SVaclav Hapla 
1170c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1171c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1172c5b5d8d5SVaclav Hapla    They can be used also in option files.
1173c5b5d8d5SVaclav Hapla 
1174811af0c4SBarry Smith    See `PetscOptionsMonitorSet()` to do monitoring programmatically.
1175c5b5d8d5SVaclav Hapla 
1176e5c89e4eSSatish Balay   Options Database Keys for Profiling:
1177a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1178811af0c4SBarry Smith + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`.
1179217044c2SLisandro Dalcin . -log_sync                                            - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1180217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1181495fc317SBarry Smith . -log_trace [filename]                                - Print traces of all PETSc calls to the screen (useful to determine where a program
1182811af0c4SBarry Smith         hangs without running in the debugger).  See `PetscLogTraceBegin()`.
1183ad2e3d55SToby Isaac . -log_view [:filename:format][,[:filename:format]...] - Prints summary of flop and timing information to screen or file, see `PetscLogView()` (up to 4 viewers)
1184811af0c4SBarry Smith . -log_view_memory                                     - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`.
1185811af0c4SBarry Smith . -log_view_gpu_time                                   - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView().
1186f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes>                  - excludes subset of object classes from logging
1187b665b14eSToby Isaac . -log [filename]                                      - Logs profiling information in a dump file, see `PetscLogDump()`.
1188b665b14eSToby Isaac . -log_all [filename]                                  - Same as `-log`.
1189c2f74817SBarry Smith . -log_mpe [filename]                                  - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
11902611ad71SToby Isaac . -log_perfstubs                                       - Starts a log handler with the perfstubs interface (which is used by TAU)
119161cc7448SToby Isaac . -log_nvtx                                            - Starts an nvtx log handler for use with Nsight
1192811af0c4SBarry Smith . -viewfromoptions on,off                              - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off
1193c2f74817SBarry Smith - -check_pointer_intensity 0,1,2                       - if pointers are checked for validity (debug version only), using 0 will result in faster code
1194495fc317SBarry Smith 
1195ffbd1cfbSBarry Smith   Options Database Keys for SAWs:
1196ffbd1cfbSBarry Smith + -saws_port <portnumber>        - port number to publish SAWs data, default is 8080
1197ffbd1cfbSBarry 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
1198ffbd1cfbSBarry Smith                                    this is useful when you are running many jobs that utilize SAWs at the same time
1199ffbd1cfbSBarry Smith . -saws_log <filename>           - save a log of all SAWs communication
1200ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1201ffbd1cfbSBarry Smith - -saws_root <directory>         - allow SAWs to have access to the given directory to search for requested resources and files
1202ffbd1cfbSBarry Smith 
1203e5c89e4eSSatish Balay   Environmental Variables:
1204811af0c4SBarry Smith +   `PETSC_TMP` - alternative tmp directory
1205811af0c4SBarry Smith .   `PETSC_SHARED_TMP` - tmp is shared by all processes
1206811af0c4SBarry Smith .   `PETSC_NOT_SHARED_TMP` - each process has its own private tmp
1207811af0c4SBarry Smith .   `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs
1208811af0c4SBarry Smith .   `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1209811af0c4SBarry Smith .   `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer
1210811af0c4SBarry Smith -   `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to
1211e5c89e4eSSatish Balay 
1212e5c89e4eSSatish Balay   Level: beginner
1213e5c89e4eSSatish Balay 
1214811af0c4SBarry Smith   Note:
1215811af0c4SBarry Smith   If for some reason you must call `MPI_Init()` separately, call
1216811af0c4SBarry Smith   it before `PetscInitialize()`.
1217e5c89e4eSSatish Balay 
1218811af0c4SBarry Smith   Fortran Notes:
1219811af0c4SBarry Smith   In Fortran this routine can be called with
1220811af0c4SBarry Smith .vb
1221811af0c4SBarry Smith        call PetscInitialize(ierr)
1222811af0c4SBarry Smith        call PetscInitialize(file,ierr) or
1223811af0c4SBarry Smith        call PetscInitialize(file,help,ierr)
1224811af0c4SBarry Smith .ve
1225e5c89e4eSSatish Balay 
1226811af0c4SBarry Smith   If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after
1227811af0c4SBarry Smith   calling `PetscInitialize()`.
1228e5c89e4eSSatish Balay 
12295dedd997SBarry Smith   Options Database Key for Developers:
12305dedd997SBarry Smith . -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form:
12315dedd997SBarry Smith     "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py)
12325dedd997SBarry Smith 
1233db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()`
1234e5c89e4eSSatish Balay @*/
1235d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[])
1236d71ae5a4SJacob Faibussowitsch {
123785649d77SJunchao Zhang   PetscMPIInt flag;
123821ef0414SBarry Smith   const char *prog = "Unknown Name", *mpienv;
1239e5c89e4eSSatish Balay 
124027104ee2SJacob Faibussowitsch   PetscFunctionBegin;
12413ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
12429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Initialized(&flag));
1243e5c89e4eSSatish Balay   if (!flag) {
124439a651e2SJacob 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");
12459566063dSJacob Faibussowitsch     PetscCall(PetscPreMPIInit_Private());
12465e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12475e765c61SJed Brown     {
1248*835d5ba2SJunchao Zhang       PetscMPIInt provided;
1249*835d5ba2SJunchao Zhang       PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE ? MPI_THREAD_FUNNELED : PETSC_MPI_THREAD_REQUIRED, &provided));
1250*835d5ba2SJunchao Zhang       PetscCheck(PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE || provided >= PETSC_MPI_THREAD_REQUIRED, PETSC_COMM_SELF, PETSC_ERR_MPI, "The MPI implementation's provided thread level is less than what you required");
1251*835d5ba2SJunchao Zhang       if (PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE) PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED; // assign it a valid value after check-up
12525e765c61SJed Brown     }
12535e765c61SJed Brown #else
12549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Init(argc, args));
12555e765c61SJed Brown #endif
125621ef0414SBarry Smith     if (PetscDefined(HAVE_MPIUNI)) {
125721ef0414SBarry Smith       mpienv = getenv("PMI_SIZE");
125821ef0414SBarry Smith       if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE");
125921ef0414SBarry Smith       if (mpienv) {
126021ef0414SBarry Smith         PetscInt isize;
126121ef0414SBarry Smith         PetscCall(PetscOptionsStringToInt(mpienv, &isize));
126221ef0414SBarry 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");
126321ef0414SBarry 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");
126421ef0414SBarry Smith       }
126521ef0414SBarry Smith     }
1266e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1267e5c89e4eSSatish Balay   }
12684dfee713SSatish Balay 
126985649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1270e5c89e4eSSatish Balay   if (argc && args) {
1271e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1272e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1273e5c89e4eSSatish Balay   }
1274811af0c4SBarry Smith   PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0));
12753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1276e5c89e4eSSatish Balay }
1277e5c89e4eSSatish Balay 
127895c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1279ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsCounts;
1280ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsMaxCounts;
128195c0884eSLisandro Dalcin PETSC_INTERN PetscBool    PetscObjectsLog;
1282e5c89e4eSSatish Balay 
1283008a6e76SBarry Smith /*
1284008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1285008a6e76SBarry Smith */
1286d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void)
1287d71ae5a4SJacob Faibussowitsch {
1288008a6e76SBarry Smith   PetscFunctionBegin;
1289613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
12909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128));
129174c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
12929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128));
1293008a6e76SBarry Smith   #endif
129474c01ffaSSatish Balay #endif
12959e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
12969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FP16));
1297008a6e76SBarry Smith #endif
1298008a6e76SBarry Smith 
1299de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
13009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_SUM));
1301613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MAX));
1302613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MIN));
13039e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
13049e517322SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128));
1305008a6e76SBarry Smith #endif
1306008a6e76SBarry Smith 
13079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR));
13089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT));
13099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT));
131040df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
13119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2INT));
1312008a6e76SBarry Smith #endif
13139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPI_4INT));
13149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_4INT));
13159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
131662e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp));
13173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1318008a6e76SBarry Smith }
1319008a6e76SBarry Smith 
1320a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1321a4af0ceeSJacob Faibussowitsch 
1322e5c89e4eSSatish Balay /*@C
1323e5c89e4eSSatish Balay   PetscFinalize - Checks for options to be called at the conclusion
1324811af0c4SBarry Smith   of the program. `MPI_Finalize()` is called only if the user had not
1325811af0c4SBarry Smith   called `MPI_Init()` before calling `PetscInitialize()`.
1326e5c89e4eSSatish Balay 
1327811af0c4SBarry Smith   Collective on `PETSC_COMM_WORLD`
1328e5c89e4eSSatish Balay 
1329e5c89e4eSSatish Balay   Options Database Keys:
1330811af0c4SBarry Smith + -options_view                    - Calls `PetscOptionsView()`
1331e5c89e4eSSatish Balay . -options_left                    - Prints unused options that remain in the database
13327eb1d149SBarry 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
1333e5c89e4eSSatish Balay . -mpidump                         - Calls PetscMPIDump()
1334811af0c4SBarry Smith . -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed
13354bd3d7f8SBarry Smith . -memory_view                     - Prints total memory usage
13364bd3d7f8SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and in what functions
1337e5c89e4eSSatish Balay 
1338e5c89e4eSSatish Balay   Level: beginner
1339e5c89e4eSSatish Balay 
1340e5c89e4eSSatish Balay   Note:
1341811af0c4SBarry Smith   See `PetscInitialize()` for other runtime options.
1342e5c89e4eSSatish Balay 
1343db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()`
1344e5c89e4eSSatish Balay @*/
1345d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void)
1346d71ae5a4SJacob Faibussowitsch {
13474bb5149bSJed Brown   PetscMPIInt rank;
1348a8d2bbe5SBarry Smith   PetscInt    nopt;
13492bf49c77SBarry Smith   PetscBool   flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE;
1350dff31646SBarry Smith   PetscBool   flg;
135110463e74SBarry Smith   char        mname[PETSC_MAX_PATH_LEN];
1352e5c89e4eSSatish Balay 
13533cbbc5ffSBarry Smith   PetscFunctionBegin;
135439a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()");
13559566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PetscFinalize() called\n"));
1356b022a5c1SBarry Smith 
1357f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1358f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd());
1359f1f2ae84SBarry Smith 
136062e5d2d2SJDBetteridge   /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */
136162e5d2d2SJDBetteridge   {
136262e5d2d2SJDBetteridge     union
136362e5d2d2SJDBetteridge     {
136462e5d2d2SJDBetteridge       MPI_Comm comm;
136562e5d2d2SJDBetteridge       void    *ptr;
136662e5d2d2SJDBetteridge     } ucomm;
136762e5d2d2SJDBetteridge     PetscMPIInt flg;
136862e5d2d2SJDBetteridge     void       *tmp;
136962e5d2d2SJDBetteridge 
137062e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
137162e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
137262e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF));
137362e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
137462e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
137562e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD));
137662e5d2d2SJDBetteridge   }
137762e5d2d2SJDBetteridge 
13789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1379a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13803ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE);
13813ba16761SJacob Faibussowitsch   PetscCallExternal(adios_finalize, rank);
1382a56f64adSBarry Smith #endif
13839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg));
1384dff31646SBarry Smith   if (flg) {
13851f817a21SBarry Smith     char *cits, filename[PETSC_MAX_PATH_LEN];
13861f817a21SBarry Smith     FILE *fd = PETSC_STDOUT;
13871f817a21SBarry Smith 
13889566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL));
138948a46eb9SPierre Jolivet     if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd));
13909566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits));
1391dff31646SBarry Smith     cits[0] = 0;
13929566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits));
13939566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n"));
13949566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
13959566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits));
13969566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
13979566063dSJacob Faibussowitsch     PetscCall(PetscFClose(PETSC_COMM_WORLD, fd));
13989566063dSJacob Faibussowitsch     PetscCall(PetscFree(cits));
1399dff31646SBarry Smith   }
14009566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&PetscCitationsList));
1401dff31646SBarry Smith 
14022d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14039566063dSJacob Faibussowitsch   PetscCall(PetscFPTDestroy());
14042d53ad75SBarry Smith #endif
14052d53ad75SBarry Smith 
1406e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1407dff31646SBarry Smith   flg = PETSC_FALSE;
14089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL));
14091baa6e33SBarry Smith   if (flg) PetscCall(PetscOptionsSAWsDestroy());
1410d5649816SBarry Smith #endif
1411d5649816SBarry Smith 
1412681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1413681455b2SBarry Smith   flg1 = PETSC_FALSE;
14149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL));
1415681455b2SBarry Smith   if (flg1) {
1416681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
14179566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL));
1418681455b2SBarry Smith   }
1419681455b2SBarry Smith #endif
1420681455b2SBarry Smith 
142167584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL));
142348a46eb9SPierre Jolivet   if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n"));
142467584ceeSBarry Smith #endif
1425e5c89e4eSSatish Balay 
14262611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
142790d69ab7SBarry Smith     flg1 = PETSC_FALSE;
14289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL));
1429e5c89e4eSSatish Balay     if (flg1) {
1430e5c89e4eSSatish Balay       PetscLogDouble flops = 0;
14319566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD));
14329566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops));
1433e5c89e4eSSatish Balay     }
14342611ad71SToby Isaac   }
1435e5c89e4eSSatish Balay 
14362611ad71SToby Isaac   if (PetscDefined(USE_LOG) && PetscDefined(HAVE_MPE)) {
1437e5c89e4eSSatish Balay     mname[0] = 0;
14389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1));
14392611ad71SToby Isaac     if (flg1) PetscCall(PetscLogMPEDump(mname[0] ? mname : NULL));
1440e5c89e4eSSatish Balay   }
1441a297a907SKarl Rupp 
14423048253cSJacob Faibussowitsch   // Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
14439566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1444dd710f27SBarry Smith 
14452611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
14469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE));
14479566063dSJacob Faibussowitsch     PetscCall(PetscLogViewFromOptions());
14489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsPopGetViewerOff());
1449473903fcSJunchao Zhang     //  It should be turned on with PetscLogGpuTime() and never turned off except in this place
1450473903fcSJunchao Zhang     PetscLogGpuTimeFlag = PETSC_FALSE;
145187c3beb6SLisandro Dalcin 
14523048253cSJacob Faibussowitsch     // Free any objects created by the last block of code.
14539566063dSJacob Faibussowitsch     PetscCall(PetscObjectRegisterDestroyAll());
1454dd710f27SBarry Smith 
1455dd710f27SBarry Smith     mname[0] = 0;
14569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1));
14579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2));
14589566063dSJacob Faibussowitsch     if (flg1 || flg2) PetscCall(PetscLogDump(mname));
14592611ad71SToby Isaac   }
146010463e74SBarry Smith 
146190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL));
14639566063dSJacob Faibussowitsch   if (!flg1) PetscCall(PetscPopSignalHandler());
146490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL));
14661baa6e33SBarry Smith   if (flg1) PetscCall(PetscMPIDump(stdout));
146790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
146890d69ab7SBarry Smith   flg2 = PETSC_FALSE;
14698bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
14709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1));
14719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
14729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL));
1473e4c476e2SSatish Balay 
1474e5c89e4eSSatish Balay   if (flg2) {
1475be56827dSJed Brown     PetscViewer viewer;
14769566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
14779566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
14789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsView(NULL, viewer));
14799566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1480e5c89e4eSSatish Balay   }
1481e5c89e4eSSatish Balay 
1482e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
14839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1));
14849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1));
1485e5c89e4eSSatish Balay 
148633fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
14879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1));
148859f199a7SJed Brown   if (!flg1) flg3 = PETSC_TRUE;
1489e5c89e4eSSatish Balay   if (flg3) {
14903de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1491be56827dSJed Brown       PetscViewer viewer;
14929566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
14939566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
14949566063dSJacob Faibussowitsch       PetscCall(PetscOptionsView(NULL, viewer));
14959566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
1496e5c89e4eSSatish Balay     }
14979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsAllUsed(NULL, &nopt));
14983de2bfdfSBarry Smith     if (nopt) {
14999566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n"));
15009566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n"));
15013de2bfdfSBarry Smith       if (nopt == 1) {
15029566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n"));
1503e5c89e4eSSatish Balay       } else {
15049566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt));
1505e5c89e4eSSatish Balay       }
15063de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15079566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n"));
1508df12ba86SBarry Smith     }
15099566063dSJacob Faibussowitsch     PetscCall(PetscOptionsLeft(NULL));
1510e5c89e4eSSatish Balay   }
1511e5c89e4eSSatish Balay 
1512e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1513d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15149566063dSJacob Faibussowitsch     PetscCall(PetscStackSAWsViewOff());
1515792fecdfSBarry Smith     PetscCallSAWs(SAWs_Finalize, ());
1516d45a07a7SBarry Smith   }
1517ec957eceSBarry Smith #endif
1518ec957eceSBarry Smith 
151910463e74SBarry Smith   /*
1520dbc8283eSBarry Smith        List all objects the user may have forgot to free
15212eff7a51SBarry Smith   */
15222611ad71SToby Isaac   if (PetscDefined(USE_LOG) && PetscObjectsLog) {
15239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
1524a64a8e02SBarry Smith     if (flg1) {
1525a64a8e02SBarry Smith       MPI_Comm local_comm;
15267eb1d149SBarry Smith       char     string[64];
1527a64a8e02SBarry Smith 
15289566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL));
1529afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
15309566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
15319566063dSJacob Faibussowitsch       PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
15329566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
15339566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
15340a1571b3SBarry Smith     }
153505df10baSBarry Smith   }
15364097062eSBarry Smith 
1537dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1538dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
15399566063dSJacob Faibussowitsch   PetscCall(PetscFree(PetscObjects));
15402eff7a51SBarry Smith 
154133f85c2fSBarry Smith   /*
154233f85c2fSBarry Smith      Destroy any packages that registered a finalize
154333f85c2fSBarry Smith   */
15449566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalizeAll());
154533f85c2fSBarry Smith 
15469566063dSJacob Faibussowitsch   PetscCall(PetscLogFinalize());
1547101409b8SToby Isaac 
154833f85c2fSBarry Smith   /*
154948dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
155048dd1dffSBarry Smith   */
15512e956fe4SStefano Zampini   if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll());
155237e93019SBarry Smith 
15534028d114SSatish Balay   if (petsc_history) {
15549566063dSJacob Faibussowitsch     PetscCall(PetscCloseHistoryFile(&petsc_history));
155502c9f0b5SLisandro Dalcin     petsc_history = NULL;
1556e5c89e4eSSatish Balay   }
15579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
15589566063dSJacob Faibussowitsch   PetscCall(PetscInfoDestroy());
1559e5c89e4eSSatish Balay 
156067584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
156192f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1562e5c89e4eSSatish Balay     char  fname[PETSC_MAX_PATH_LEN];
156392f119d6SBarry Smith     char  sname[PETSC_MAX_PATH_LEN];
1564e5c89e4eSSatish Balay     FILE *fd;
1565ed9cf6e9SBarry Smith     int   err;
1566e5c89e4eSSatish Balay 
1567dc92acbaSJed Brown     flg2 = PETSC_FALSE;
156892f119d6SBarry Smith     flg3 = PETSC_FALSE;
15699566063dSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL));
15709566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL));
157192f119d6SBarry Smith     fname[0] = 0;
15729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1));
1573e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
15743ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
15759371c9d4SSatish Balay       fd = fopen(sname, "w");
15769371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
15779566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(fd));
1578ed9cf6e9SBarry Smith       err = fclose(fd);
157928b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
158092f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1581e5c89e4eSSatish Balay       MPI_Comm local_comm;
1582e5c89e4eSSatish Balay 
1583afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
15849566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
15859566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(stdout));
15869566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
15879566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1588e5c89e4eSSatish Balay     }
1589e5c89e4eSSatish Balay     fname[0] = 0;
15909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1));
1591e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
15923ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
15939371c9d4SSatish Balay       fd = fopen(sname, "w");
15949371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
15959566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(fd));
1596ed9cf6e9SBarry Smith       err = fclose(fd);
159728b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
159892f119d6SBarry Smith     } else if (flg1) {
159992f119d6SBarry Smith       MPI_Comm local_comm;
160092f119d6SBarry Smith 
1601afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16029566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16039566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(stdout));
16049566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16059566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1606e5c89e4eSSatish Balay     }
1607e5c89e4eSSatish Balay   }
160867584ceeSBarry Smith #endif
160920e2c332SMatthew G. Knepley 
16105486ca60SMatthew G. Knepley   /*
16115486ca60SMatthew G. Knepley      Close any open dynamic libraries
16125486ca60SMatthew G. Knepley   */
16139566063dSJacob Faibussowitsch   PetscCall(PetscFinalize_DynamicLibraries());
16145486ca60SMatthew G. Knepley 
1615e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroyDefault());
1617e5c89e4eSSatish Balay 
1618e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
161902c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1620e5c89e4eSSatish Balay 
1621c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1622c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
16239566063dSJacob Faibussowitsch     PetscCall(PetscKokkosFinalize_Private());
1624c2b86a48SJunchao Zhang     PetscBeganKokkos       = PETSC_FALSE;
162545639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1626c2b86a48SJunchao Zhang   }
1627c2b86a48SJunchao Zhang #endif
1628c2b86a48SJunchao Zhang 
162971438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
163071438e86SJunchao Zhang   if (PetscBeganNvshmem) {
16319566063dSJacob Faibussowitsch     PetscCall(PetscNvshmemFinalize());
163271438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
163371438e86SJunchao Zhang   }
163471438e86SJunchao Zhang #endif
1635a0e72f99SJunchao Zhang 
16369566063dSJacob Faibussowitsch   PetscCall(PetscFreeMPIResources());
1637e5c89e4eSSatish Balay 
1638dbc8283eSBarry Smith   /*
1639efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1640efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1641efb80d3cSBarry Smith 
1642efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1643efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1644dbc8283eSBarry Smith  */
1645b770b1f6SSatish Balay   {
1646dbc8283eSBarry Smith     PetscCommCounter *counter;
1647dbc8283eSBarry Smith     PetscMPIInt       flg;
1648dbc8283eSBarry Smith     MPI_Comm          icomm;
16499371c9d4SSatish Balay     union
16509371c9d4SSatish Balay     {
16519371c9d4SSatish Balay       MPI_Comm comm;
16529371c9d4SSatish Balay       void    *ptr;
16539371c9d4SSatish Balay     } ucomm;
16549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
1655dbc8283eSBarry Smith     if (flg) {
1656265f3f35SJed Brown       icomm = ucomm.comm;
16579566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
165828b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1659dbc8283eSBarry Smith 
16609566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval));
16619566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
16629566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1663dbc8283eSBarry Smith     }
16649566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
1665dbc8283eSBarry Smith     if (flg) {
1666265f3f35SJed Brown       icomm = ucomm.comm;
16679566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
166828b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1669dbc8283eSBarry Smith 
16709566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval));
16719566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
16729566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1673dbc8283eSBarry Smith     }
1674b770b1f6SSatish Balay   }
1675dbc8283eSBarry Smith 
16769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
16779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
16789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
16799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
168062e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval));
168162e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval));
1682480cf27aSJed Brown 
16835ea2b939SDuncan Campbell   // Free keyvals which may be silently created by some routines
16845ea2b939SDuncan Campbell   if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval));
16855ea2b939SDuncan Campbell   if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval));
16865ea2b939SDuncan Campbell 
16879566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
16889566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
16899566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
16909566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock));
1691ef19f930SBarry Smith 
1692e5c89e4eSSatish Balay   if (PetscBeganMPI) {
169399b1327fSBarry Smith     PetscMPIInt flag;
16949566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalized(&flag));
169539a651e2SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
169639a651e2SJacob Faibussowitsch     /* wait until the very last moment to disable error handling */
169739a651e2SJacob Faibussowitsch     PetscErrorHandlingInitialized = PETSC_FALSE;
16989566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalize());
169939a651e2SJacob Faibussowitsch   } else PetscErrorHandlingInitialized = PETSC_FALSE;
170039a651e2SJacob Faibussowitsch 
1701e5c89e4eSSatish Balay   /*
1702e5c89e4eSSatish Balay 
1703e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1704e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1705e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1706e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1707e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17080e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1709e5c89e4eSSatish Balay    memory was not freed.
1710e5c89e4eSSatish Balay 
1711e5c89e4eSSatish Balay */
17129566063dSJacob Faibussowitsch   PetscCall(PetscMallocClear());
17139566063dSJacob Faibussowitsch   PetscCall(PetscStackReset());
1714a297a907SKarl Rupp 
1715e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1716e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
17177ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE)
171856883f60SBarry Smith   /*
171956883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
172056883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
172156883f60SBarry Smith    */
172256883f60SBarry Smith   __gcov_flush();
172356883f60SBarry Smith #endif
17241724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
17251724198aSStefano Zampini   PetscStackClearTop;
17263ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1727e5c89e4eSSatish Balay }
1728e5c89e4eSSatish Balay 
172943db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
1730d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b)
1731d71ae5a4SJacob Faibussowitsch {
173243db4dbbSBarry Smith   if (*a == *b) return 1;
173343db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
173443db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
173543db4dbbSBarry Smith   return 0;
173643db4dbbSBarry Smith }
1737a70650f6SBarry Smith #endif
173843db4dbbSBarry Smith 
173943db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
1740d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b)
1741d71ae5a4SJacob Faibussowitsch {
174243db4dbbSBarry Smith   if (*a == *b) return 1;
174343db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
174443db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
174543db4dbbSBarry Smith   return 0;
174643db4dbbSBarry Smith }
174743db4dbbSBarry Smith #endif
1748