xref: /petsc/src/sys/objects/pinit.c (revision 9a2cd68e563fb69dfcb2578ce58f7db064c2e0fb)
1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I  "petscsys.h"   I*/
6665c2dedSJed Brown #include <petscviewer.h>
762e5d2d2SJDBetteridge #include <petsc/private/garbagecollector.h>
8a0e72f99SJunchao Zhang 
96edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
106edef35eSSatish Balay   #include <petsc/private/valgrind/valgrind.h>
116edef35eSSatish Balay #endif
126edef35eSSatish Balay 
1385649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
1485649d77SJunchao Zhang   #include <petsc/private/fortranimpl.h>
1585649d77SJunchao Zhang #endif
1685649d77SJunchao Zhang 
177ce81a4bSJacob Faibussowitsch #if PetscDefined(USE_COVERAGE)
1856883f60SBarry Smith EXTERN_C_BEGIN
19aaf3846cSSatish Balay   #if defined(PETSC_HAVE___GCOV_DUMP)
20aaf3846cSSatish Balay     #define __gcov_flush(x) __gcov_dump(x)
21aaf3846cSSatish Balay   #endif
2256883f60SBarry Smith void __gcov_flush(void);
2356883f60SBarry Smith EXTERN_C_END
2456883f60SBarry Smith #endif
258101f56cSMatthew Knepley 
262d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
2795c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
282d53ad75SBarry Smith PetscFPT              PetscFPTData = 0;
292d53ad75SBarry Smith #endif
302d53ad75SBarry Smith 
3127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
3216ad0300SBarry Smith   #include <petscviewersaws.h>
33a6790183SBarry Smith #endif
34eb02082bSJunchao Zhang 
3595c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
36e5c89e4eSSatish Balay 
3795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm, int);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm, int);
4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE **);
420069ddf5SShri Abhyankar 
436de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
44e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4527104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
466de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
476de5d289SStefano Zampini #else
486de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
496de5d289SStefano Zampini #endif
50e5c89e4eSSatish Balay 
51480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval      = MPI_KEYVAL_INVALID;
52480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval    = MPI_KEYVAL_INVALID;
53480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval    = MPI_KEYVAL_INVALID;
545f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval      = MPI_KEYVAL_INVALID;
5562e5d2d2SJDBetteridge PetscMPIInt Petsc_CreationIdx_keyval  = MPI_KEYVAL_INVALID;
5662e5d2d2SJDBetteridge PetscMPIInt Petsc_Garbage_HMap_keyval = MPI_KEYVAL_INVALID;
57480cf27aSJed Brown 
585ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedWD_keyval  = MPI_KEYVAL_INVALID;
595ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedTmp_keyval = MPI_KEYVAL_INVALID;
605ea2b939SDuncan Campbell 
61e5c89e4eSSatish Balay /*
62e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
63e5c89e4eSSatish Balay */
6402c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL};
6502c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL};
66e5c89e4eSSatish Balay 
67ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
68ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
690f8e0872SSatish Balay 
70a2f94806SJed Brown PetscInt PetscHotRegionDepth;
71a2f94806SJed Brown 
726edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
736edef35eSSatish Balay 
74b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
75b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
76b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
77b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
78b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
79b22622e2STadeu Manoel #endif
80b22622e2STadeu Manoel 
81e5c89e4eSSatish Balay /*
82945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8372a42c3cSBarry Smith 
8472a42c3cSBarry Smith    Collective
8572a42c3cSBarry Smith 
8672a42c3cSBarry Smith    Level: advanced
8772a42c3cSBarry Smith 
8895452b02SPatrick Sanan     Notes:
89a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
900f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
91a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
920f11a792SBarry Smith 
93a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
941ea5a559SBarry Smith 
95db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()`
960f11a792SBarry Smith */
97d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoPointers(int argc, char **args, const char *filename, const char *help)
98d71ae5a4SJacob Faibussowitsch {
9972a42c3cSBarry Smith   int    myargc = argc;
10072a42c3cSBarry Smith   char **myargs = args;
10172a42c3cSBarry Smith 
10272a42c3cSBarry Smith   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&myargc, &myargs, filename, help));
1049566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
105df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10772a42c3cSBarry Smith }
10872a42c3cSBarry Smith 
109f0865b08SBarry Smith /*
110a56f64adSBarry Smith       Used by Julia interface to get communicator
111f0865b08SBarry Smith */
112d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
113d71ae5a4SJacob Faibussowitsch {
114f0865b08SBarry Smith   PetscFunctionBegin;
11527104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidPointer(comm, 1);
116f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118f0865b08SBarry Smith }
119f0865b08SBarry Smith 
120e5c89e4eSSatish Balay /*@C
121811af0c4SBarry Smith       PetscInitializeNoArguments - Calls `PetscInitialize()` from C/C++ without
122e5c89e4eSSatish Balay         the command line arguments.
123e5c89e4eSSatish Balay 
124e5c89e4eSSatish Balay    Collective
125e5c89e4eSSatish Balay 
126e5c89e4eSSatish Balay    Level: advanced
127e5c89e4eSSatish Balay 
128db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`
129e5c89e4eSSatish Balay @*/
130d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoArguments(void)
131d71ae5a4SJacob Faibussowitsch {
132e5c89e4eSSatish Balay   int    argc = 0;
13302c9f0b5SLisandro Dalcin   char **args = NULL;
134e5c89e4eSSatish Balay 
135e5c89e4eSSatish Balay   PetscFunctionBegin;
1369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, NULL, NULL));
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138e5c89e4eSSatish Balay }
139e5c89e4eSSatish Balay 
140e5c89e4eSSatish Balay /*@
141e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
142e5c89e4eSSatish Balay 
14393b6d2d1SJed Brown    Level: beginner
144e5c89e4eSSatish Balay 
145db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
146e5c89e4eSSatish Balay @*/
147d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialized(PetscBool *isInitialized)
148d71ae5a4SJacob Faibussowitsch {
14927104ee2SJacob Faibussowitsch   PetscFunctionBegin;
15027104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized, 1);
151e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153e5c89e4eSSatish Balay }
154e5c89e4eSSatish Balay 
155e5c89e4eSSatish Balay /*@
156811af0c4SBarry Smith       PetscFinalized - Determine whether `PetscFinalize()` has been called yet
157e5c89e4eSSatish Balay 
158e5c89e4eSSatish Balay    Level: developer
159e5c89e4eSSatish Balay 
160db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
161e5c89e4eSSatish Balay @*/
162d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalized(PetscBool *isFinalized)
163d71ae5a4SJacob Faibussowitsch {
16427104ee2SJacob Faibussowitsch   PetscFunctionBegin;
16527104ee2SJacob Faibussowitsch   if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized, 1);
166e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168e5c89e4eSSatish Balay }
169e5c89e4eSSatish Balay 
17057171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]);
171e5c89e4eSSatish Balay 
172e5c89e4eSSatish Balay /*
173e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
174e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
175e5c89e4eSSatish Balay */
176367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP               = 0;
17762e5d2d2SJDBetteridge MPI_Op Petsc_Garbage_SetIntersectOp = 0;
178e5c89e4eSSatish Balay 
179d71ae5a4SJacob Faibussowitsch PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, int *cnt, MPI_Datatype *datatype)
180d71ae5a4SJacob Faibussowitsch {
181e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out, i, count = *cnt;
182e5c89e4eSSatish Balay 
183e5c89e4eSSatish Balay   PetscFunctionBegin;
184e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
1853ba16761SJacob Faibussowitsch     PetscErrorCode ierr = (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
1863ba16761SJacob Faibussowitsch     (void)ierr;
18741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
188e5c89e4eSSatish Balay   }
189e5c89e4eSSatish Balay 
190e5c89e4eSSatish Balay   for (i = 0; i < count; i++) {
191e5c89e4eSSatish Balay     xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]);
192e5c89e4eSSatish Balay     xout[2 * i + 1] += xin[2 * i + 1];
193e5c89e4eSSatish Balay   }
194812af9f3SBarry Smith   PetscFunctionReturnVoid();
195e5c89e4eSSatish Balay }
196e5c89e4eSSatish Balay 
197e5c89e4eSSatish Balay /*
198e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
199e5c89e4eSSatish Balay sum of the second entry.
200b693b147SBarry Smith 
20176f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
202367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
203b693b147SBarry Smith there would be no place to store the both needed results.
204e5c89e4eSSatish Balay */
205d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt sizes[], PetscInt *max, PetscInt *sum)
206d71ae5a4SJacob Faibussowitsch {
207e5c89e4eSSatish Balay   PetscFunctionBegin;
208d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
209d6e4c47cSJed Brown   {
2109371c9d4SSatish Balay     struct {
2119371c9d4SSatish Balay       PetscInt max, sum;
2129371c9d4SSatish Balay     } work;
2139566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce_scatter_block((void *)sizes, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm));
214d6e4c47cSJed Brown     *max = work.max;
215d6e4c47cSJed Brown     *sum = work.sum;
216d6e4c47cSJed Brown   }
217d6e4c47cSJed Brown #else
218d6e4c47cSJed Brown   {
219d6e4c47cSJed Brown     PetscMPIInt size, rank;
2209371c9d4SSatish Balay     struct {
2219371c9d4SSatish Balay       PetscInt max, sum;
2229371c9d4SSatish Balay     } *work;
2239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
2249566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
2259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &work));
2261c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce((void *)sizes, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm));
2276ac3741eSJed Brown     *max = work[rank].max;
2286ac3741eSJed Brown     *sum = work[rank].sum;
2299566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
230d6e4c47cSJed Brown   }
231d6e4c47cSJed Brown #endif
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
233e5c89e4eSSatish Balay }
234e5c89e4eSSatish Balay 
2359e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
236613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
237613bf2b2SPierre Jolivet     #include <quadmath.h>
238613bf2b2SPierre Jolivet   #endif
2399e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0;
240de272c7aSSatish Balay   #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
24106a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
242613bf2b2SPierre Jolivet   #endif
243e5c89e4eSSatish Balay 
244d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
245d71ae5a4SJacob Faibussowitsch {
246e5c89e4eSSatish Balay   PetscInt i, count = *cnt;
247e5c89e4eSSatish Balay 
248e5c89e4eSSatish Balay   PetscFunctionBegin;
2497c2de775SJed Brown   if (*datatype == MPIU_REAL) {
250e2e03761SBarry Smith     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
251a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2527c2de775SJed Brown   }
2537c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
254c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2557c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
256a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2577c2de775SJed Brown   }
2587c2de775SJed Brown   #endif
259613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
260613bf2b2SPierre Jolivet   else if (*datatype == MPIU___FLOAT128) {
261613bf2b2SPierre Jolivet     __float128 *xin = (__float128 *)in, *xout = (__float128 *)out;
262613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
26374c01ffaSSatish Balay     #if defined(PETSC_HAVE_COMPLEX)
2649371c9d4SSatish Balay   } else if (*datatype == MPIU___COMPLEX128) {
265613bf2b2SPierre Jolivet     __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out;
266613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
26774c01ffaSSatish Balay     #endif
268613bf2b2SPierre Jolivet   }
269613bf2b2SPierre Jolivet   #endif
2709e517322SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FP16)
2719e517322SPierre Jolivet   else if (*datatype == MPIU___FP16) {
2729e517322SPierre Jolivet     __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out;
2739e517322SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2749e517322SPierre Jolivet   }
2759e517322SPierre Jolivet   #endif
2767c2de775SJed Brown   else {
2779e517322SPierre Jolivet   #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16)
2783ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SElF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
2799e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FP16)
2803ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types"));
2819e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FLOAT128)
2823ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types"));
2839e517322SPierre Jolivet   #else
2843ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types"));
285613bf2b2SPierre Jolivet   #endif
28641e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
287e2e03761SBarry Smith   }
288812af9f3SBarry Smith   PetscFunctionReturnVoid();
289e5c89e4eSSatish Balay }
290e5c89e4eSSatish Balay #endif
291e5c89e4eSSatish Balay 
292570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
293d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
294d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
295d9822059SBarry Smith 
296d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
297d71ae5a4SJacob Faibussowitsch {
298d9822059SBarry Smith   PetscInt i, count = *cnt;
299d9822059SBarry Smith 
300d9822059SBarry Smith   PetscFunctionBegin;
3017c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3028c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
303a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]);
3047c2de775SJed Brown   }
3057c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3067c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3077c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
308ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3097c2de775SJed Brown   }
3107c2de775SJed Brown   #endif
3117c2de775SJed Brown   else {
3123ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
31341e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3148c764dc5SJose Roman   }
315d9822059SBarry Smith   PetscFunctionReturnVoid();
316d9822059SBarry Smith }
317d9822059SBarry Smith 
318d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
319d71ae5a4SJacob Faibussowitsch {
320d9822059SBarry Smith   PetscInt i, count = *cnt;
321d9822059SBarry Smith 
322d9822059SBarry Smith   PetscFunctionBegin;
3237c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3248c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
325a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]);
3267c2de775SJed Brown   }
3277c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3287c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3297c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
330ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3317c2de775SJed Brown   }
3327c2de775SJed Brown   #endif
3337c2de775SJed Brown   else {
3343ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"));
33541e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3368c764dc5SJose Roman   }
337d9822059SBarry Smith   PetscFunctionReturnVoid();
338d9822059SBarry Smith }
339d9822059SBarry Smith #endif
340d9822059SBarry Smith 
341480cf27aSJed Brown /*
342480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
343480cf27aSJed Brown 
344ff0e51ddSBarry 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.
345480cf27aSJed Brown 
34612801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
347480cf27aSJed Brown 
348480cf27aSJed Brown */
349d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state)
350d71ae5a4SJacob Faibussowitsch {
35105342407SJunchao Zhang   PetscCommCounter      *counter = (PetscCommCounter *)count_val;
35257f21012SBarry Smith   struct PetscCommStash *comms   = counter->comms, *pcomm;
353480cf27aSJed Brown 
354480cf27aSJed Brown   PetscFunctionBegin;
3559566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm));
3569566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter->iflags));
35757f21012SBarry Smith   while (comms) {
3589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&comms->comm));
35957f21012SBarry Smith     pcomm = comms;
36057f21012SBarry Smith     comms = comms->next;
3619566063dSJacob Faibussowitsch     PetscCall(PetscFree(pcomm));
36257f21012SBarry Smith   }
3639566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter));
364480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
365480cf27aSJed Brown }
366480cf27aSJed Brown 
367480cf27aSJed Brown /*
36847435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
369da3039f7SJed Brown   calls MPI_Comm_free().
370da3039f7SJed Brown 
371da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
372480cf27aSJed Brown 
373ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
374480cf27aSJed Brown 
37512801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
376480cf27aSJed Brown 
377480cf27aSJed Brown */
378d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
379d71ae5a4SJacob Faibussowitsch {
3809371c9d4SSatish Balay   union
381480cf27aSJed Brown   {
3829371c9d4SSatish Balay     MPI_Comm comm;
3839371c9d4SSatish Balay     void    *ptr;
3849371c9d4SSatish Balay   } icomm;
385480cf27aSJed Brown 
386480cf27aSJed Brown   PetscFunctionBegin;
38712801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval");
388ec4fadc2SJed Brown   icomm.ptr = attr_val;
38976bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
39033779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
39133779a13SJunchao Zhang     PetscMPIInt flg;
3929371c9d4SSatish Balay     union
3939371c9d4SSatish Balay     {
3949371c9d4SSatish Balay       MPI_Comm comm;
3959371c9d4SSatish Balay       void    *ptr;
3969371c9d4SSatish Balay     } ocomm;
3979566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg));
39833779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute");
39933779a13SJunchao 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");
40033779a13SJunchao Zhang   }
4019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval));
4029566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm));
403da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
404b89831e5SBarry Smith }
405da3039f7SJed Brown 
406da3039f7SJed Brown /*
40733779a13SJunchao 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.
408da3039f7SJed Brown  */
409d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
410d71ae5a4SJacob Faibussowitsch {
411da3039f7SJed Brown   PetscFunctionBegin;
4129566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm));
413480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
414480cf27aSJed Brown }
415480cf27aSJed Brown 
41633779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *);
41742218b76SBarry Smith 
418951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
4198cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *);
4208cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
4218cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
422e39fd77fSBarry Smith #endif
423e39fd77fSBarry Smith 
4240f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE;
42512801b39SBarry Smith 
426eb27c7c8SSatish Balay PETSC_INTERN int    PetscGlobalArgc;
427eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4286ae9a8a6SBarry Smith int                 PetscGlobalArgc = 0;
42902c9f0b5SLisandro Dalcin char              **PetscGlobalArgs = NULL;
430dff31646SBarry Smith PetscSegBuffer      PetscCitationsList;
431e5c89e4eSSatish Balay 
432d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void)
433d71ae5a4SJacob Faibussowitsch {
434051e4cf2SJed Brown   PetscFunctionBegin;
4359566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList));
43630c35bf2SSatish Balay 
43730c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\
43830c35bf2SSatish Balay   Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\
43930c35bf2SSatish Balay     and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\
4403f81df79SBarry Smith     and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\
44130c35bf2SSatish Balay     and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\
44230c35bf2SSatish Balay     and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\
44330c35bf2SSatish Balay     and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\
44430c35bf2SSatish Balay     and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\
44530c35bf2SSatish Balay   Title = {{PETSc/TAO} Users Manual},\n\
446c4dc7257SSatish Balay   Number = {ANL-21/39 - Revision 3.19},\n\
447*9a2cd68eSMatthew Knepley   Doi = {10.2172/1968587},\n\
44830c35bf2SSatish Balay   Institution = {Argonne National Laboratory},\n\
449c4dc7257SSatish Balay   Year = {2023}\n}\n",
4509371c9d4SSatish Balay                                    NULL));
45130c35bf2SSatish Balay 
45230c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\
45330c35bf2SSatish Balay   Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\
45430c35bf2SSatish Balay   Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\
45530c35bf2SSatish Balay   Booktitle = {Modern Software Tools in Scientific Computing},\n\
45630c35bf2SSatish Balay   Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\
45730c35bf2SSatish Balay   Pages = {163--202},\n\
45830c35bf2SSatish Balay   Publisher = {Birkh{\\\"{a}}user Press},\n\
4599371c9d4SSatish Balay   Year = {1997}\n}\n",
4609371c9d4SSatish Balay                                    NULL));
46130c35bf2SSatish Balay 
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
463051e4cf2SJed Brown }
464e5c89e4eSSatish Balay 
4652d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4662d747510SLisandro Dalcin 
467d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[])
468d71ae5a4SJacob Faibussowitsch {
4692d747510SLisandro Dalcin   PetscFunctionBegin;
4709566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(programname, name, sizeof(programname)));
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4722d747510SLisandro Dalcin }
4732d747510SLisandro Dalcin 
4742d747510SLisandro Dalcin /*@C
4752d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4762d747510SLisandro Dalcin 
4772d747510SLisandro Dalcin     Not Collective
4782d747510SLisandro Dalcin 
4792d747510SLisandro Dalcin     Input Parameter:
4802d747510SLisandro Dalcin .   len - length of the string name
4812d747510SLisandro Dalcin 
4822d747510SLisandro Dalcin     Output Parameter:
483811af0c4SBarry Smith .   name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN`
4842d747510SLisandro Dalcin 
4852d747510SLisandro Dalcin    Level: advanced
4862d747510SLisandro Dalcin 
48721532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()`
4882d747510SLisandro Dalcin @*/
489d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProgramName(char name[], size_t len)
490d71ae5a4SJacob Faibussowitsch {
4912d747510SLisandro Dalcin   PetscFunctionBegin;
4929566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(name, programname, len));
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4942d747510SLisandro Dalcin }
4952d747510SLisandro Dalcin 
496e5c89e4eSSatish Balay /*@C
497e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
498811af0c4SBarry Smith      after PetscInitialize() is called but before `PetscFinalize()`.
499e5c89e4eSSatish Balay 
500e5c89e4eSSatish Balay    Not Collective
501e5c89e4eSSatish Balay 
502e5c89e4eSSatish Balay    Output Parameters:
503e5c89e4eSSatish Balay +  argc - count of number of command line arguments
504e5c89e4eSSatish Balay -  args - the command line arguments
505e5c89e4eSSatish Balay 
506e5c89e4eSSatish Balay    Level: intermediate
507e5c89e4eSSatish Balay 
508e5c89e4eSSatish Balay    Notes:
509e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
510e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
511e5c89e4eSSatish Balay 
512f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
513f177e3b1SBarry Smith 
51421532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()`
515e5c89e4eSSatish Balay @*/
516d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args)
517d71ae5a4SJacob Faibussowitsch {
518e5c89e4eSSatish Balay   PetscFunctionBegin;
51939a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
520e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
521e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
523e5c89e4eSSatish Balay }
524e5c89e4eSSatish Balay 
525793721a6SBarry Smith /*@C
526793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
527811af0c4SBarry Smith      after `PetscInitialize()` is called but before `PetscFinalize()`.
528793721a6SBarry Smith 
529793721a6SBarry Smith    Not Collective
530793721a6SBarry Smith 
5312fe279fdSBarry Smith    Output Parameter:
532793721a6SBarry Smith .  args - the command line arguments
533793721a6SBarry Smith 
534793721a6SBarry Smith    Level: intermediate
535793721a6SBarry Smith 
53621532e8aSBarry Smith    Note:
53721532e8aSBarry Smith       This does NOT start with the program name and IS `NULL` terminated (final arg is void)
538793721a6SBarry Smith 
53921532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`, `PetscInitialize()`
540793721a6SBarry Smith @*/
541d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args)
542d71ae5a4SJacob Faibussowitsch {
543793721a6SBarry Smith   PetscInt i, argc = PetscGlobalArgc;
544793721a6SBarry Smith 
545793721a6SBarry Smith   PetscFunctionBegin;
54639a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
5479371c9d4SSatish Balay   if (!argc) {
5489371c9d4SSatish Balay     *args = NULL;
5493ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5509371c9d4SSatish Balay   }
5519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(argc, args));
5529566063dSJacob Faibussowitsch   for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i]));
5532d747510SLisandro Dalcin   (*args)[argc - 1] = NULL;
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
555793721a6SBarry Smith }
556793721a6SBarry Smith 
557793721a6SBarry Smith /*@C
558811af0c4SBarry Smith    PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()`
559793721a6SBarry Smith 
560793721a6SBarry Smith    Not Collective
561793721a6SBarry Smith 
5622fe279fdSBarry Smith    Output Parameter:
563793721a6SBarry Smith .  args - the command line arguments
564793721a6SBarry Smith 
565793721a6SBarry Smith    Level: intermediate
566793721a6SBarry Smith 
567db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()`
568793721a6SBarry Smith @*/
569d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args)
570d71ae5a4SJacob Faibussowitsch {
57139a651e2SJacob Faibussowitsch   PetscFunctionBegin;
57239a651e2SJacob Faibussowitsch   if (args) {
573793721a6SBarry Smith     PetscInt i = 0;
574793721a6SBarry Smith 
5759566063dSJacob Faibussowitsch     while (args[i]) PetscCall(PetscFree(args[i++]));
5769566063dSJacob Faibussowitsch     PetscCall(PetscFree(args));
57739a651e2SJacob Faibussowitsch   }
5783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
579793721a6SBarry Smith }
580793721a6SBarry Smith 
58127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
58230befbd2SBarry Smith   #include <petscconfiginfo.h>
58330befbd2SBarry Smith 
584d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
585d71ae5a4SJacob Faibussowitsch {
58627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
58711525c0dSBarry Smith   if (!PetscGlobalRank) {
58830befbd2SBarry Smith     char      cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64];
58911525c0dSBarry Smith     int       port;
590ffbd1cfbSBarry Smith     PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE;
59111525c0dSBarry Smith     size_t    applinelen, introlen;
592ffbd1cfbSBarry Smith     char      sawsurl[256];
59311525c0dSBarry Smith 
5949566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg));
59511525c0dSBarry Smith     if (flg) {
59611525c0dSBarry Smith       char sawslog[PETSC_MAX_PATH_LEN];
59711525c0dSBarry Smith 
5989566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL));
59911525c0dSBarry Smith       if (sawslog[0]) {
600792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog));
60111525c0dSBarry Smith       } else {
602792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL));
60311525c0dSBarry Smith       }
60411525c0dSBarry Smith     }
6059566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg));
60648a46eb9SPierre Jolivet     if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert));
6079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL));
608ffbd1cfbSBarry Smith     if (selectport) {
609792fecdfSBarry Smith       PetscCallSAWs(SAWs_Get_Available_Port, (&port));
610792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Port, (port));
611ffbd1cfbSBarry Smith     } else {
6129566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg));
61348a46eb9SPierre Jolivet       if (flg) PetscCallSAWs(SAWs_Set_Port, (port));
614ffbd1cfbSBarry Smith     }
6159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg));
61611525c0dSBarry Smith     if (flg) {
617792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Document_Root, (root));
6189566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(root, ".", &rootlocal));
6199c1e0ce8SBarry Smith     } else {
6209566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg));
6219c1e0ce8SBarry Smith       if (flg) {
6229566063dSJacob Faibussowitsch         PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root)));
623792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Document_Root, (root));
6249c1e0ce8SBarry Smith       }
62511525c0dSBarry Smith     }
6269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2));
62711525c0dSBarry Smith     if (flg2) {
62811525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
62928b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option");
6309566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root));
6319566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(jsdir, 'r', &flg));
63228b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory");
633792fecdfSBarry Smith       PetscCallSAWs(SAWs_Push_Local_Header, ());
63411525c0dSBarry Smith     }
6359566063dSJacob Faibussowitsch     PetscCall(PetscGetProgramName(programname, sizeof(programname)));
6369566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(help, &applinelen));
63711525c0dSBarry Smith     introlen = 4096 + applinelen;
63830a8c9c0SSurtai Han     applinelen += 1024;
6399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(applinelen, &appline));
6409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(introlen, &intro));
64111525c0dSBarry Smith 
64211525c0dSBarry Smith     if (rootlocal) {
6439566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname));
6449566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(appline, 'r', &rootlocal));
64511525c0dSBarry Smith     }
6469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetAll(NULL, &options));
64711525c0dSBarry Smith     if (rootlocal && help) {
6489566063dSJacob 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));
64911525c0dSBarry Smith     } else if (help) {
6509566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help));
65111525c0dSBarry Smith     } else {
6529566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options));
65311525c0dSBarry Smith     }
6549566063dSJacob Faibussowitsch     PetscCall(PetscFree(options));
6559566063dSJacob Faibussowitsch     PetscCall(PetscGetVersion(version, sizeof(version)));
6569371c9d4SSatish Balay     PetscCall(PetscSNPrintf(intro, introlen,
6579371c9d4SSatish Balay                             "<body>\n"
658a17b96a8SKyle 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"
659df62c222SBarry 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"
6609371c9d4SSatish Balay                             "%s",
6619371c9d4SSatish Balay                             version, petscconfigureoptions, appline));
662792fecdfSBarry Smith     PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro));
6639566063dSJacob Faibussowitsch     PetscCall(PetscFree(intro));
6649566063dSJacob Faibussowitsch     PetscCall(PetscFree(appline));
665ffbd1cfbSBarry Smith     if (selectport) {
666aa573868SBarry Smith       PetscBool silent;
6677d812c46SBarry Smith 
6687d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
66939a651e2SJacob Faibussowitsch       while (SAWs_Initialize()) {
670792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_Available_Port, (&port));
671792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port, (port));
6727d812c46SBarry Smith       }
6737d812c46SBarry Smith 
6749566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL));
675aa573868SBarry Smith       if (!silent) {
676792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl));
6779566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl));
678ffbd1cfbSBarry Smith       }
6797d812c46SBarry Smith     } else {
680792fecdfSBarry Smith       PetscCallSAWs(SAWs_Initialize, ());
681aa573868SBarry Smith     }
6829566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister("@TechReport{ saws,\n"
6830af79b04SBarry Smith                                      "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6840af79b04SBarry Smith                                      "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6850af79b04SBarry Smith                                      "  Institution = {Argonne National Laboratory},\n"
6869371c9d4SSatish Balay                                      "  Year   = 2013\n}\n",
6879371c9d4SSatish Balay                                      NULL));
68811525c0dSBarry Smith   }
6893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69011525c0dSBarry Smith }
69111525c0dSBarry Smith #endif
69211525c0dSBarry Smith 
6934dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
694d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
695d71ae5a4SJacob Faibussowitsch {
6964dfee713SSatish Balay   PetscFunctionBegin;
6974dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6984dfee713SSatish Balay   /* see MPI.py for details on this bug */
6994dfee713SSatish Balay   (void)setenv("HWLOC_COMPONENTS", "-x86", 1);
7004dfee713SSatish Balay #endif
7013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7024dfee713SSatish Balay }
7034dfee713SSatish Balay 
704a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
705a56f64adSBarry Smith   #include <adios.h>
70622580e64SBarry Smith   #include <adios_read.h>
7077b56e58cSSatish Balay int64_t Petsc_adios_group;
708a56f64adSBarry Smith #endif
709a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
710cd1b99a6SBarry Smith   #include <omp.h>
711f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
712cd1b99a6SBarry Smith #endif
713a56f64adSBarry Smith 
714a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
715a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA)
7160e6b6b59SJacob Faibussowitsch   #include <petscdevice_cuda.h>
717a4af0ceeSJacob Faibussowitsch // REMOVE ME
718a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
719a4af0ceeSJacob Faibussowitsch #endif
720a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP)
7210e6b6b59SJacob Faibussowitsch   #include <petscdevice_hip.h>
722a4af0ceeSJacob Faibussowitsch // REMOVE ME
723a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
724a4af0ceeSJacob Faibussowitsch #endif
725a4af0ceeSJacob Faibussowitsch 
72627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
727bc8a24ecSBarry Smith   #include <dlfcn.h>
728bc8a24ecSBarry Smith #endif
729a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
730a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
731a4af0ceeSJacob Faibussowitsch #endif
732a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
7333274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscViennaCLInit(void);
734a4af0ceeSJacob Faibussowitsch PetscBool                   PetscViennaCLSynchronize = PETSC_FALSE;
735a4af0ceeSJacob Faibussowitsch #endif
736bc8a24ecSBarry Smith 
737660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE;
738660278c0SBarry Smith 
73985649d77SJunchao Zhang /*
74085649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
74185649d77SJunchao Zhang 
74285649d77SJunchao Zhang   prog:     program name
74302101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
74485649d77SJunchao Zhang   help:     program help message
745da81f932SPierre Jolivet   ftn:      is it called from Fortran initialization (petscinitializef_)?
74685649d77SJunchao Zhang   readarguments,len: used when fortran is true
74785649d77SJunchao Zhang */
748d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len)
749d71ae5a4SJacob Faibussowitsch {
75085649d77SJunchao Zhang   PetscMPIInt size;
75185649d77SJunchao Zhang   PetscBool   flg = PETSC_TRUE;
75285649d77SJunchao Zhang   char        hostname[256];
75385649d77SJunchao Zhang 
75427104ee2SJacob Faibussowitsch   PetscFunctionBegin;
7553ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
75639a651e2SJacob Faibussowitsch   /* these must be initialized in a routine, not as a constant declaration */
75739a651e2SJacob Faibussowitsch   PETSC_STDOUT = stdout;
75839a651e2SJacob Faibussowitsch   PETSC_STDERR = stderr;
75939a651e2SJacob Faibussowitsch 
7609566063dSJacob Faibussowitsch   /* PetscCall can be used from now */
76139a651e2SJacob Faibussowitsch   PetscErrorHandlingInitialized = PETSC_TRUE;
76239a651e2SJacob Faibussowitsch 
76385649d77SJunchao Zhang   /*
76485649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
76585649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
76685649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
76785649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
76885649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
76985649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
77085649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
77185649d77SJunchao Zhang       listed above and since that time are compatible.
77285649d77SJunchao Zhang 
77385649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
77485649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
77585649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
77685649d77SJunchao Zhang       to perform the checking.
77785649d77SJunchao Zhang 
77885649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
77985649d77SJunchao Zhang 
78085649d77SJunchao Zhang       Questions:
78185649d77SJunchao Zhang 
78285649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
78385649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
78485649d77SJunchao Zhang   */
78585649d77SJunchao Zhang 
78685649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
78785649d77SJunchao Zhang   {
78885649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
78985649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
79039a651e2SJacob Faibussowitsch 
7919566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength));
79285649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
79385649d77SJunchao Zhang   #if defined(MPICH_VERSION)
79485649d77SJunchao Zhang     #if MPICH_NUMVERSION < 30100000
79585649d77SJunchao Zhang     {
79685649d77SJunchao Zhang       char     *ver, *lf;
79785649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
79839a651e2SJacob Faibussowitsch 
7999566063dSJacob Faibussowitsch       PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver));
80039a651e2SJacob Faibussowitsch       if (ver) {
8019566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(ver, '\n', &lf));
80239a651e2SJacob Faibussowitsch         if (lf) {
80385649d77SJunchao Zhang           *lf = 0;
8049566063dSJacob Faibussowitsch           PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg));
80585649d77SJunchao Zhang         }
80685649d77SJunchao Zhang       }
80785649d77SJunchao Zhang       if (!flg) {
808d5b396e8SSatish Balay         PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION));
80985649d77SJunchao Zhang         flg = PETSC_TRUE;
81085649d77SJunchao Zhang       }
81185649d77SJunchao Zhang     }
81285649d77SJunchao Zhang     #endif
81385649d77SJunchao 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?) */
81485649d77SJunchao Zhang   #elif defined(OMPI_MAJOR_VERSION)
81585649d77SJunchao Zhang     {
81685649d77SJunchao Zhang       char     *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf;
81785649d77SJunchao Zhang       PetscBool flg                                              = PETSC_FALSE;
81885649d77SJunchao Zhang     #define PSTRSZ 2
81985649d77SJunchao Zhang       char      ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"};
82085649d77SJunchao Zhang       char      ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "};
82185649d77SJunchao Zhang       int       i;
82285649d77SJunchao Zhang       for (i = 0; i < PSTRSZ; i++) {
8239566063dSJacob Faibussowitsch         PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver));
82439a651e2SJacob Faibussowitsch         if (ver) {
8259566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION));
8269566063dSJacob Faibussowitsch           PetscCall(PetscStrstr(ver, bs, &bsf));
82739a651e2SJacob Faibussowitsch           if (bsf) flg = PETSC_TRUE;
82885649d77SJunchao Zhang           break;
82985649d77SJunchao Zhang         }
83085649d77SJunchao Zhang       }
83185649d77SJunchao Zhang       if (!flg) {
8323ba16761SJacob Faibussowitsch         PetscCall(PetscInfo(NULL, "PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n", mpilibraryversion, OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION));
83385649d77SJunchao Zhang         flg = PETSC_TRUE;
83485649d77SJunchao Zhang       }
83585649d77SJunchao Zhang     }
83685649d77SJunchao Zhang   #endif
83785649d77SJunchao Zhang   }
83885649d77SJunchao Zhang #endif
83985649d77SJunchao Zhang 
840e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__))
84185649d77SJunchao 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 */
84239a651e2SJacob 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");
84385649d77SJunchao Zhang #endif
84485649d77SJunchao Zhang 
84585649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
84685649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
84785649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
84885649d77SJunchao Zhang #endif
84985649d77SJunchao Zhang 
8509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
85185649d77SJunchao Zhang 
85285649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
85385649d77SJunchao Zhang 
8549566063dSJacob Faibussowitsch   PetscCall(PetscSetProgramName(prog));
8559566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen));
8569566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout));
8579566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr));
8589566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscCommSpinLock));
85985649d77SJunchao Zhang 
86085649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
8619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN));
86285649d77SJunchao Zhang 
86385649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
8649566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS));
8659566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE));
86685649d77SJunchao Zhang   }
86785649d77SJunchao Zhang 
86885649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
8699566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
87085649d77SJunchao Zhang 
8719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank));
8729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize));
87385649d77SJunchao Zhang 
87485649d77SJunchao Zhang   MPIU_BOOL        = MPI_INT;
87585649d77SJunchao Zhang   MPIU_ENUM        = MPI_INT;
87685649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64;
87785649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
87885649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
87985649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
88085649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
88185649d77SJunchao Zhang #endif
88239a651e2SJacob Faibussowitsch   else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t");
88385649d77SJunchao Zhang 
88485649d77SJunchao Zhang     /*
88585649d77SJunchao Zhang      Initialized the global complex variable; this is because with
88685649d77SJunchao Zhang      shared libraries the constructors for global variables
88785649d77SJunchao Zhang      are not called; at least on IRIX.
88885649d77SJunchao Zhang   */
88985649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
89085649d77SJunchao Zhang   {
89185649d77SJunchao Zhang   #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
89285649d77SJunchao Zhang     PetscComplex ic(0.0, 1.0);
89385649d77SJunchao Zhang     PETSC_i = ic;
89485649d77SJunchao Zhang   #else
89585649d77SJunchao Zhang     PETSC_i = _Complex_I;
89685649d77SJunchao Zhang   #endif
89785649d77SJunchao Zhang   }
89885649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
89985649d77SJunchao Zhang 
90085649d77SJunchao Zhang   /*
90185649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
90285649d77SJunchao Zhang      half of the entries and maxes the second half.
90385649d77SJunchao Zhang   */
9049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP));
90585649d77SJunchao Zhang 
906613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
9079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128));
9089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128));
90974c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
9109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128));
9119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128));
91285649d77SJunchao Zhang   #endif
91374c01ffaSSatish Balay #endif
9149e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
9159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16));
9169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FP16));
91785649d77SJunchao Zhang #endif
91885649d77SJunchao Zhang 
91985649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
9209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM));
921613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX));
922613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN));
9239e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
9249e517322SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128));
92585649d77SJunchao Zhang #endif
92685649d77SJunchao Zhang 
9279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR));
92862e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp));
9299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR));
93085649d77SJunchao Zhang 
93185649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
93285649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
93385649d77SJunchao Zhang   {
93485649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
93593d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)};
93685649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_REAL, MPIU_INT}, tmpStruct;
93785649d77SJunchao Zhang 
9389566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
93993d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT));
9409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT));
94285649d77SJunchao Zhang   }
94385649d77SJunchao Zhang   {
94485649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
94593d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)};
94685649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_SCALAR, MPIU_INT}, tmpStruct;
94785649d77SJunchao Zhang 
9489566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
94993d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT));
9509566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
95285649d77SJunchao Zhang   }
95385649d77SJunchao Zhang #endif
95485649d77SJunchao Zhang 
95585649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT));
9579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2INT));
95885649d77SJunchao Zhang #endif
9599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT));
9609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPI_4INT));
9619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT));
9629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_4INT));
96385649d77SJunchao Zhang 
96485649d77SJunchao Zhang   /*
96585649d77SJunchao Zhang      Attributes to be set on PETSc communicators
96685649d77SJunchao Zhang   */
9679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0));
9689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0));
9699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0));
9709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0));
97162e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, (void *)0));
97262e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, (void *)0));
97385649d77SJunchao Zhang 
97485649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
9759566063dSJacob Faibussowitsch   if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len));
97685649d77SJunchao Zhang   else
97785649d77SJunchao Zhang #endif
9789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file));
97985649d77SJunchao Zhang 
98085649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9819566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
98285649d77SJunchao Zhang 
98385649d77SJunchao Zhang   /*
98485649d77SJunchao Zhang      Check system options and print help
98585649d77SJunchao Zhang   */
9869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCheckInitial_Private(help));
98785649d77SJunchao Zhang 
988a4af0ceeSJacob Faibussowitsch   /*
9890e6b6b59SJacob Faibussowitsch     Creates the logging data structures; this is enabled even if logging is not turned on
9900e6b6b59SJacob Faibussowitsch     This is the last thing we do before returning to the user code to prevent having the
9910e6b6b59SJacob Faibussowitsch     logging numbers contaminated by any startup time associated with MPI
9920e6b6b59SJacob Faibussowitsch   */
9930e6b6b59SJacob Faibussowitsch #if defined(PETSC_USE_LOG)
9940e6b6b59SJacob Faibussowitsch   PetscCall(PetscLogInitialize());
9950e6b6b59SJacob Faibussowitsch #endif
9960e6b6b59SJacob Faibussowitsch 
9970e6b6b59SJacob Faibussowitsch   /*
998a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
999a4af0ceeSJacob Faibussowitsch 
1000a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
1001a4af0ceeSJacob Faibussowitsch    1. MPI initialized
1002a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
10030e6b6b59SJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up
10040e6b6b59SJacob Faibussowitsch       its own SIGSEV handler via the push/pop interface.
10050e6b6b59SJacob Faibussowitsch    4. Logging initialized
1006a4af0ceeSJacob Faibussowitsch   */
10079566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
1008a4af0ceeSJacob Faibussowitsch 
1009a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
1010a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
10119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_summary", &flg));
10129566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg));
10139566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL));
1014a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
10159566063dSJacob Faibussowitsch   PetscCall(PetscViennaCLInit());
1016a4af0ceeSJacob Faibussowitsch #endif
1017a4af0ceeSJacob Faibussowitsch 
10189566063dSJacob Faibussowitsch   PetscCall(PetscCitationsInitialize());
101985649d77SJunchao Zhang 
102085649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
10219566063dSJacob Faibussowitsch   PetscCall(PetscInitializeSAWs(ftn ? NULL : help));
102227104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
10239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg));
10249566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscStackViewSAWs());
102585649d77SJunchao Zhang #endif
102685649d77SJunchao Zhang 
102785649d77SJunchao Zhang   /*
102885649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
102985649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
103085649d77SJunchao Zhang   */
10319566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_DynamicLibraries());
103285649d77SJunchao Zhang 
10339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10349566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size));
103596dcfd6cSBarry Smith   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
10369566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname));
103785649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
103885649d77SJunchao Zhang   {
103985649d77SJunchao Zhang     PetscBool omp_view_flag;
104085649d77SJunchao Zhang     char     *threads = getenv("OMP_NUM_THREADS");
104185649d77SJunchao Zhang 
104285649d77SJunchao Zhang     if (threads) {
10439566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads));
104485649d77SJunchao Zhang       (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads);
104585649d77SJunchao Zhang     } else {
10462f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt)omp_get_max_threads();
10479566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads));
104885649d77SJunchao Zhang     }
1049d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys");
10509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg));
10519566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag));
1052d0609cedSBarry Smith     PetscOptionsEnd();
105385649d77SJunchao Zhang     if (flg) {
10549566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads));
105585649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
105685649d77SJunchao Zhang     }
105748a46eb9SPierre Jolivet     if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads));
105885649d77SJunchao Zhang   }
105985649d77SJunchao Zhang #endif
106085649d77SJunchao Zhang 
106185649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
106285649d77SJunchao Zhang   /*
106385649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
106485649d77SJunchao Zhang 
106585649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
106685649d77SJunchao Zhang   */
10679566063dSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL));
106885649d77SJunchao Zhang #endif
106985649d77SJunchao Zhang 
107085649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10719566063dSJacob Faibussowitsch   PetscCall(PetscFPTCreate(10000));
107285649d77SJunchao Zhang #endif
107385649d77SJunchao Zhang 
107485649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
107585649d77SJunchao Zhang   {
107685649d77SJunchao Zhang     PetscViewer viewer;
10779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg));
107885649d77SJunchao Zhang     if (flg) {
10799566063dSJacob Faibussowitsch       PetscCall(PetscProcessPlacementView(viewer));
10809566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
108185649d77SJunchao Zhang     }
108285649d77SJunchao Zhang   }
108385649d77SJunchao Zhang #endif
108485649d77SJunchao Zhang 
108585649d77SJunchao Zhang   flg = PETSC_TRUE;
10869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL));
10879566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE));
108885649d77SJunchao Zhang 
108985649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10903ba16761SJacob Faibussowitsch   PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD);
10913ba16761SJacob Faibussowitsch   PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default);
10923ba16761SJacob Faibussowitsch   PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", "");
10933ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, "");
109485649d77SJunchao Zhang #endif
109585649d77SJunchao Zhang 
109685649d77SJunchao Zhang #if defined(__VALGRIND_H)
109785649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE;
109885649d77SJunchao Zhang   #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
10999566063dSJacob 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"));
110085649d77SJunchao Zhang   #endif
110185649d77SJunchao Zhang #endif
110285649d77SJunchao Zhang   /*
110385649d77SJunchao Zhang       Set flag that we are completely initialized
110485649d77SJunchao Zhang   */
110585649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
110685649d77SJunchao Zhang 
11079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg));
11089566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPythonInitialize(NULL, NULL));
1109f1f2ae84SBarry Smith 
1110f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1111f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin());
1112f1f2ae84SBarry 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");
11133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111485649d77SJunchao Zhang }
111585649d77SJunchao Zhang 
1116e5c89e4eSSatish Balay /*@C
1117e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1118811af0c4SBarry Smith    `PetscInitialize()` calls MPI_Init() if that has yet to be called,
1119e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1120e5c89e4eSSatish Balay    your program -- usually the very first line!
1121e5c89e4eSSatish Balay 
1122811af0c4SBarry Smith    Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set
1123e5c89e4eSSatish Balay 
1124e5c89e4eSSatish Balay    Input Parameters:
1125e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1126e5c89e4eSSatish Balay .  args - the command line arguments
1127be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1128be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1129be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1130c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
11310298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1132e5c89e4eSSatish Balay 
1133811af0c4SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that
1134811af0c4SBarry Smith    communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a
1135811af0c4SBarry Smith    four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not,
1136811af0c4SBarry Smith    then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even
113705827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1138e5c89e4eSSatish Balay 
1139e5c89e4eSSatish Balay    Options Database Keys:
11407ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
11417ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1142e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11438a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
1144811af0c4SBarry Smith .  -on_error_abort - calls `abort()` when error detected (no traceback)
1145811af0c4SBarry Smith .  -on_error_mpiabort - calls `MPI_abort()` when error detected
11461af3601dSBarry Smith .  -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr
11478a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1148bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1149e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1150e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1151e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
115279dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
115379dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
1154811af0c4SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see `PetscMallocSetDebug()`
1155aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
115692f119d6SBarry 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
1157811af0c4SBarry Smith .  -malloc_view - show a list of all allocated memory during `PetscFinalize()`
115892f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1159608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
116092f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1161e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1162e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1163e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1164e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1165e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
11660841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1167e5c89e4eSSatish Balay 
1168c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1169c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1170c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1171811af0c4SBarry Smith -  -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()`
1172c5b5d8d5SVaclav Hapla 
1173c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1174c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1175c5b5d8d5SVaclav Hapla    They can be used also in option files.
1176c5b5d8d5SVaclav Hapla 
1177811af0c4SBarry Smith    See `PetscOptionsMonitorSet()` to do monitoring programmatically.
1178c5b5d8d5SVaclav Hapla 
1179e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1180a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1181811af0c4SBarry Smith +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`.
1182217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1183217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1184495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1185811af0c4SBarry Smith         hangs without running in the debugger).  See `PetscLogTraceBegin()`.
1186811af0c4SBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see `PetscLogView()`.
1187811af0c4SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`.
1188811af0c4SBarry Smith .  -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView().
11899a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1190495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1191f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1192811af0c4SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See `PetscLogDump()`.
1193811af0c4SBarry Smith .  -log [filename] - Logs basic profiline information  See `PetscLogDump()`.
1194c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
1195811af0c4SBarry Smith .  -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off
1196c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1197495fc317SBarry Smith 
11986a77f485SPierre Jolivet     Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time
1199e5c89e4eSSatish Balay 
1200ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1201ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1202ffbd1cfbSBarry 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
1203ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1204ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1205ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1206ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1207ffbd1cfbSBarry Smith 
1208e5c89e4eSSatish Balay    Environmental Variables:
1209811af0c4SBarry Smith +   `PETSC_TMP` - alternative tmp directory
1210811af0c4SBarry Smith .   `PETSC_SHARED_TMP` - tmp is shared by all processes
1211811af0c4SBarry Smith .   `PETSC_NOT_SHARED_TMP` - each process has its own private tmp
1212811af0c4SBarry Smith .   `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs
1213811af0c4SBarry Smith .   `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1214811af0c4SBarry Smith .   `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer
1215811af0c4SBarry Smith -   `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to
1216e5c89e4eSSatish Balay 
1217e5c89e4eSSatish Balay    Level: beginner
1218e5c89e4eSSatish Balay 
1219811af0c4SBarry Smith    Note:
1220811af0c4SBarry Smith    If for some reason you must call `MPI_Init()` separately, call
1221811af0c4SBarry Smith    it before `PetscInitialize()`.
1222e5c89e4eSSatish Balay 
1223811af0c4SBarry Smith    Fortran Notes:
1224811af0c4SBarry Smith    In Fortran this routine can be called with
1225811af0c4SBarry Smith .vb
1226811af0c4SBarry Smith        call PetscInitialize(ierr)
1227811af0c4SBarry Smith        call PetscInitialize(file,ierr) or
1228811af0c4SBarry Smith        call PetscInitialize(file,help,ierr)
1229811af0c4SBarry Smith .ve
1230e5c89e4eSSatish Balay 
1231811af0c4SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after
1232811af0c4SBarry Smith    calling `PetscInitialize()`.
1233e5c89e4eSSatish Balay 
12345dedd997SBarry Smith    Options Database Key for Developers:
12355dedd997SBarry Smith .  -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form:
12365dedd997SBarry Smith     "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py)
12375dedd997SBarry Smith 
1238db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()`
1239e5c89e4eSSatish Balay @*/
1240d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[])
1241d71ae5a4SJacob Faibussowitsch {
124285649d77SJunchao Zhang   PetscMPIInt flag;
124321ef0414SBarry Smith   const char *prog = "Unknown Name", *mpienv;
1244e5c89e4eSSatish Balay 
124527104ee2SJacob Faibussowitsch   PetscFunctionBegin;
12463ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
12479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Initialized(&flag));
1248e5c89e4eSSatish Balay   if (!flag) {
124939a651e2SJacob 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");
12509566063dSJacob Faibussowitsch     PetscCall(PetscPreMPIInit_Private());
12515e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12525e765c61SJed Brown     {
125339a651e2SJacob Faibussowitsch       PetscMPIInt PETSC_UNUSED provided;
12549566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED, &provided));
12555e765c61SJed Brown     }
12565e765c61SJed Brown #else
12579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Init(argc, args));
12585e765c61SJed Brown #endif
125921ef0414SBarry Smith     if (PetscDefined(HAVE_MPIUNI)) {
126021ef0414SBarry Smith       mpienv = getenv("PMI_SIZE");
126121ef0414SBarry Smith       if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE");
126221ef0414SBarry Smith       if (mpienv) {
126321ef0414SBarry Smith         PetscInt isize;
126421ef0414SBarry Smith         PetscCall(PetscOptionsStringToInt(mpienv, &isize));
126521ef0414SBarry 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");
126621ef0414SBarry 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");
126721ef0414SBarry Smith       }
126821ef0414SBarry Smith     }
1269e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1270e5c89e4eSSatish Balay   }
12714dfee713SSatish Balay 
127285649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1273e5c89e4eSSatish Balay   if (argc && args) {
1274e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1275e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1276e5c89e4eSSatish Balay   }
1277811af0c4SBarry Smith   PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0));
12783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1279e5c89e4eSSatish Balay }
1280e5c89e4eSSatish Balay 
1281a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
128295c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1283ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsCounts;
1284ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsMaxCounts;
128595c0884eSLisandro Dalcin PETSC_INTERN PetscBool    PetscObjectsLog;
12864097062eSBarry Smith #endif
1287e5c89e4eSSatish Balay 
1288008a6e76SBarry Smith /*
1289008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1290008a6e76SBarry Smith */
1291d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void)
1292d71ae5a4SJacob Faibussowitsch {
1293008a6e76SBarry Smith   PetscFunctionBegin;
1294613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
12959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128));
129674c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
12979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128));
1298008a6e76SBarry Smith   #endif
129974c01ffaSSatish Balay #endif
13009e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
13019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FP16));
1302008a6e76SBarry Smith #endif
1303008a6e76SBarry Smith 
1304de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
13059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_SUM));
1306613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MAX));
1307613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MIN));
13089e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
13099e517322SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128));
1310008a6e76SBarry Smith #endif
1311008a6e76SBarry Smith 
13129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR));
13139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT));
13149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT));
131540df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
13169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2INT));
1317008a6e76SBarry Smith #endif
13189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPI_4INT));
13199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_4INT));
13209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
132162e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp));
13223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1323008a6e76SBarry Smith }
1324008a6e76SBarry Smith 
1325a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
1326a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1327a4af0ceeSJacob Faibussowitsch #endif
1328a4af0ceeSJacob Faibussowitsch 
1329e5c89e4eSSatish Balay /*@C
1330e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1331811af0c4SBarry Smith    of the program. `MPI_Finalize()` is called only if the user had not
1332811af0c4SBarry Smith    called `MPI_Init()` before calling `PetscInitialize()`.
1333e5c89e4eSSatish Balay 
1334811af0c4SBarry Smith    Collective on `PETSC_COMM_WORLD`
1335e5c89e4eSSatish Balay 
1336e5c89e4eSSatish Balay    Options Database Keys:
1337811af0c4SBarry Smith +  -options_view - Calls `PetscOptionsView()`
1338e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
13397eb1d149SBarry Smith .  -objects_dump [all] - Prints list of objects allocated by the user that have not been freed, the option all cause all outstanding objects to be listed
1340e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1341811af0c4SBarry Smith .  -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed
1342e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
134392f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1344e5c89e4eSSatish Balay 
1345e5c89e4eSSatish Balay    Level: beginner
1346e5c89e4eSSatish Balay 
1347e5c89e4eSSatish Balay    Note:
1348811af0c4SBarry Smith    See `PetscInitialize()` for other runtime options.
1349e5c89e4eSSatish Balay 
1350db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()`
1351e5c89e4eSSatish Balay @*/
1352d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void)
1353d71ae5a4SJacob Faibussowitsch {
13544bb5149bSJed Brown   PetscMPIInt rank;
1355a8d2bbe5SBarry Smith   PetscInt    nopt;
13562bf49c77SBarry Smith   PetscBool   flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE;
1357dff31646SBarry Smith   PetscBool   flg;
135810463e74SBarry Smith #if defined(PETSC_USE_LOG)
135910463e74SBarry Smith   char mname[PETSC_MAX_PATH_LEN];
136010463e74SBarry Smith #endif
1361e5c89e4eSSatish Balay 
13623cbbc5ffSBarry Smith   PetscFunctionBegin;
136339a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()");
13649566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PetscFinalize() called\n"));
1365b022a5c1SBarry Smith 
1366f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1367f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd());
1368f1f2ae84SBarry Smith 
136962e5d2d2SJDBetteridge   /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */
137062e5d2d2SJDBetteridge   {
137162e5d2d2SJDBetteridge     union
137262e5d2d2SJDBetteridge     {
137362e5d2d2SJDBetteridge       MPI_Comm comm;
137462e5d2d2SJDBetteridge       void    *ptr;
137562e5d2d2SJDBetteridge     } ucomm;
137662e5d2d2SJDBetteridge     PetscMPIInt flg;
137762e5d2d2SJDBetteridge     void       *tmp;
137862e5d2d2SJDBetteridge 
137962e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
138062e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
138162e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF));
138262e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
138362e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
138462e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD));
138562e5d2d2SJDBetteridge   }
138662e5d2d2SJDBetteridge 
13879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1388a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13893ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE);
13903ba16761SJacob Faibussowitsch   PetscCallExternal(adios_finalize, rank);
1391a56f64adSBarry Smith #endif
13929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg));
1393dff31646SBarry Smith   if (flg) {
13941f817a21SBarry Smith     char *cits, filename[PETSC_MAX_PATH_LEN];
13951f817a21SBarry Smith     FILE *fd = PETSC_STDOUT;
13961f817a21SBarry Smith 
13979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL));
139848a46eb9SPierre Jolivet     if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd));
13999566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits));
1400dff31646SBarry Smith     cits[0] = 0;
14019566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits));
14029566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n"));
14039566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14049566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits));
14059566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14069566063dSJacob Faibussowitsch     PetscCall(PetscFClose(PETSC_COMM_WORLD, fd));
14079566063dSJacob Faibussowitsch     PetscCall(PetscFree(cits));
1408dff31646SBarry Smith   }
14099566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&PetscCitationsList));
1410dff31646SBarry Smith 
1411c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
141204102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
141304102261SBarry Smith   {
141404102261SBarry Smith     PetscInt nmax = 2;
141504102261SBarry Smith     char   **buffs;
14169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
14179566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-textbelt", buffs, &nmax, &flg1));
141804102261SBarry Smith     if (flg1) {
141928b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-textbelt requires either the phone number or number,\"message\"");
142004102261SBarry Smith       if (nmax == 1) {
1421c6a7a370SJeremy L Thompson         size_t len = 128;
1422c6a7a370SJeremy L Thompson         PetscCall(PetscMalloc1(len, &buffs[1]));
14239566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
1424c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(buffs[1], " has completed", len));
142504102261SBarry Smith       }
14269566063dSJacob Faibussowitsch       PetscCall(PetscTextBelt(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
14279566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
14289566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
142904102261SBarry Smith     }
14309566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
143104102261SBarry Smith   }
1432763ec7b1SBarry Smith   {
1433763ec7b1SBarry Smith     PetscInt nmax = 2;
1434763ec7b1SBarry Smith     char   **buffs;
14359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
14369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-tellmycell", buffs, &nmax, &flg1));
1437763ec7b1SBarry Smith     if (flg1) {
143828b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-tellmycell requires either the phone number or number,\"message\"");
1439763ec7b1SBarry Smith       if (nmax == 1) {
1440c6a7a370SJeremy L Thompson         size_t len = 128;
1441c6a7a370SJeremy L Thompson         PetscCall(PetscMalloc1(len, &buffs[1]));
14429566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
1443c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(buffs[1], " has completed", len));
1444763ec7b1SBarry Smith       }
14459566063dSJacob Faibussowitsch       PetscCall(PetscTellMyCell(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
14469566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
14479566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
1448763ec7b1SBarry Smith     }
14499566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
1450763ec7b1SBarry Smith   }
145104102261SBarry Smith #endif
145204102261SBarry Smith 
14532d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14549566063dSJacob Faibussowitsch   PetscCall(PetscFPTDestroy());
14552d53ad75SBarry Smith #endif
14562d53ad75SBarry Smith 
1457e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1458dff31646SBarry Smith   flg = PETSC_FALSE;
14599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL));
14601baa6e33SBarry Smith   if (flg) PetscCall(PetscOptionsSAWsDestroy());
1461d5649816SBarry Smith #endif
1462d5649816SBarry Smith 
1463681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1464681455b2SBarry Smith   flg1 = PETSC_FALSE;
14659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL));
1466681455b2SBarry Smith   if (flg1) {
1467681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
14689566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL));
1469681455b2SBarry Smith   }
1470681455b2SBarry Smith #endif
1471681455b2SBarry Smith 
147267584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_info", &flg2, NULL));
1474e5c89e4eSSatish Balay   if (!flg2) {
147590d69ab7SBarry Smith     flg2 = PETSC_FALSE;
14769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL));
1477e5c89e4eSSatish Balay   }
147848a46eb9SPierre Jolivet   if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n"));
147967584ceeSBarry Smith #endif
1480e5c89e4eSSatish Balay 
1481e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
148290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL));
1484e5c89e4eSSatish Balay   if (flg1) {
1485e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
14869566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD));
14879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops));
1488e5c89e4eSSatish Balay   }
1489e5c89e4eSSatish Balay #endif
1490e5c89e4eSSatish Balay 
1491e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1492e5c89e4eSSatish Balay   #if defined(PETSC_HAVE_MPE)
1493e5c89e4eSSatish Balay   mname[0] = 0;
14949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1));
1495e5c89e4eSSatish Balay   if (flg1) {
14969566063dSJacob Faibussowitsch     if (mname[0]) PetscCall(PetscLogMPEDump(mname));
14979566063dSJacob Faibussowitsch     else PetscCall(PetscLogMPEDump(0));
1498e5c89e4eSSatish Balay   }
1499e5c89e4eSSatish Balay   #endif
1500356e5837SBarry Smith #endif
1501a297a907SKarl Rupp 
1502dd710f27SBarry Smith   /*
1503dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1504dd710f27SBarry Smith   */
15059566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1506dd710f27SBarry Smith 
1507356e5837SBarry Smith #if defined(PETSC_USE_LOG)
15089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE));
15099566063dSJacob Faibussowitsch   PetscCall(PetscLogViewFromOptions());
15109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPopGetViewerOff());
151187c3beb6SLisandro Dalcin 
1512356e5837SBarry Smith   mname[0] = 0;
15139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_summary", mname, sizeof(mname), &flg1));
1514e5c89e4eSSatish Balay   if (flg1) {
151591eabc43SBarry Smith     PetscViewer viewer;
15169566063dSJacob Faibussowitsch     PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD, "\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n"));
151791eabc43SBarry Smith     if (mname[0]) {
15189566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, mname, &viewer));
15199566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
15209566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
152133f85c2fSBarry Smith     } else {
152233f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
15239566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
15249566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
15259566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
152633f85c2fSBarry Smith     }
1527e5c89e4eSSatish Balay   }
1528a297a907SKarl Rupp 
1529dd710f27SBarry Smith   /*
1530dd710f27SBarry Smith      Free any objects created by the last block of code.
1531dd710f27SBarry Smith   */
15329566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1533dd710f27SBarry Smith 
1534dd710f27SBarry Smith   mname[0] = 0;
15359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1));
15369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2));
15379566063dSJacob Faibussowitsch   if (flg1 || flg2) PetscCall(PetscLogDump(mname));
1538e5c89e4eSSatish Balay #endif
153910463e74SBarry Smith 
154090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL));
15429566063dSJacob Faibussowitsch   if (!flg1) PetscCall(PetscPopSignalHandler());
154390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL));
15451baa6e33SBarry Smith   if (flg1) PetscCall(PetscMPIDump(stdout));
154690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
154790d69ab7SBarry Smith   flg2 = PETSC_FALSE;
15488bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
15499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1));
15509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
15519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL));
1552e4c476e2SSatish Balay 
1553e5c89e4eSSatish Balay   if (flg2) {
1554be56827dSJed Brown     PetscViewer viewer;
15559566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15569566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsView(NULL, viewer));
15589566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1559e5c89e4eSSatish Balay   }
1560e5c89e4eSSatish Balay 
1561e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
15629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1));
15639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1));
1564e5c89e4eSSatish Balay 
156533fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
15669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1));
15679245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1568e5c89e4eSSatish Balay   if (flg3) {
15693de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1570be56827dSJed Brown       PetscViewer viewer;
15719566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15729566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15739566063dSJacob Faibussowitsch       PetscCall(PetscOptionsView(NULL, viewer));
15749566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
1575e5c89e4eSSatish Balay     }
15769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsAllUsed(NULL, &nopt));
15773de2bfdfSBarry Smith     if (nopt) {
15789566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n"));
15799566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n"));
15803de2bfdfSBarry Smith       if (nopt == 1) {
15819566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n"));
1582e5c89e4eSSatish Balay       } else {
15839566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt));
1584e5c89e4eSSatish Balay       }
15853de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15869566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n"));
1587df12ba86SBarry Smith     }
15889566063dSJacob Faibussowitsch     PetscCall(PetscOptionsLeft(NULL));
1589e5c89e4eSSatish Balay   }
1590e5c89e4eSSatish Balay 
1591e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1592d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15939566063dSJacob Faibussowitsch     PetscCall(PetscStackSAWsViewOff());
1594792fecdfSBarry Smith     PetscCallSAWs(SAWs_Finalize, ());
1595d45a07a7SBarry Smith   }
1596ec957eceSBarry Smith #endif
1597ec957eceSBarry Smith 
15984097062eSBarry Smith #if defined(PETSC_USE_LOG)
159910463e74SBarry Smith   /*
1600dbc8283eSBarry Smith        List all objects the user may have forgot to free
16012eff7a51SBarry Smith   */
160205df10baSBarry Smith   if (PetscObjectsLog) {
16039566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
1604a64a8e02SBarry Smith     if (flg1) {
1605a64a8e02SBarry Smith       MPI_Comm local_comm;
16067eb1d149SBarry Smith       char     string[64];
1607a64a8e02SBarry Smith 
16089566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL));
1609afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16109566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16119566063dSJacob Faibussowitsch       PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
16129566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16139566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
16140a1571b3SBarry Smith     }
161505df10baSBarry Smith   }
16164097062eSBarry Smith #endif
16174097062eSBarry Smith 
16184097062eSBarry Smith #if defined(PETSC_USE_LOG)
1619dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1620dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
16219566063dSJacob Faibussowitsch   PetscCall(PetscFree(PetscObjects));
16224097062eSBarry Smith #endif
16232eff7a51SBarry Smith 
162433f85c2fSBarry Smith   /*
162533f85c2fSBarry Smith      Destroy any packages that registered a finalize
162633f85c2fSBarry Smith   */
16279566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalizeAll());
162833f85c2fSBarry Smith 
1629101409b8SToby Isaac #if defined(PETSC_USE_LOG)
16309566063dSJacob Faibussowitsch   PetscCall(PetscLogFinalize());
1631101409b8SToby Isaac #endif
1632101409b8SToby Isaac 
163333f85c2fSBarry Smith   /*
163448dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
163548dd1dffSBarry Smith   */
16362e956fe4SStefano Zampini   if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll());
163737e93019SBarry Smith 
16384028d114SSatish Balay   if (petsc_history) {
16399566063dSJacob Faibussowitsch     PetscCall(PetscCloseHistoryFile(&petsc_history));
164002c9f0b5SLisandro Dalcin     petsc_history = NULL;
1641e5c89e4eSSatish Balay   }
16429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
16439566063dSJacob Faibussowitsch   PetscCall(PetscInfoDestroy());
1644e5c89e4eSSatish Balay 
164567584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
164692f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1647e5c89e4eSSatish Balay     char  fname[PETSC_MAX_PATH_LEN];
164892f119d6SBarry Smith     char  sname[PETSC_MAX_PATH_LEN];
1649e5c89e4eSSatish Balay     FILE *fd;
1650ed9cf6e9SBarry Smith     int   err;
1651e5c89e4eSSatish Balay 
1652dc92acbaSJed Brown     flg2 = PETSC_FALSE;
165392f119d6SBarry Smith     flg3 = PETSC_FALSE;
16549566063dSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL));
16559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL));
165692f119d6SBarry Smith     fname[0] = 0;
16579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1));
1658e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
16593ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
16609371c9d4SSatish Balay       fd = fopen(sname, "w");
16619371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16629566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(fd));
1663ed9cf6e9SBarry Smith       err = fclose(fd);
166428b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
166592f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1666e5c89e4eSSatish Balay       MPI_Comm local_comm;
1667e5c89e4eSSatish Balay 
1668afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16699566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16709566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(stdout));
16719566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16729566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1673e5c89e4eSSatish Balay     }
1674e5c89e4eSSatish Balay     fname[0] = 0;
16759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1));
1676e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
16773ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
16789371c9d4SSatish Balay       fd = fopen(sname, "w");
16799371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16809566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(fd));
1681ed9cf6e9SBarry Smith       err = fclose(fd);
168228b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
168392f119d6SBarry Smith     } else if (flg1) {
168492f119d6SBarry Smith       MPI_Comm local_comm;
168592f119d6SBarry Smith 
1686afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16879566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16889566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(stdout));
16899566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16909566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1691e5c89e4eSSatish Balay     }
1692e5c89e4eSSatish Balay   }
169367584ceeSBarry Smith #endif
169420e2c332SMatthew G. Knepley 
16955486ca60SMatthew G. Knepley   /*
16965486ca60SMatthew G. Knepley      Close any open dynamic libraries
16975486ca60SMatthew G. Knepley   */
16989566063dSJacob Faibussowitsch   PetscCall(PetscFinalize_DynamicLibraries());
16995486ca60SMatthew G. Knepley 
1700e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
17019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroyDefault());
1702e5c89e4eSSatish Balay 
1703e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
170402c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1705e5c89e4eSSatish Balay 
1706c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1707c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
17089566063dSJacob Faibussowitsch     PetscCall(PetscKokkosFinalize_Private());
1709c2b86a48SJunchao Zhang     PetscBeganKokkos       = PETSC_FALSE;
171045639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1711c2b86a48SJunchao Zhang   }
1712c2b86a48SJunchao Zhang #endif
1713c2b86a48SJunchao Zhang 
171471438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
171571438e86SJunchao Zhang   if (PetscBeganNvshmem) {
17169566063dSJacob Faibussowitsch     PetscCall(PetscNvshmemFinalize());
171771438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
171871438e86SJunchao Zhang   }
171971438e86SJunchao Zhang #endif
1720a0e72f99SJunchao Zhang 
17219566063dSJacob Faibussowitsch   PetscCall(PetscFreeMPIResources());
1722e5c89e4eSSatish Balay 
1723dbc8283eSBarry Smith   /*
1724efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1725efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1726efb80d3cSBarry Smith 
1727efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1728efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1729dbc8283eSBarry Smith  */
1730b770b1f6SSatish Balay   {
1731dbc8283eSBarry Smith     PetscCommCounter *counter;
1732dbc8283eSBarry Smith     PetscMPIInt       flg;
1733dbc8283eSBarry Smith     MPI_Comm          icomm;
17349371c9d4SSatish Balay     union
17359371c9d4SSatish Balay     {
17369371c9d4SSatish Balay       MPI_Comm comm;
17379371c9d4SSatish Balay       void    *ptr;
17389371c9d4SSatish Balay     } ucomm;
17399566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
1740dbc8283eSBarry Smith     if (flg) {
1741265f3f35SJed Brown       icomm = ucomm.comm;
17429566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
174328b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1744dbc8283eSBarry Smith 
17459566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval));
17469566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17479566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1748dbc8283eSBarry Smith     }
17499566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
1750dbc8283eSBarry Smith     if (flg) {
1751265f3f35SJed Brown       icomm = ucomm.comm;
17529566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
175328b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1754dbc8283eSBarry Smith 
17559566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval));
17569566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17579566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1758dbc8283eSBarry Smith     }
1759b770b1f6SSatish Balay   }
1760dbc8283eSBarry Smith 
17619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
17629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
17639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
17649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
176562e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval));
176662e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval));
1767480cf27aSJed Brown 
17685ea2b939SDuncan Campbell   // Free keyvals which may be silently created by some routines
17695ea2b939SDuncan Campbell   if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval));
17705ea2b939SDuncan Campbell   if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval));
17715ea2b939SDuncan Campbell 
17729566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
17739566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
17749566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
17759566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock));
1776ef19f930SBarry Smith 
1777e5c89e4eSSatish Balay   if (PetscBeganMPI) {
177899b1327fSBarry Smith     PetscMPIInt flag;
17799566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalized(&flag));
178039a651e2SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
178139a651e2SJacob Faibussowitsch     /* wait until the very last moment to disable error handling */
178239a651e2SJacob Faibussowitsch     PetscErrorHandlingInitialized = PETSC_FALSE;
17839566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalize());
178439a651e2SJacob Faibussowitsch   } else PetscErrorHandlingInitialized = PETSC_FALSE;
178539a651e2SJacob Faibussowitsch 
1786e5c89e4eSSatish Balay   /*
1787e5c89e4eSSatish Balay 
1788e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1789e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1790e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1791e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1792e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17930e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1794e5c89e4eSSatish Balay    memory was not freed.
1795e5c89e4eSSatish Balay 
1796e5c89e4eSSatish Balay */
17979566063dSJacob Faibussowitsch   PetscCall(PetscMallocClear());
17989566063dSJacob Faibussowitsch   PetscCall(PetscStackReset());
1799a297a907SKarl Rupp 
1800e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1801e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
18027ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE)
180356883f60SBarry Smith   /*
180456883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
180556883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
180656883f60SBarry Smith    */
180756883f60SBarry Smith   __gcov_flush();
180856883f60SBarry Smith #endif
18091724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
18101724198aSStefano Zampini   PetscStackClearTop;
18113ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1812e5c89e4eSSatish Balay }
1813e5c89e4eSSatish Balay 
181443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
1815d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b)
1816d71ae5a4SJacob Faibussowitsch {
181743db4dbbSBarry Smith   if (*a == *b) return 1;
181843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
181943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
182043db4dbbSBarry Smith   return 0;
182143db4dbbSBarry Smith }
1822a70650f6SBarry Smith #endif
182343db4dbbSBarry Smith 
182443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
1825d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b)
1826d71ae5a4SJacob Faibussowitsch {
182743db4dbbSBarry Smith   if (*a == *b) return 1;
182843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
182943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
183043db4dbbSBarry Smith   return 0;
183143db4dbbSBarry Smith }
183243db4dbbSBarry Smith #endif
1833