xref: /petsc/src/sys/objects/pinit.c (revision e8953f014f11b2d2ab249b0e93ff2da3d35372f1)
1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I  "petscsys.h"   I*/
6665c2dedSJed Brown #include <petscviewer.h>
7a0e72f99SJunchao Zhang 
86edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
96edef35eSSatish Balay   #include <petsc/private/valgrind/valgrind.h>
106edef35eSSatish Balay #endif
116edef35eSSatish Balay 
1285649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
1385649d77SJunchao Zhang   #include <petsc/private/fortranimpl.h>
1485649d77SJunchao Zhang #endif
1585649d77SJunchao Zhang 
1656883f60SBarry Smith #if defined(PETSC_USE_GCOV)
1756883f60SBarry Smith EXTERN_C_BEGIN
18aaf3846cSSatish Balay   #if defined(PETSC_HAVE___GCOV_DUMP)
19aaf3846cSSatish Balay     #define __gcov_flush(x) __gcov_dump(x)
20aaf3846cSSatish Balay   #endif
2156883f60SBarry Smith void __gcov_flush(void);
2256883f60SBarry Smith EXTERN_C_END
2356883f60SBarry Smith #endif
248101f56cSMatthew Knepley 
252d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
2695c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
272d53ad75SBarry Smith PetscFPT              PetscFPTData = 0;
282d53ad75SBarry Smith #endif
292d53ad75SBarry Smith 
3027104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
3116ad0300SBarry Smith   #include <petscviewersaws.h>
32a6790183SBarry Smith #endif
33eb02082bSJunchao Zhang 
34e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
35e5c89e4eSSatish Balay 
3695c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
37e5c89e4eSSatish Balay 
3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm, int);
4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm, int);
4295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE **);
430069ddf5SShri Abhyankar 
446de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
45e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
476de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
486de5d289SStefano Zampini #else
496de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
506de5d289SStefano Zampini #endif
51e5c89e4eSSatish Balay 
52480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
53480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
54480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
555f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
56480cf27aSJed Brown 
57e5c89e4eSSatish Balay /*
58e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
59e5c89e4eSSatish Balay */
6002c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL};
6102c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL};
62e5c89e4eSSatish Balay 
63ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
64ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
650f8e0872SSatish Balay 
66a2f94806SJed Brown PetscInt PetscHotRegionDepth;
67a2f94806SJed Brown 
686edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
696edef35eSSatish Balay 
70b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
71b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
72b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
73b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
74b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
75b22622e2STadeu Manoel #endif
76b22622e2STadeu Manoel 
77e5c89e4eSSatish Balay /*
78945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
7972a42c3cSBarry Smith 
8072a42c3cSBarry Smith    Collective
8172a42c3cSBarry Smith 
8272a42c3cSBarry Smith    Level: advanced
8372a42c3cSBarry Smith 
8495452b02SPatrick Sanan     Notes:
85a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
860f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
87a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
880f11a792SBarry Smith 
89a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
901ea5a559SBarry Smith 
91db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()`
920f11a792SBarry Smith */
93d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoPointers(int argc, char **args, const char *filename, const char *help)
94d71ae5a4SJacob Faibussowitsch {
9572a42c3cSBarry Smith   int    myargc = argc;
9672a42c3cSBarry Smith   char **myargs = args;
9772a42c3cSBarry Smith 
9872a42c3cSBarry Smith   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&myargc, &myargs, filename, help));
1009566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
101df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
10227104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
10372a42c3cSBarry Smith }
10472a42c3cSBarry Smith 
105f0865b08SBarry Smith /*
106a56f64adSBarry Smith       Used by Julia interface to get communicator
107f0865b08SBarry Smith */
108d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
109d71ae5a4SJacob Faibussowitsch {
110f0865b08SBarry Smith   PetscFunctionBegin;
11127104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidPointer(comm, 1);
112f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
113f0865b08SBarry Smith   PetscFunctionReturn(0);
114f0865b08SBarry Smith }
115f0865b08SBarry Smith 
116e5c89e4eSSatish Balay /*@C
117811af0c4SBarry Smith       PetscInitializeNoArguments - Calls `PetscInitialize()` from C/C++ without
118e5c89e4eSSatish Balay         the command line arguments.
119e5c89e4eSSatish Balay 
120e5c89e4eSSatish Balay    Collective
121e5c89e4eSSatish Balay 
122e5c89e4eSSatish Balay    Level: advanced
123e5c89e4eSSatish Balay 
124db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`
125e5c89e4eSSatish Balay @*/
126d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoArguments(void)
127d71ae5a4SJacob Faibussowitsch {
128e5c89e4eSSatish Balay   int    argc = 0;
12902c9f0b5SLisandro Dalcin   char **args = NULL;
130e5c89e4eSSatish Balay 
131e5c89e4eSSatish Balay   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, NULL, NULL));
13339a651e2SJacob Faibussowitsch   PetscFunctionReturn(0);
134e5c89e4eSSatish Balay }
135e5c89e4eSSatish Balay 
136e5c89e4eSSatish Balay /*@
137e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
138e5c89e4eSSatish Balay 
13993b6d2d1SJed Brown    Level: beginner
140e5c89e4eSSatish Balay 
141db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
142e5c89e4eSSatish Balay @*/
143d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialized(PetscBool *isInitialized)
144d71ae5a4SJacob Faibussowitsch {
14527104ee2SJacob Faibussowitsch   PetscFunctionBegin;
14627104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized, 1);
147e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
14827104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
149e5c89e4eSSatish Balay }
150e5c89e4eSSatish Balay 
151e5c89e4eSSatish Balay /*@
152811af0c4SBarry Smith       PetscFinalized - Determine whether `PetscFinalize()` has been called yet
153e5c89e4eSSatish Balay 
154e5c89e4eSSatish Balay    Level: developer
155e5c89e4eSSatish Balay 
156db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
157e5c89e4eSSatish Balay @*/
158d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalized(PetscBool *isFinalized)
159d71ae5a4SJacob Faibussowitsch {
16027104ee2SJacob Faibussowitsch   PetscFunctionBegin;
16127104ee2SJacob Faibussowitsch   if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized, 1);
162e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
16327104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
164e5c89e4eSSatish Balay }
165e5c89e4eSSatish Balay 
16657171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]);
167e5c89e4eSSatish Balay 
168e5c89e4eSSatish Balay /*
169e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
170e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
171e5c89e4eSSatish Balay */
172367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
173e5c89e4eSSatish Balay 
174d71ae5a4SJacob Faibussowitsch PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, int *cnt, MPI_Datatype *datatype)
175d71ae5a4SJacob Faibussowitsch {
176e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out, i, count = *cnt;
177e5c89e4eSSatish Balay 
178e5c89e4eSSatish Balay   PetscFunctionBegin;
179e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
180e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
18141e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
182e5c89e4eSSatish Balay   }
183e5c89e4eSSatish Balay 
184e5c89e4eSSatish Balay   for (i = 0; i < count; i++) {
185e5c89e4eSSatish Balay     xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]);
186e5c89e4eSSatish Balay     xout[2 * i + 1] += xin[2 * i + 1];
187e5c89e4eSSatish Balay   }
188812af9f3SBarry Smith   PetscFunctionReturnVoid();
189e5c89e4eSSatish Balay }
190e5c89e4eSSatish Balay 
191e5c89e4eSSatish Balay /*
192e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
193e5c89e4eSSatish Balay sum of the second entry.
194b693b147SBarry Smith 
19576f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
196367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
197b693b147SBarry Smith there would be no place to store the both needed results.
198e5c89e4eSSatish Balay */
199d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt sizes[], PetscInt *max, PetscInt *sum)
200d71ae5a4SJacob Faibussowitsch {
201e5c89e4eSSatish Balay   PetscFunctionBegin;
202d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
203d6e4c47cSJed Brown   {
2049371c9d4SSatish Balay     struct {
2059371c9d4SSatish Balay       PetscInt max, sum;
2069371c9d4SSatish Balay     } work;
2079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce_scatter_block((void *)sizes, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm));
208d6e4c47cSJed Brown     *max = work.max;
209d6e4c47cSJed Brown     *sum = work.sum;
210d6e4c47cSJed Brown   }
211d6e4c47cSJed Brown #else
212d6e4c47cSJed Brown   {
213d6e4c47cSJed Brown     PetscMPIInt size, rank;
2149371c9d4SSatish Balay     struct {
2159371c9d4SSatish Balay       PetscInt max, sum;
2169371c9d4SSatish Balay     } *work;
2179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
2189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
2199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &work));
2201c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce((void *)sizes, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm));
2216ac3741eSJed Brown     *max = work[rank].max;
2226ac3741eSJed Brown     *sum = work[rank].sum;
2239566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
224d6e4c47cSJed Brown   }
225d6e4c47cSJed Brown #endif
226e5c89e4eSSatish Balay   PetscFunctionReturn(0);
227e5c89e4eSSatish Balay }
228e5c89e4eSSatish Balay 
229e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
230e5c89e4eSSatish Balay 
2319e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
232613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
233613bf2b2SPierre Jolivet     #include <quadmath.h>
234613bf2b2SPierre Jolivet   #endif
2359e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0;
236de272c7aSSatish Balay   #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
23706a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
238613bf2b2SPierre Jolivet   #endif
239e5c89e4eSSatish Balay 
240d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
241d71ae5a4SJacob Faibussowitsch {
242e5c89e4eSSatish Balay   PetscInt i, count = *cnt;
243e5c89e4eSSatish Balay 
244e5c89e4eSSatish Balay   PetscFunctionBegin;
2457c2de775SJed Brown   if (*datatype == MPIU_REAL) {
246e2e03761SBarry Smith     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
247a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2487c2de775SJed Brown   }
2497c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
250c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2517c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
252a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2537c2de775SJed Brown   }
2547c2de775SJed Brown   #endif
255613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
256613bf2b2SPierre Jolivet   else if (*datatype == MPIU___FLOAT128) {
257613bf2b2SPierre Jolivet     __float128 *xin = (__float128 *)in, *xout = (__float128 *)out;
258613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2599371c9d4SSatish Balay   } else if (*datatype == MPIU___COMPLEX128) {
260613bf2b2SPierre Jolivet     __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out;
261613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
262613bf2b2SPierre Jolivet   }
263613bf2b2SPierre Jolivet   #endif
2649e517322SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FP16)
2659e517322SPierre Jolivet   else if (*datatype == MPIU___FP16) {
2669e517322SPierre Jolivet     __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out;
2679e517322SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2689e517322SPierre Jolivet   }
2699e517322SPierre Jolivet   #endif
2707c2de775SJed Brown   else {
2719e517322SPierre Jolivet   #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16)
2727c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
2739e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FP16)
274613bf2b2SPierre Jolivet     (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types");
2759e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FLOAT128)
2769e517322SPierre Jolivet     (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types");
2779e517322SPierre Jolivet   #else
2789e517322SPierre Jolivet     (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types");
279613bf2b2SPierre Jolivet   #endif
28041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
281e2e03761SBarry Smith   }
282812af9f3SBarry Smith   PetscFunctionReturnVoid();
283e5c89e4eSSatish Balay }
284e5c89e4eSSatish Balay #endif
285e5c89e4eSSatish Balay 
286570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
287d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
288d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
289d9822059SBarry Smith 
290d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
291d71ae5a4SJacob Faibussowitsch {
292d9822059SBarry Smith   PetscInt i, count = *cnt;
293d9822059SBarry Smith 
294d9822059SBarry Smith   PetscFunctionBegin;
2957c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2968c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
297a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]);
2987c2de775SJed Brown   }
2997c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3007c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3017c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
302ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3037c2de775SJed Brown   }
3047c2de775SJed Brown   #endif
3057c2de775SJed Brown   else {
3067c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
30741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3088c764dc5SJose Roman   }
309d9822059SBarry Smith   PetscFunctionReturnVoid();
310d9822059SBarry Smith }
311d9822059SBarry Smith 
312d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
313d71ae5a4SJacob Faibussowitsch {
314d9822059SBarry Smith   PetscInt i, count = *cnt;
315d9822059SBarry Smith 
316d9822059SBarry Smith   PetscFunctionBegin;
3177c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3188c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
319a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]);
3207c2de775SJed Brown   }
3217c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3227c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3237c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
324ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3257c2de775SJed Brown   }
3267c2de775SJed Brown   #endif
3277c2de775SJed Brown   else {
3288c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
32941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3308c764dc5SJose Roman   }
331d9822059SBarry Smith   PetscFunctionReturnVoid();
332d9822059SBarry Smith }
333d9822059SBarry Smith #endif
334d9822059SBarry Smith 
335480cf27aSJed Brown /*
336480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
337480cf27aSJed Brown 
338ff0e51ddSBarry 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.
339480cf27aSJed Brown 
34012801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
341480cf27aSJed Brown 
342480cf27aSJed Brown */
343d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state)
344d71ae5a4SJacob Faibussowitsch {
34505342407SJunchao Zhang   PetscCommCounter      *counter = (PetscCommCounter *)count_val;
34657f21012SBarry Smith   struct PetscCommStash *comms   = counter->comms, *pcomm;
347480cf27aSJed Brown 
348480cf27aSJed Brown   PetscFunctionBegin;
3499566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm));
3509566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter->iflags));
35157f21012SBarry Smith   while (comms) {
3529566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&comms->comm));
35357f21012SBarry Smith     pcomm = comms;
35457f21012SBarry Smith     comms = comms->next;
3559566063dSJacob Faibussowitsch     PetscCall(PetscFree(pcomm));
35657f21012SBarry Smith   }
3579566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter));
358480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
359480cf27aSJed Brown }
360480cf27aSJed Brown 
361480cf27aSJed Brown /*
36247435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
363da3039f7SJed Brown   calls MPI_Comm_free().
364da3039f7SJed Brown 
365da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
366480cf27aSJed Brown 
367ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
368480cf27aSJed Brown 
36912801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
370480cf27aSJed Brown 
371480cf27aSJed Brown */
372d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
373d71ae5a4SJacob Faibussowitsch {
3749371c9d4SSatish Balay   union
375480cf27aSJed Brown   {
3769371c9d4SSatish Balay     MPI_Comm comm;
3779371c9d4SSatish Balay     void    *ptr;
3789371c9d4SSatish Balay   } icomm;
379480cf27aSJed Brown 
380480cf27aSJed Brown   PetscFunctionBegin;
38112801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval");
382ec4fadc2SJed Brown   icomm.ptr = attr_val;
38376bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
38433779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
38533779a13SJunchao Zhang     PetscMPIInt flg;
3869371c9d4SSatish Balay     union
3879371c9d4SSatish Balay     {
3889371c9d4SSatish Balay       MPI_Comm comm;
3899371c9d4SSatish Balay       void    *ptr;
3909371c9d4SSatish Balay     } ocomm;
3919566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg));
39233779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute");
39333779a13SJunchao 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");
39433779a13SJunchao Zhang   }
3959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval));
3969566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm));
397da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
398b89831e5SBarry Smith }
399da3039f7SJed Brown 
400da3039f7SJed Brown /*
40133779a13SJunchao 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.
402da3039f7SJed Brown  */
403d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
404d71ae5a4SJacob Faibussowitsch {
405da3039f7SJed Brown   PetscFunctionBegin;
4069566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm));
407480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
408480cf27aSJed Brown }
409480cf27aSJed Brown 
41033779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *);
41142218b76SBarry Smith 
412951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
4138cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *);
4148cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
4158cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
416e39fd77fSBarry Smith #endif
417e39fd77fSBarry Smith 
4180f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE;
41912801b39SBarry Smith 
420eb27c7c8SSatish Balay PETSC_INTERN int    PetscGlobalArgc;
421eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4226ae9a8a6SBarry Smith int                 PetscGlobalArgc = 0;
42302c9f0b5SLisandro Dalcin char              **PetscGlobalArgs = NULL;
424dff31646SBarry Smith PetscSegBuffer      PetscCitationsList;
425e5c89e4eSSatish Balay 
426d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void)
427d71ae5a4SJacob Faibussowitsch {
428051e4cf2SJed Brown   PetscFunctionBegin;
4299566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList));
43030c35bf2SSatish Balay 
43130c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\
43230c35bf2SSatish Balay   Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\
43330c35bf2SSatish Balay     and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\
4343f81df79SBarry Smith     and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\
43530c35bf2SSatish Balay     and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\
43630c35bf2SSatish Balay     and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\
43730c35bf2SSatish Balay     and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\
43830c35bf2SSatish Balay     and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\
43930c35bf2SSatish Balay   Title = {{PETSc/TAO} Users Manual},\n\
44030c35bf2SSatish Balay   Number = {ANL-21/39 - Revision 3.17},\n\
44130c35bf2SSatish Balay   Institution = {Argonne National Laboratory},\n\
4429371c9d4SSatish Balay   Year = {2022}\n}\n",
4439371c9d4SSatish Balay                                    NULL));
44430c35bf2SSatish Balay 
44530c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\
44630c35bf2SSatish Balay   Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\
44730c35bf2SSatish Balay   Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\
44830c35bf2SSatish Balay   Booktitle = {Modern Software Tools in Scientific Computing},\n\
44930c35bf2SSatish Balay   Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\
45030c35bf2SSatish Balay   Pages = {163--202},\n\
45130c35bf2SSatish Balay   Publisher = {Birkh{\\\"{a}}user Press},\n\
4529371c9d4SSatish Balay   Year = {1997}\n}\n",
4539371c9d4SSatish Balay                                    NULL));
45430c35bf2SSatish Balay 
455051e4cf2SJed Brown   PetscFunctionReturn(0);
456051e4cf2SJed Brown }
457e5c89e4eSSatish Balay 
4582d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4592d747510SLisandro Dalcin 
460d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[])
461d71ae5a4SJacob Faibussowitsch {
4622d747510SLisandro Dalcin   PetscFunctionBegin;
4639566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(programname, name, sizeof(programname)));
4642d747510SLisandro Dalcin   PetscFunctionReturn(0);
4652d747510SLisandro Dalcin }
4662d747510SLisandro Dalcin 
4672d747510SLisandro Dalcin /*@C
4682d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4692d747510SLisandro Dalcin 
4702d747510SLisandro Dalcin     Not Collective
4712d747510SLisandro Dalcin 
4722d747510SLisandro Dalcin     Input Parameter:
4732d747510SLisandro Dalcin .   len - length of the string name
4742d747510SLisandro Dalcin 
4752d747510SLisandro Dalcin     Output Parameter:
476811af0c4SBarry Smith .   name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN`
4772d747510SLisandro Dalcin 
4782d747510SLisandro Dalcin    Level: advanced
4792d747510SLisandro Dalcin 
4802d747510SLisandro Dalcin @*/
481d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProgramName(char name[], size_t len)
482d71ae5a4SJacob Faibussowitsch {
4832d747510SLisandro Dalcin   PetscFunctionBegin;
4849566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(name, programname, len));
4852d747510SLisandro Dalcin   PetscFunctionReturn(0);
4862d747510SLisandro Dalcin }
4872d747510SLisandro Dalcin 
488e5c89e4eSSatish Balay /*@C
489e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
490811af0c4SBarry Smith      after PetscInitialize() is called but before `PetscFinalize()`.
491e5c89e4eSSatish Balay 
492e5c89e4eSSatish Balay    Not Collective
493e5c89e4eSSatish Balay 
494e5c89e4eSSatish Balay    Output Parameters:
495e5c89e4eSSatish Balay +  argc - count of number of command line arguments
496e5c89e4eSSatish Balay -  args - the command line arguments
497e5c89e4eSSatish Balay 
498e5c89e4eSSatish Balay    Level: intermediate
499e5c89e4eSSatish Balay 
500e5c89e4eSSatish Balay    Notes:
501e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
502e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
503e5c89e4eSSatish Balay 
504f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
505f177e3b1SBarry Smith 
506db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`
507e5c89e4eSSatish Balay @*/
508d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args)
509d71ae5a4SJacob Faibussowitsch {
510e5c89e4eSSatish Balay   PetscFunctionBegin;
51139a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
512e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
513e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
514e5c89e4eSSatish Balay   PetscFunctionReturn(0);
515e5c89e4eSSatish Balay }
516e5c89e4eSSatish Balay 
517793721a6SBarry Smith /*@C
518793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
519811af0c4SBarry Smith      after `PetscInitialize()` is called but before `PetscFinalize()`.
520793721a6SBarry Smith 
521793721a6SBarry Smith    Not Collective
522793721a6SBarry Smith 
523793721a6SBarry Smith    Output Parameters:
524793721a6SBarry Smith .  args - the command line arguments
525793721a6SBarry Smith 
526793721a6SBarry Smith    Level: intermediate
527793721a6SBarry Smith 
528793721a6SBarry Smith    Notes:
529793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
530793721a6SBarry Smith 
531db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`
532793721a6SBarry Smith @*/
533d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args)
534d71ae5a4SJacob Faibussowitsch {
535793721a6SBarry Smith   PetscInt i, argc = PetscGlobalArgc;
536793721a6SBarry Smith 
537793721a6SBarry Smith   PetscFunctionBegin;
53839a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
5399371c9d4SSatish Balay   if (!argc) {
5409371c9d4SSatish Balay     *args = NULL;
5419371c9d4SSatish Balay     PetscFunctionReturn(0);
5429371c9d4SSatish Balay   }
5439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(argc, args));
5449566063dSJacob Faibussowitsch   for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i]));
5452d747510SLisandro Dalcin   (*args)[argc - 1] = NULL;
546793721a6SBarry Smith   PetscFunctionReturn(0);
547793721a6SBarry Smith }
548793721a6SBarry Smith 
549793721a6SBarry Smith /*@C
550811af0c4SBarry Smith    PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()`
551793721a6SBarry Smith 
552793721a6SBarry Smith    Not Collective
553793721a6SBarry Smith 
554793721a6SBarry Smith    Output Parameters:
555793721a6SBarry Smith .  args - the command line arguments
556793721a6SBarry Smith 
557793721a6SBarry Smith    Level: intermediate
558793721a6SBarry Smith 
559db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()`
560793721a6SBarry Smith @*/
561d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args)
562d71ae5a4SJacob Faibussowitsch {
56339a651e2SJacob Faibussowitsch   PetscFunctionBegin;
56439a651e2SJacob Faibussowitsch   if (args) {
565793721a6SBarry Smith     PetscInt i = 0;
566793721a6SBarry Smith 
5679566063dSJacob Faibussowitsch     while (args[i]) PetscCall(PetscFree(args[i++]));
5689566063dSJacob Faibussowitsch     PetscCall(PetscFree(args));
56939a651e2SJacob Faibussowitsch   }
570793721a6SBarry Smith   PetscFunctionReturn(0);
571793721a6SBarry Smith }
572793721a6SBarry Smith 
57327104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
57430befbd2SBarry Smith   #include <petscconfiginfo.h>
57530befbd2SBarry Smith 
576d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
577d71ae5a4SJacob Faibussowitsch {
57827104ee2SJacob Faibussowitsch   PetscFunctionBegin;
57911525c0dSBarry Smith   if (!PetscGlobalRank) {
58030befbd2SBarry Smith     char      cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64];
58111525c0dSBarry Smith     int       port;
582ffbd1cfbSBarry Smith     PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE;
58311525c0dSBarry Smith     size_t    applinelen, introlen;
584ffbd1cfbSBarry Smith     char      sawsurl[256];
58511525c0dSBarry Smith 
5869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg));
58711525c0dSBarry Smith     if (flg) {
58811525c0dSBarry Smith       char sawslog[PETSC_MAX_PATH_LEN];
58911525c0dSBarry Smith 
5909566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL));
59111525c0dSBarry Smith       if (sawslog[0]) {
592792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog));
59311525c0dSBarry Smith       } else {
594792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL));
59511525c0dSBarry Smith       }
59611525c0dSBarry Smith     }
5979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg));
59848a46eb9SPierre Jolivet     if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert));
5999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL));
600ffbd1cfbSBarry Smith     if (selectport) {
601792fecdfSBarry Smith       PetscCallSAWs(SAWs_Get_Available_Port, (&port));
602792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Port, (port));
603ffbd1cfbSBarry Smith     } else {
6049566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg));
60548a46eb9SPierre Jolivet       if (flg) PetscCallSAWs(SAWs_Set_Port, (port));
606ffbd1cfbSBarry Smith     }
6079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg));
60811525c0dSBarry Smith     if (flg) {
609792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Document_Root, (root));
6109566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(root, ".", &rootlocal));
6119c1e0ce8SBarry Smith     } else {
6129566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg));
6139c1e0ce8SBarry Smith       if (flg) {
6149566063dSJacob Faibussowitsch         PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root)));
615792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Document_Root, (root));
6169c1e0ce8SBarry Smith       }
61711525c0dSBarry Smith     }
6189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2));
61911525c0dSBarry Smith     if (flg2) {
62011525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
62128b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option");
6229566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root));
6239566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(jsdir, 'r', &flg));
62428b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory");
625792fecdfSBarry Smith       PetscCallSAWs(SAWs_Push_Local_Header, ());
62611525c0dSBarry Smith     }
6279566063dSJacob Faibussowitsch     PetscCall(PetscGetProgramName(programname, sizeof(programname)));
6289566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(help, &applinelen));
62911525c0dSBarry Smith     introlen = 4096 + applinelen;
63030a8c9c0SSurtai Han     applinelen += 1024;
6319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(applinelen, &appline));
6329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(introlen, &intro));
63311525c0dSBarry Smith 
63411525c0dSBarry Smith     if (rootlocal) {
6359566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname));
6369566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(appline, 'r', &rootlocal));
63711525c0dSBarry Smith     }
6389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetAll(NULL, &options));
63911525c0dSBarry Smith     if (rootlocal && help) {
6409566063dSJacob 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));
64111525c0dSBarry Smith     } else if (help) {
6429566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help));
64311525c0dSBarry Smith     } else {
6449566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options));
64511525c0dSBarry Smith     }
6469566063dSJacob Faibussowitsch     PetscCall(PetscFree(options));
6479566063dSJacob Faibussowitsch     PetscCall(PetscGetVersion(version, sizeof(version)));
6489371c9d4SSatish Balay     PetscCall(PetscSNPrintf(intro, introlen,
6499371c9d4SSatish Balay                             "<body>\n"
650a17b96a8SKyle 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"
651df62c222SBarry 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"
6529371c9d4SSatish Balay                             "%s",
6539371c9d4SSatish Balay                             version, petscconfigureoptions, appline));
654792fecdfSBarry Smith     PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro));
6559566063dSJacob Faibussowitsch     PetscCall(PetscFree(intro));
6569566063dSJacob Faibussowitsch     PetscCall(PetscFree(appline));
657ffbd1cfbSBarry Smith     if (selectport) {
658aa573868SBarry Smith       PetscBool silent;
6597d812c46SBarry Smith 
6607d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
66139a651e2SJacob Faibussowitsch       while (SAWs_Initialize()) {
662792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_Available_Port, (&port));
663792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port, (port));
6647d812c46SBarry Smith       }
6657d812c46SBarry Smith 
6669566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL));
667aa573868SBarry Smith       if (!silent) {
668792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl));
6699566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl));
670ffbd1cfbSBarry Smith       }
6717d812c46SBarry Smith     } else {
672792fecdfSBarry Smith       PetscCallSAWs(SAWs_Initialize, ());
673aa573868SBarry Smith     }
6749566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister("@TechReport{ saws,\n"
6750af79b04SBarry Smith                                      "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6760af79b04SBarry Smith                                      "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6770af79b04SBarry Smith                                      "  Institution = {Argonne National Laboratory},\n"
6789371c9d4SSatish Balay                                      "  Year   = 2013\n}\n",
6799371c9d4SSatish Balay                                      NULL));
68011525c0dSBarry Smith   }
68111525c0dSBarry Smith   PetscFunctionReturn(0);
68211525c0dSBarry Smith }
68311525c0dSBarry Smith #endif
68411525c0dSBarry Smith 
6854dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
686d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
687d71ae5a4SJacob Faibussowitsch {
6884dfee713SSatish Balay   PetscFunctionBegin;
6894dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6904dfee713SSatish Balay   /* see MPI.py for details on this bug */
6914dfee713SSatish Balay   (void)setenv("HWLOC_COMPONENTS", "-x86", 1);
6924dfee713SSatish Balay #endif
6934dfee713SSatish Balay   PetscFunctionReturn(0);
6944dfee713SSatish Balay }
6954dfee713SSatish Balay 
696a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
697a56f64adSBarry Smith   #include <adios.h>
69822580e64SBarry Smith   #include <adios_read.h>
6997b56e58cSSatish Balay int64_t Petsc_adios_group;
700a56f64adSBarry Smith #endif
701a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
702cd1b99a6SBarry Smith   #include <omp.h>
703f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
704cd1b99a6SBarry Smith #endif
705a56f64adSBarry Smith 
706a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
707a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA)
7080e6b6b59SJacob Faibussowitsch   #include <petscdevice_cuda.h>
709a4af0ceeSJacob Faibussowitsch // REMOVE ME
710a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
711a4af0ceeSJacob Faibussowitsch #endif
712a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP)
7130e6b6b59SJacob Faibussowitsch   #include <petscdevice_hip.h>
714a4af0ceeSJacob Faibussowitsch // REMOVE ME
715a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
716a4af0ceeSJacob Faibussowitsch #endif
717a4af0ceeSJacob Faibussowitsch 
71827104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
719bc8a24ecSBarry Smith   #include <dlfcn.h>
720bc8a24ecSBarry Smith #endif
721a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
722a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
723a4af0ceeSJacob Faibussowitsch #endif
724a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
725a4af0ceeSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViennaCLInit();
726a4af0ceeSJacob Faibussowitsch PetscBool                   PetscViennaCLSynchronize = PETSC_FALSE;
727a4af0ceeSJacob Faibussowitsch #endif
728bc8a24ecSBarry Smith 
729660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE;
730660278c0SBarry Smith 
73185649d77SJunchao Zhang /*
73285649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
73385649d77SJunchao Zhang 
73485649d77SJunchao Zhang   prog:     program name
73502101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
73685649d77SJunchao Zhang   help:     program help message
73702101c96SSatish Balay   ftn:      is it called from Fortran initilization (petscinitializef_)?
73885649d77SJunchao Zhang   readarguments,len: used when fortran is true
73985649d77SJunchao Zhang */
740d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len)
741d71ae5a4SJacob Faibussowitsch {
74285649d77SJunchao Zhang   PetscMPIInt size;
74385649d77SJunchao Zhang   PetscBool   flg = PETSC_TRUE;
74485649d77SJunchao Zhang   char        hostname[256];
74585649d77SJunchao Zhang 
74627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
74727104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
74839a651e2SJacob Faibussowitsch   /* these must be initialized in a routine, not as a constant declaration */
74939a651e2SJacob Faibussowitsch   PETSC_STDOUT = stdout;
75039a651e2SJacob Faibussowitsch   PETSC_STDERR = stderr;
75139a651e2SJacob Faibussowitsch 
7529566063dSJacob Faibussowitsch   /* PetscCall can be used from now */
75339a651e2SJacob Faibussowitsch   PetscErrorHandlingInitialized = PETSC_TRUE;
75439a651e2SJacob Faibussowitsch 
75585649d77SJunchao Zhang   /*
75685649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
75785649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
75885649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
75985649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
76085649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
76185649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
76285649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
76385649d77SJunchao Zhang       listed above and since that time are compatible.
76485649d77SJunchao Zhang 
76585649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
76685649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
76785649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
76885649d77SJunchao Zhang       to perform the checking.
76985649d77SJunchao Zhang 
77085649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
77185649d77SJunchao Zhang 
77285649d77SJunchao Zhang       Questions:
77385649d77SJunchao Zhang 
77485649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
77585649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
77685649d77SJunchao Zhang   */
77785649d77SJunchao Zhang 
77885649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
77985649d77SJunchao Zhang   {
78085649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
78185649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
78239a651e2SJacob Faibussowitsch 
7839566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength));
78485649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
78585649d77SJunchao Zhang   #if defined(MPICH_VERSION)
78685649d77SJunchao Zhang     #if MPICH_NUMVERSION < 30100000
78785649d77SJunchao Zhang     {
78885649d77SJunchao Zhang       char     *ver, *lf;
78985649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
79039a651e2SJacob Faibussowitsch 
7919566063dSJacob Faibussowitsch       PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver));
79239a651e2SJacob Faibussowitsch       if (ver) {
7939566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(ver, '\n', &lf));
79439a651e2SJacob Faibussowitsch         if (lf) {
79585649d77SJunchao Zhang           *lf = 0;
7969566063dSJacob Faibussowitsch           PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg));
79785649d77SJunchao Zhang         }
79885649d77SJunchao Zhang       }
79985649d77SJunchao Zhang       if (!flg) {
800d5b396e8SSatish Balay         PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION));
80185649d77SJunchao Zhang         flg = PETSC_TRUE;
80285649d77SJunchao Zhang       }
80385649d77SJunchao Zhang     }
80485649d77SJunchao Zhang     #endif
80585649d77SJunchao 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?) */
80685649d77SJunchao Zhang   #elif defined(OMPI_MAJOR_VERSION)
80785649d77SJunchao Zhang     {
80885649d77SJunchao Zhang       char     *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf;
80985649d77SJunchao Zhang       PetscBool flg                                              = PETSC_FALSE;
81085649d77SJunchao Zhang     #define PSTRSZ 2
81185649d77SJunchao Zhang       char      ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"};
81285649d77SJunchao Zhang       char      ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "};
81385649d77SJunchao Zhang       int       i;
81485649d77SJunchao Zhang       for (i = 0; i < PSTRSZ; i++) {
8159566063dSJacob Faibussowitsch         PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver));
81639a651e2SJacob Faibussowitsch         if (ver) {
8179566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION));
8189566063dSJacob Faibussowitsch           PetscCall(PetscStrstr(ver, bs, &bsf));
81939a651e2SJacob Faibussowitsch           if (bsf) flg = PETSC_TRUE;
82085649d77SJunchao Zhang           break;
82185649d77SJunchao Zhang         }
82285649d77SJunchao Zhang       }
82385649d77SJunchao Zhang       if (!flg) {
8247d3de750SJacob Faibussowitsch         PetscInfo(NULL, "PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n", mpilibraryversion, OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION);
82585649d77SJunchao Zhang         flg = PETSC_TRUE;
82685649d77SJunchao Zhang       }
82785649d77SJunchao Zhang     }
82885649d77SJunchao Zhang   #endif
82985649d77SJunchao Zhang   }
83085649d77SJunchao Zhang #endif
83185649d77SJunchao Zhang 
832*e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__))
83385649d77SJunchao 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 */
83439a651e2SJacob 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");
83585649d77SJunchao Zhang #endif
83685649d77SJunchao Zhang 
83785649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
83885649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
83985649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
84085649d77SJunchao Zhang #endif
84185649d77SJunchao Zhang 
8429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
84385649d77SJunchao Zhang 
84485649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
84585649d77SJunchao Zhang 
8469566063dSJacob Faibussowitsch   PetscCall(PetscSetProgramName(prog));
8479566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen));
8489566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout));
8499566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr));
8509566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscCommSpinLock));
85185649d77SJunchao Zhang 
85285649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
8539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN));
85485649d77SJunchao Zhang 
85585649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
8569566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS));
8579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE));
85885649d77SJunchao Zhang   }
85985649d77SJunchao Zhang 
86085649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
8619566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
86285649d77SJunchao Zhang 
8639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank));
8649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize));
86585649d77SJunchao Zhang 
86685649d77SJunchao Zhang   MPIU_BOOL        = MPI_INT;
86785649d77SJunchao Zhang   MPIU_ENUM        = MPI_INT;
86885649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64;
86985649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
87085649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
87185649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
87285649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
87385649d77SJunchao Zhang #endif
87439a651e2SJacob Faibussowitsch   else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t");
87585649d77SJunchao Zhang 
87685649d77SJunchao Zhang     /*
87785649d77SJunchao Zhang      Initialized the global complex variable; this is because with
87885649d77SJunchao Zhang      shared libraries the constructors for global variables
87985649d77SJunchao Zhang      are not called; at least on IRIX.
88085649d77SJunchao Zhang   */
88185649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
88285649d77SJunchao Zhang   {
88385649d77SJunchao Zhang   #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
88485649d77SJunchao Zhang     PetscComplex ic(0.0, 1.0);
88585649d77SJunchao Zhang     PETSC_i = ic;
88685649d77SJunchao Zhang   #else
88785649d77SJunchao Zhang     PETSC_i = _Complex_I;
88885649d77SJunchao Zhang   #endif
88985649d77SJunchao Zhang   }
89085649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
89185649d77SJunchao Zhang 
89285649d77SJunchao Zhang   /*
89385649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
89485649d77SJunchao Zhang      half of the entries and maxes the second half.
89585649d77SJunchao Zhang   */
8969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP));
89785649d77SJunchao Zhang 
898613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
8999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128));
9009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128));
9019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128));
9029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128));
90385649d77SJunchao Zhang #endif
9049e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
9059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16));
9069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FP16));
90785649d77SJunchao Zhang #endif
90885649d77SJunchao Zhang 
90985649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
9109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM));
911613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX));
912613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN));
9139e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
9149e517322SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128));
91585649d77SJunchao Zhang #endif
91685649d77SJunchao Zhang 
9179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR));
9189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR));
91985649d77SJunchao Zhang 
92085649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
92185649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
92285649d77SJunchao Zhang   {
92385649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
92493d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)};
92585649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_REAL, MPIU_INT}, tmpStruct;
92685649d77SJunchao Zhang 
9279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
92893d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT));
9299566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9309566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT));
93185649d77SJunchao Zhang   }
93285649d77SJunchao Zhang   {
93385649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
93493d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)};
93585649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_SCALAR, MPIU_INT}, tmpStruct;
93685649d77SJunchao Zhang 
9379566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
93893d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT));
9399566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
94185649d77SJunchao Zhang   }
94285649d77SJunchao Zhang #endif
94385649d77SJunchao Zhang 
94485649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT));
9469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2INT));
94785649d77SJunchao Zhang #endif
9489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT));
9499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPI_4INT));
9509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT));
9519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_4INT));
95285649d77SJunchao Zhang 
95385649d77SJunchao Zhang   /*
95485649d77SJunchao Zhang      Attributes to be set on PETSc communicators
95585649d77SJunchao Zhang   */
9569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0));
9579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0));
9589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0));
9599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0));
96085649d77SJunchao Zhang 
96185649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
9629566063dSJacob Faibussowitsch   if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len));
96385649d77SJunchao Zhang   else
96485649d77SJunchao Zhang #endif
9659566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file));
96685649d77SJunchao Zhang 
96785649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9689566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
96985649d77SJunchao Zhang 
97085649d77SJunchao Zhang   /*
97185649d77SJunchao Zhang      Check system options and print help
97285649d77SJunchao Zhang   */
9739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCheckInitial_Private(help));
97485649d77SJunchao Zhang 
975a4af0ceeSJacob Faibussowitsch   /*
9760e6b6b59SJacob Faibussowitsch     Creates the logging data structures; this is enabled even if logging is not turned on
9770e6b6b59SJacob Faibussowitsch     This is the last thing we do before returning to the user code to prevent having the
9780e6b6b59SJacob Faibussowitsch     logging numbers contaminated by any startup time associated with MPI
9790e6b6b59SJacob Faibussowitsch   */
9800e6b6b59SJacob Faibussowitsch #if defined(PETSC_USE_LOG)
9810e6b6b59SJacob Faibussowitsch   PetscCall(PetscLogInitialize());
9820e6b6b59SJacob Faibussowitsch #endif
9830e6b6b59SJacob Faibussowitsch 
9840e6b6b59SJacob Faibussowitsch   /*
985a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
986a4af0ceeSJacob Faibussowitsch 
987a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
988a4af0ceeSJacob Faibussowitsch    1. MPI initialized
989a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
9900e6b6b59SJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up
9910e6b6b59SJacob Faibussowitsch       its own SIGSEV handler via the push/pop interface.
9920e6b6b59SJacob Faibussowitsch    4. Logging initialized
993a4af0ceeSJacob Faibussowitsch   */
9949566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
995a4af0ceeSJacob Faibussowitsch 
996a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
997a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
9989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_summary", &flg));
9999566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg));
10009566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL));
1001a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
10029566063dSJacob Faibussowitsch   PetscCall(PetscViennaCLInit());
1003a4af0ceeSJacob Faibussowitsch #endif
1004a4af0ceeSJacob Faibussowitsch 
10059566063dSJacob Faibussowitsch   PetscCall(PetscCitationsInitialize());
100685649d77SJunchao Zhang 
100785649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
10089566063dSJacob Faibussowitsch   PetscCall(PetscInitializeSAWs(ftn ? NULL : help));
100927104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
10109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg));
10119566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscStackViewSAWs());
101285649d77SJunchao Zhang #endif
101385649d77SJunchao Zhang 
101485649d77SJunchao Zhang   /*
101585649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
101685649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
101785649d77SJunchao Zhang   */
10189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_DynamicLibraries());
101985649d77SJunchao Zhang 
10209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10219566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size));
10229566063dSJacob Faibussowitsch   PetscCall(PetscGetHostName(hostname, 256));
10239566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname));
102485649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
102585649d77SJunchao Zhang   {
102685649d77SJunchao Zhang     PetscBool omp_view_flag;
102785649d77SJunchao Zhang     char     *threads = getenv("OMP_NUM_THREADS");
102885649d77SJunchao Zhang 
102985649d77SJunchao Zhang     if (threads) {
10309566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads));
103185649d77SJunchao Zhang       (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads);
103285649d77SJunchao Zhang     } else {
10332f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt)omp_get_max_threads();
10349566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads));
103585649d77SJunchao Zhang     }
1036d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys");
10379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg));
10389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag));
1039d0609cedSBarry Smith     PetscOptionsEnd();
104085649d77SJunchao Zhang     if (flg) {
10419566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads));
104285649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
104385649d77SJunchao Zhang     }
104448a46eb9SPierre Jolivet     if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads));
104585649d77SJunchao Zhang   }
104685649d77SJunchao Zhang #endif
104785649d77SJunchao Zhang 
104885649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
104985649d77SJunchao Zhang   /*
105085649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
105185649d77SJunchao Zhang 
105285649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
105385649d77SJunchao Zhang   */
10549566063dSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL));
105585649d77SJunchao Zhang #endif
105685649d77SJunchao Zhang 
105785649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10589566063dSJacob Faibussowitsch   PetscCall(PetscFPTCreate(10000));
105985649d77SJunchao Zhang #endif
106085649d77SJunchao Zhang 
106185649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
106285649d77SJunchao Zhang   {
106385649d77SJunchao Zhang     PetscViewer viewer;
10649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg));
106585649d77SJunchao Zhang     if (flg) {
10669566063dSJacob Faibussowitsch       PetscCall(PetscProcessPlacementView(viewer));
10679566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
106885649d77SJunchao Zhang     }
106985649d77SJunchao Zhang   }
107085649d77SJunchao Zhang #endif
107185649d77SJunchao Zhang 
107285649d77SJunchao Zhang   flg = PETSC_TRUE;
10739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL));
10749566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE));
107585649d77SJunchao Zhang 
107685649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10779566063dSJacob Faibussowitsch   PetscCall(adios_init_noxml(PETSC_COMM_WORLD));
10789566063dSJacob Faibussowitsch   PetscCall(adios_declare_group(&Petsc_adios_group, "PETSc", "", adios_stat_default));
10799566063dSJacob Faibussowitsch   PetscCall(adios_select_method(Petsc_adios_group, "MPI", "", ""));
10809566063dSJacob Faibussowitsch   PetscCall(adios_read_init_method(ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, ""));
108185649d77SJunchao Zhang #endif
108285649d77SJunchao Zhang 
108385649d77SJunchao Zhang #if defined(__VALGRIND_H)
108485649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE;
108585649d77SJunchao Zhang   #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
10869566063dSJacob 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"));
108785649d77SJunchao Zhang   #endif
108885649d77SJunchao Zhang #endif
108985649d77SJunchao Zhang   /*
109085649d77SJunchao Zhang       Set flag that we are completely initialized
109185649d77SJunchao Zhang   */
109285649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
109385649d77SJunchao Zhang 
10949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg));
10959566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPythonInitialize(NULL, NULL));
1096f1f2ae84SBarry Smith 
1097f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1098f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin());
1099f1f2ae84SBarry 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");
110027104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
110185649d77SJunchao Zhang }
110285649d77SJunchao Zhang 
1103e5c89e4eSSatish Balay /*@C
1104e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1105811af0c4SBarry Smith    `PetscInitialize()` calls MPI_Init() if that has yet to be called,
1106e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1107e5c89e4eSSatish Balay    your program -- usually the very first line!
1108e5c89e4eSSatish Balay 
1109811af0c4SBarry Smith    Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set
1110e5c89e4eSSatish Balay 
1111e5c89e4eSSatish Balay    Input Parameters:
1112e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1113e5c89e4eSSatish Balay .  args - the command line arguments
1114be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1115be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1116be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1117c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
11180298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1119e5c89e4eSSatish Balay 
1120811af0c4SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that
1121811af0c4SBarry Smith    communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a
1122811af0c4SBarry Smith    four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not,
1123811af0c4SBarry Smith    then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even
112405827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1125e5c89e4eSSatish Balay 
1126e5c89e4eSSatish Balay    Options Database Keys:
11277ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
11287ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1129e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11308a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
1131811af0c4SBarry Smith .  -on_error_abort - calls `abort()` when error detected (no traceback)
1132811af0c4SBarry Smith .  -on_error_mpiabort - calls `MPI_abort()` when error detected
11331af3601dSBarry Smith .  -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr
11348a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1135bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1136e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1137e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1138e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
113979dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
114079dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
1141811af0c4SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see `PetscMallocSetDebug()`
1142aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
114392f119d6SBarry 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
1144811af0c4SBarry Smith .  -malloc_view - show a list of all allocated memory during `PetscFinalize()`
114592f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1146608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
114792f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1148e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1149e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1150e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1151e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1152e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
11530841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1154e5c89e4eSSatish Balay 
1155c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1156c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1157c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1158811af0c4SBarry Smith -  -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()`
1159c5b5d8d5SVaclav Hapla 
1160c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1161c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1162c5b5d8d5SVaclav Hapla    They can be used also in option files.
1163c5b5d8d5SVaclav Hapla 
1164811af0c4SBarry Smith    See `PetscOptionsMonitorSet()` to do monitoring programmatically.
1165c5b5d8d5SVaclav Hapla 
1166e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1167a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1168811af0c4SBarry Smith +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`.
1169217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1170217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1171495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1172811af0c4SBarry Smith         hangs without running in the debugger).  See `PetscLogTraceBegin()`.
1173811af0c4SBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see `PetscLogView()`.
1174811af0c4SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`.
1175811af0c4SBarry Smith .  -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView().
11769a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1177495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1178f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1179811af0c4SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See `PetscLogDump()`.
1180811af0c4SBarry Smith .  -log [filename] - Logs basic profiline information  See `PetscLogDump()`.
1181c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
1182811af0c4SBarry Smith .  -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off
1183c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1184495fc317SBarry Smith 
11856a77f485SPierre Jolivet     Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time
1186e5c89e4eSSatish Balay 
1187ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1188ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1189ffbd1cfbSBarry 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
1190ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1191ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1192ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1193ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1194ffbd1cfbSBarry Smith 
1195e5c89e4eSSatish Balay    Environmental Variables:
1196811af0c4SBarry Smith +   `PETSC_TMP` - alternative tmp directory
1197811af0c4SBarry Smith .   `PETSC_SHARED_TMP` - tmp is shared by all processes
1198811af0c4SBarry Smith .   `PETSC_NOT_SHARED_TMP` - each process has its own private tmp
1199811af0c4SBarry Smith .   `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs
1200811af0c4SBarry Smith .   `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1201811af0c4SBarry Smith .   `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer
1202811af0c4SBarry Smith -   `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to
1203e5c89e4eSSatish Balay 
1204e5c89e4eSSatish Balay    Level: beginner
1205e5c89e4eSSatish Balay 
1206811af0c4SBarry Smith    Note:
1207811af0c4SBarry Smith    If for some reason you must call `MPI_Init()` separately, call
1208811af0c4SBarry Smith    it before `PetscInitialize()`.
1209e5c89e4eSSatish Balay 
1210811af0c4SBarry Smith    Fortran Notes:
1211811af0c4SBarry Smith    In Fortran this routine can be called with
1212811af0c4SBarry Smith .vb
1213811af0c4SBarry Smith        call PetscInitialize(ierr)
1214811af0c4SBarry Smith        call PetscInitialize(file,ierr) or
1215811af0c4SBarry Smith        call PetscInitialize(file,help,ierr)
1216811af0c4SBarry Smith .ve
1217e5c89e4eSSatish Balay 
1218811af0c4SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after
1219811af0c4SBarry Smith    calling `PetscInitialize()`.
1220e5c89e4eSSatish Balay 
1221db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()`
1222e5c89e4eSSatish Balay @*/
1223d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[])
1224d71ae5a4SJacob Faibussowitsch {
122585649d77SJunchao Zhang   PetscMPIInt flag;
122621ef0414SBarry Smith   const char *prog = "Unknown Name", *mpienv;
1227e5c89e4eSSatish Balay 
122827104ee2SJacob Faibussowitsch   PetscFunctionBegin;
122927104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
12309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Initialized(&flag));
1231e5c89e4eSSatish Balay   if (!flag) {
123239a651e2SJacob 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");
12339566063dSJacob Faibussowitsch     PetscCall(PetscPreMPIInit_Private());
12345e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12355e765c61SJed Brown     {
123639a651e2SJacob Faibussowitsch       PetscMPIInt PETSC_UNUSED provided;
12379566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED, &provided));
12385e765c61SJed Brown     }
12395e765c61SJed Brown #else
12409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Init(argc, args));
12415e765c61SJed Brown #endif
124221ef0414SBarry Smith     if (PetscDefined(HAVE_MPIUNI)) {
124321ef0414SBarry Smith       mpienv = getenv("PMI_SIZE");
124421ef0414SBarry Smith       if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE");
124521ef0414SBarry Smith       if (mpienv) {
124621ef0414SBarry Smith         PetscInt isize;
124721ef0414SBarry Smith         PetscCall(PetscOptionsStringToInt(mpienv, &isize));
124821ef0414SBarry 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");
124921ef0414SBarry 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");
125021ef0414SBarry Smith       }
125121ef0414SBarry Smith     }
1252e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1253e5c89e4eSSatish Balay   }
12544dfee713SSatish Balay 
125585649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1256e5c89e4eSSatish Balay   if (argc && args) {
1257e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1258e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1259e5c89e4eSSatish Balay   }
1260811af0c4SBarry Smith   PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0));
126127104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
1262e5c89e4eSSatish Balay }
1263e5c89e4eSSatish Balay 
1264a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
126595c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1266ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsCounts;
1267ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsMaxCounts;
126895c0884eSLisandro Dalcin PETSC_INTERN PetscBool    PetscObjectsLog;
12694097062eSBarry Smith #endif
1270e5c89e4eSSatish Balay 
1271008a6e76SBarry Smith /*
1272008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1273008a6e76SBarry Smith */
1274d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void)
1275d71ae5a4SJacob Faibussowitsch {
1276008a6e76SBarry Smith   PetscFunctionBegin;
1277613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
12789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128));
12799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128));
1280008a6e76SBarry Smith #endif
12819e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
12829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FP16));
1283008a6e76SBarry Smith #endif
1284008a6e76SBarry Smith 
1285de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
12869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_SUM));
1287613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MAX));
1288613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MIN));
12899e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
12909e517322SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128));
1291008a6e76SBarry Smith #endif
1292008a6e76SBarry Smith 
12939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR));
12949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT));
12959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT));
129640df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
12979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2INT));
1298008a6e76SBarry Smith #endif
12999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPI_4INT));
13009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_4INT));
13019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
1302008a6e76SBarry Smith   PetscFunctionReturn(0);
1303008a6e76SBarry Smith }
1304008a6e76SBarry Smith 
1305a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
1306a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1307a4af0ceeSJacob Faibussowitsch #endif
1308a4af0ceeSJacob Faibussowitsch 
1309e5c89e4eSSatish Balay /*@C
1310e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1311811af0c4SBarry Smith    of the program. `MPI_Finalize()` is called only if the user had not
1312811af0c4SBarry Smith    called `MPI_Init()` before calling `PetscInitialize()`.
1313e5c89e4eSSatish Balay 
1314811af0c4SBarry Smith    Collective on `PETSC_COMM_WORLD`
1315e5c89e4eSSatish Balay 
1316e5c89e4eSSatish Balay    Options Database Keys:
1317811af0c4SBarry Smith +  -options_view - Calls `PetscOptionsView()`
1318e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
13197eb1d149SBarry 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
1320e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1321811af0c4SBarry Smith .  -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed
1322e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
132392f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1324e5c89e4eSSatish Balay 
1325e5c89e4eSSatish Balay    Level: beginner
1326e5c89e4eSSatish Balay 
1327e5c89e4eSSatish Balay    Note:
1328811af0c4SBarry Smith    See `PetscInitialize()` for other runtime options.
1329e5c89e4eSSatish Balay 
1330db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()`
1331e5c89e4eSSatish Balay @*/
1332d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void)
1333d71ae5a4SJacob Faibussowitsch {
13344bb5149bSJed Brown   PetscMPIInt rank;
1335a8d2bbe5SBarry Smith   PetscInt    nopt;
13362bf49c77SBarry Smith   PetscBool   flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE;
1337dff31646SBarry Smith   PetscBool   flg;
133810463e74SBarry Smith #if defined(PETSC_USE_LOG)
133910463e74SBarry Smith   char mname[PETSC_MAX_PATH_LEN];
134010463e74SBarry Smith #endif
1341e5c89e4eSSatish Balay 
13423cbbc5ffSBarry Smith   PetscFunctionBegin;
134339a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()");
13449566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PetscFinalize() called\n"));
1345b022a5c1SBarry Smith 
1346f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1347f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd());
1348f1f2ae84SBarry Smith 
13499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1350a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13519566063dSJacob Faibussowitsch   PetscCall(adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE));
13529566063dSJacob Faibussowitsch   PetscCall(adios_finalize(rank));
1353a56f64adSBarry Smith #endif
13549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg));
1355dff31646SBarry Smith   if (flg) {
13561f817a21SBarry Smith     char *cits, filename[PETSC_MAX_PATH_LEN];
13571f817a21SBarry Smith     FILE *fd = PETSC_STDOUT;
13581f817a21SBarry Smith 
13599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL));
136048a46eb9SPierre Jolivet     if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd));
13619566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits));
1362dff31646SBarry Smith     cits[0] = 0;
13639566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits));
13649566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n"));
13659566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
13669566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits));
13679566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
13689566063dSJacob Faibussowitsch     PetscCall(PetscFClose(PETSC_COMM_WORLD, fd));
13699566063dSJacob Faibussowitsch     PetscCall(PetscFree(cits));
1370dff31646SBarry Smith   }
13719566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&PetscCitationsList));
1372dff31646SBarry Smith 
1373c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
137404102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
137504102261SBarry Smith   {
137604102261SBarry Smith     PetscInt nmax = 2;
137704102261SBarry Smith     char   **buffs;
13789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
13799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-textbelt", buffs, &nmax, &flg1));
138004102261SBarry Smith     if (flg1) {
138128b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-textbelt requires either the phone number or number,\"message\"");
138204102261SBarry Smith       if (nmax == 1) {
13839566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(128, &buffs[1]));
13849566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
13859566063dSJacob Faibussowitsch         PetscCall(PetscStrcat(buffs[1], " has completed"));
138604102261SBarry Smith       }
13879566063dSJacob Faibussowitsch       PetscCall(PetscTextBelt(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
13889566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
13899566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
139004102261SBarry Smith     }
13919566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
139204102261SBarry Smith   }
1393763ec7b1SBarry Smith   {
1394763ec7b1SBarry Smith     PetscInt nmax = 2;
1395763ec7b1SBarry Smith     char   **buffs;
13969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
13979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-tellmycell", buffs, &nmax, &flg1));
1398763ec7b1SBarry Smith     if (flg1) {
139928b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-tellmycell requires either the phone number or number,\"message\"");
1400763ec7b1SBarry Smith       if (nmax == 1) {
14019566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(128, &buffs[1]));
14029566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
14039566063dSJacob Faibussowitsch         PetscCall(PetscStrcat(buffs[1], " has completed"));
1404763ec7b1SBarry Smith       }
14059566063dSJacob Faibussowitsch       PetscCall(PetscTellMyCell(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
14069566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
14079566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
1408763ec7b1SBarry Smith     }
14099566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
1410763ec7b1SBarry Smith   }
141104102261SBarry Smith #endif
141204102261SBarry Smith 
14132d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14149566063dSJacob Faibussowitsch   PetscCall(PetscFPTDestroy());
14152d53ad75SBarry Smith #endif
14162d53ad75SBarry Smith 
1417e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1418dff31646SBarry Smith   flg = PETSC_FALSE;
14199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL));
14201baa6e33SBarry Smith   if (flg) PetscCall(PetscOptionsSAWsDestroy());
1421d5649816SBarry Smith #endif
1422d5649816SBarry Smith 
1423681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1424681455b2SBarry Smith   flg1 = PETSC_FALSE;
14259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL));
1426681455b2SBarry Smith   if (flg1) {
1427681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
14289566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL));
1429681455b2SBarry Smith   }
1430681455b2SBarry Smith #endif
1431681455b2SBarry Smith 
143267584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_info", &flg2, NULL));
1434e5c89e4eSSatish Balay   if (!flg2) {
143590d69ab7SBarry Smith     flg2 = PETSC_FALSE;
14369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL));
1437e5c89e4eSSatish Balay   }
143848a46eb9SPierre Jolivet   if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n"));
143967584ceeSBarry Smith #endif
1440e5c89e4eSSatish Balay 
1441e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
144290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL));
1444e5c89e4eSSatish Balay   if (flg1) {
1445e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
14469566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD));
14479566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops));
1448e5c89e4eSSatish Balay   }
1449e5c89e4eSSatish Balay #endif
1450e5c89e4eSSatish Balay 
1451e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1452e5c89e4eSSatish Balay   #if defined(PETSC_HAVE_MPE)
1453e5c89e4eSSatish Balay   mname[0] = 0;
14549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1));
1455e5c89e4eSSatish Balay   if (flg1) {
14569566063dSJacob Faibussowitsch     if (mname[0]) PetscCall(PetscLogMPEDump(mname));
14579566063dSJacob Faibussowitsch     else PetscCall(PetscLogMPEDump(0));
1458e5c89e4eSSatish Balay   }
1459e5c89e4eSSatish Balay   #endif
1460356e5837SBarry Smith #endif
1461a297a907SKarl Rupp 
1462dd710f27SBarry Smith   /*
1463dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1464dd710f27SBarry Smith   */
14659566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1466dd710f27SBarry Smith 
1467356e5837SBarry Smith #if defined(PETSC_USE_LOG)
14689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE));
14699566063dSJacob Faibussowitsch   PetscCall(PetscLogViewFromOptions());
14709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPopGetViewerOff());
147187c3beb6SLisandro Dalcin 
1472356e5837SBarry Smith   mname[0] = 0;
14739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_summary", mname, sizeof(mname), &flg1));
1474e5c89e4eSSatish Balay   if (flg1) {
147591eabc43SBarry Smith     PetscViewer viewer;
14769566063dSJacob Faibussowitsch     PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD, "\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n"));
147791eabc43SBarry Smith     if (mname[0]) {
14789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, mname, &viewer));
14799566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
14809566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
148133f85c2fSBarry Smith     } else {
148233f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
14839566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
14849566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
14859566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
148633f85c2fSBarry Smith     }
1487e5c89e4eSSatish Balay   }
1488a297a907SKarl Rupp 
1489dd710f27SBarry Smith   /*
1490dd710f27SBarry Smith      Free any objects created by the last block of code.
1491dd710f27SBarry Smith   */
14929566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1493dd710f27SBarry Smith 
1494dd710f27SBarry Smith   mname[0] = 0;
14959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1));
14969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2));
14979566063dSJacob Faibussowitsch   if (flg1 || flg2) PetscCall(PetscLogDump(mname));
1498e5c89e4eSSatish Balay #endif
149910463e74SBarry Smith 
150090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL));
15029566063dSJacob Faibussowitsch   if (!flg1) PetscCall(PetscPopSignalHandler());
150390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL));
15051baa6e33SBarry Smith   if (flg1) PetscCall(PetscMPIDump(stdout));
150690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
150790d69ab7SBarry Smith   flg2 = PETSC_FALSE;
15088bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
15099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1));
15109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
15119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL));
1512e4c476e2SSatish Balay 
1513e5c89e4eSSatish Balay   if (flg2) {
1514be56827dSJed Brown     PetscViewer viewer;
15159566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15169566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15179566063dSJacob Faibussowitsch     PetscCall(PetscOptionsView(NULL, viewer));
15189566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1519e5c89e4eSSatish Balay   }
1520e5c89e4eSSatish Balay 
1521e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
15229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1));
15239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1));
1524e5c89e4eSSatish Balay 
152533fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
15269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1));
15279245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1528e5c89e4eSSatish Balay   if (flg3) {
15293de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1530be56827dSJed Brown       PetscViewer viewer;
15319566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15329566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15339566063dSJacob Faibussowitsch       PetscCall(PetscOptionsView(NULL, viewer));
15349566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
1535e5c89e4eSSatish Balay     }
15369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsAllUsed(NULL, &nopt));
15373de2bfdfSBarry Smith     if (nopt) {
15389566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n"));
15399566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n"));
15403de2bfdfSBarry Smith       if (nopt == 1) {
15419566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n"));
1542e5c89e4eSSatish Balay       } else {
15439566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt));
1544e5c89e4eSSatish Balay       }
15453de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15469566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n"));
1547df12ba86SBarry Smith     }
15489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsLeft(NULL));
1549e5c89e4eSSatish Balay   }
1550e5c89e4eSSatish Balay 
1551e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1552d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15539566063dSJacob Faibussowitsch     PetscCall(PetscStackSAWsViewOff());
1554792fecdfSBarry Smith     PetscCallSAWs(SAWs_Finalize, ());
1555d45a07a7SBarry Smith   }
1556ec957eceSBarry Smith #endif
1557ec957eceSBarry Smith 
15584097062eSBarry Smith #if defined(PETSC_USE_LOG)
155910463e74SBarry Smith   /*
1560dbc8283eSBarry Smith        List all objects the user may have forgot to free
15612eff7a51SBarry Smith   */
156205df10baSBarry Smith   if (PetscObjectsLog) {
15639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
1564a64a8e02SBarry Smith     if (flg1) {
1565a64a8e02SBarry Smith       MPI_Comm local_comm;
15667eb1d149SBarry Smith       char     string[64];
1567a64a8e02SBarry Smith 
15689566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL));
1569afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
15709566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
15719566063dSJacob Faibussowitsch       PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
15729566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
15739566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
15740a1571b3SBarry Smith     }
157505df10baSBarry Smith   }
15764097062eSBarry Smith #endif
15774097062eSBarry Smith 
15784097062eSBarry Smith #if defined(PETSC_USE_LOG)
1579dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1580dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
15819566063dSJacob Faibussowitsch   PetscCall(PetscFree(PetscObjects));
15824097062eSBarry Smith #endif
15832eff7a51SBarry Smith 
158433f85c2fSBarry Smith   /*
158533f85c2fSBarry Smith      Destroy any packages that registered a finalize
158633f85c2fSBarry Smith   */
15879566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalizeAll());
158833f85c2fSBarry Smith 
1589101409b8SToby Isaac #if defined(PETSC_USE_LOG)
15909566063dSJacob Faibussowitsch   PetscCall(PetscLogFinalize());
1591101409b8SToby Isaac #endif
1592101409b8SToby Isaac 
159333f85c2fSBarry Smith   /*
159448dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
159548dd1dffSBarry Smith   */
15962e956fe4SStefano Zampini   if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll());
159737e93019SBarry Smith 
15984028d114SSatish Balay   if (petsc_history) {
15999566063dSJacob Faibussowitsch     PetscCall(PetscCloseHistoryFile(&petsc_history));
160002c9f0b5SLisandro Dalcin     petsc_history = NULL;
1601e5c89e4eSSatish Balay   }
16029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
16039566063dSJacob Faibussowitsch   PetscCall(PetscInfoDestroy());
1604e5c89e4eSSatish Balay 
160567584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
160692f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1607e5c89e4eSSatish Balay     char  fname[PETSC_MAX_PATH_LEN];
160892f119d6SBarry Smith     char  sname[PETSC_MAX_PATH_LEN];
1609e5c89e4eSSatish Balay     FILE *fd;
1610ed9cf6e9SBarry Smith     int   err;
1611e5c89e4eSSatish Balay 
1612dc92acbaSJed Brown     flg2 = PETSC_FALSE;
161392f119d6SBarry Smith     flg3 = PETSC_FALSE;
16149566063dSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL));
16159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL));
161692f119d6SBarry Smith     fname[0] = 0;
16179566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1));
1618e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1619589a23caSBarry Smith       PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank);
16209371c9d4SSatish Balay       fd = fopen(sname, "w");
16219371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16229566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(fd));
1623ed9cf6e9SBarry Smith       err = fclose(fd);
162428b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
162592f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1626e5c89e4eSSatish Balay       MPI_Comm local_comm;
1627e5c89e4eSSatish Balay 
1628afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16299566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16309566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(stdout));
16319566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16329566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1633e5c89e4eSSatish Balay     }
1634e5c89e4eSSatish Balay     fname[0] = 0;
16359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1));
1636e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1637589a23caSBarry Smith       PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank);
16389371c9d4SSatish Balay       fd = fopen(sname, "w");
16399371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16409566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(fd));
1641ed9cf6e9SBarry Smith       err = fclose(fd);
164228b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
164392f119d6SBarry Smith     } else if (flg1) {
164492f119d6SBarry Smith       MPI_Comm local_comm;
164592f119d6SBarry Smith 
1646afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16479566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16489566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(stdout));
16499566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16509566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1651e5c89e4eSSatish Balay     }
1652e5c89e4eSSatish Balay   }
165367584ceeSBarry Smith #endif
165420e2c332SMatthew G. Knepley 
16555486ca60SMatthew G. Knepley   /*
16565486ca60SMatthew G. Knepley      Close any open dynamic libraries
16575486ca60SMatthew G. Knepley   */
16589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize_DynamicLibraries());
16595486ca60SMatthew G. Knepley 
1660e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroyDefault());
1662e5c89e4eSSatish Balay 
1663e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
166402c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1665e5c89e4eSSatish Balay 
1666c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1667c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
16689566063dSJacob Faibussowitsch     PetscCall(PetscKokkosFinalize_Private());
1669c2b86a48SJunchao Zhang     PetscBeganKokkos       = PETSC_FALSE;
167045639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1671c2b86a48SJunchao Zhang   }
1672c2b86a48SJunchao Zhang #endif
1673c2b86a48SJunchao Zhang 
167471438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
167571438e86SJunchao Zhang   if (PetscBeganNvshmem) {
16769566063dSJacob Faibussowitsch     PetscCall(PetscNvshmemFinalize());
167771438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
167871438e86SJunchao Zhang   }
167971438e86SJunchao Zhang #endif
1680a0e72f99SJunchao Zhang 
16819566063dSJacob Faibussowitsch   PetscCall(PetscFreeMPIResources());
1682e5c89e4eSSatish Balay 
1683dbc8283eSBarry Smith   /*
1684efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1685efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1686efb80d3cSBarry Smith 
1687efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1688efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1689dbc8283eSBarry Smith  */
1690b770b1f6SSatish Balay   {
1691dbc8283eSBarry Smith     PetscCommCounter *counter;
1692dbc8283eSBarry Smith     PetscMPIInt       flg;
1693dbc8283eSBarry Smith     MPI_Comm          icomm;
16949371c9d4SSatish Balay     union
16959371c9d4SSatish Balay     {
16969371c9d4SSatish Balay       MPI_Comm comm;
16979371c9d4SSatish Balay       void    *ptr;
16989371c9d4SSatish Balay     } ucomm;
16999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
1700dbc8283eSBarry Smith     if (flg) {
1701265f3f35SJed Brown       icomm = ucomm.comm;
17029566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
170328b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1704dbc8283eSBarry Smith 
17059566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval));
17069566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17079566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1708dbc8283eSBarry Smith     }
17099566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
1710dbc8283eSBarry Smith     if (flg) {
1711265f3f35SJed Brown       icomm = ucomm.comm;
17129566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
171328b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1714dbc8283eSBarry Smith 
17159566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval));
17169566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17179566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1718dbc8283eSBarry Smith     }
1719b770b1f6SSatish Balay   }
1720dbc8283eSBarry Smith 
17219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
17229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
17239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
17249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
1725480cf27aSJed Brown 
17269566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
17279566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
17289566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
17299566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock));
1730ef19f930SBarry Smith 
1731e5c89e4eSSatish Balay   if (PetscBeganMPI) {
173299b1327fSBarry Smith     PetscMPIInt flag;
17339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalized(&flag));
173439a651e2SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
173539a651e2SJacob Faibussowitsch     /* wait until the very last moment to disable error handling */
173639a651e2SJacob Faibussowitsch     PetscErrorHandlingInitialized = PETSC_FALSE;
17379566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalize());
173839a651e2SJacob Faibussowitsch   } else PetscErrorHandlingInitialized = PETSC_FALSE;
173939a651e2SJacob Faibussowitsch 
1740e5c89e4eSSatish Balay   /*
1741e5c89e4eSSatish Balay 
1742e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1743e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1744e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1745e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1746e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17470e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1748e5c89e4eSSatish Balay    memory was not freed.
1749e5c89e4eSSatish Balay 
1750e5c89e4eSSatish Balay */
17519566063dSJacob Faibussowitsch   PetscCall(PetscMallocClear());
17529566063dSJacob Faibussowitsch   PetscCall(PetscStackReset());
1753a297a907SKarl Rupp 
1754e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1755e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
175656883f60SBarry Smith #if defined(PETSC_USE_GCOV)
175756883f60SBarry Smith   /*
175856883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
175956883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
176056883f60SBarry Smith    */
176156883f60SBarry Smith   __gcov_flush();
176256883f60SBarry Smith #endif
17631724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
17641724198aSStefano Zampini   PetscStackClearTop;
176527104ee2SJacob Faibussowitsch   return 0;
1766e5c89e4eSSatish Balay }
1767e5c89e4eSSatish Balay 
176843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
1769d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b)
1770d71ae5a4SJacob Faibussowitsch {
177143db4dbbSBarry Smith   if (*a == *b) return 1;
177243db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
177343db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
177443db4dbbSBarry Smith   return 0;
177543db4dbbSBarry Smith }
1776a70650f6SBarry Smith #endif
177743db4dbbSBarry Smith 
177843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
1779d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b)
1780d71ae5a4SJacob Faibussowitsch {
178143db4dbbSBarry Smith   if (*a == *b) return 1;
178243db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
178343db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
178443db4dbbSBarry Smith   return 0;
178543db4dbbSBarry Smith }
178643db4dbbSBarry Smith #endif
1787