xref: /petsc/src/sys/objects/pinit.c (revision 5ea2b939da046e2ffbe5d369773bb579d11fd712) !
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 
35e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
36e5c89e4eSSatish Balay 
3795c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
38e5c89e4eSSatish Balay 
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm, int);
4295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm, int);
4395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE **);
440069ddf5SShri Abhyankar 
456de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
46e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4727104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
486de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
496de5d289SStefano Zampini #else
506de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
516de5d289SStefano Zampini #endif
52e5c89e4eSSatish Balay 
53480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval      = MPI_KEYVAL_INVALID;
54480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval    = MPI_KEYVAL_INVALID;
55480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval    = MPI_KEYVAL_INVALID;
565f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval      = MPI_KEYVAL_INVALID;
5762e5d2d2SJDBetteridge PetscMPIInt Petsc_CreationIdx_keyval  = MPI_KEYVAL_INVALID;
5862e5d2d2SJDBetteridge PetscMPIInt Petsc_Garbage_HMap_keyval = MPI_KEYVAL_INVALID;
59480cf27aSJed Brown 
60*5ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedWD_keyval  = MPI_KEYVAL_INVALID;
61*5ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedTmp_keyval = MPI_KEYVAL_INVALID;
62*5ea2b939SDuncan Campbell 
63e5c89e4eSSatish Balay /*
64e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
65e5c89e4eSSatish Balay */
6602c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL};
6702c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL};
68e5c89e4eSSatish Balay 
69ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
70ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
710f8e0872SSatish Balay 
72a2f94806SJed Brown PetscInt PetscHotRegionDepth;
73a2f94806SJed Brown 
746edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
756edef35eSSatish Balay 
76b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
77b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
78b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
79b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
80b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
81b22622e2STadeu Manoel #endif
82b22622e2STadeu Manoel 
83e5c89e4eSSatish Balay /*
84945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8572a42c3cSBarry Smith 
8672a42c3cSBarry Smith    Collective
8772a42c3cSBarry Smith 
8872a42c3cSBarry Smith    Level: advanced
8972a42c3cSBarry Smith 
9095452b02SPatrick Sanan     Notes:
91a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
920f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
93a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
940f11a792SBarry Smith 
95a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
961ea5a559SBarry Smith 
97db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()`
980f11a792SBarry Smith */
99d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoPointers(int argc, char **args, const char *filename, const char *help)
100d71ae5a4SJacob Faibussowitsch {
10172a42c3cSBarry Smith   int    myargc = argc;
10272a42c3cSBarry Smith   char **myargs = args;
10372a42c3cSBarry Smith 
10472a42c3cSBarry Smith   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&myargc, &myargs, filename, help));
1069566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
107df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10972a42c3cSBarry Smith }
11072a42c3cSBarry Smith 
111f0865b08SBarry Smith /*
112a56f64adSBarry Smith       Used by Julia interface to get communicator
113f0865b08SBarry Smith */
114d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
115d71ae5a4SJacob Faibussowitsch {
116f0865b08SBarry Smith   PetscFunctionBegin;
11727104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidPointer(comm, 1);
118f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120f0865b08SBarry Smith }
121f0865b08SBarry Smith 
122e5c89e4eSSatish Balay /*@C
123811af0c4SBarry Smith       PetscInitializeNoArguments - Calls `PetscInitialize()` from C/C++ without
124e5c89e4eSSatish Balay         the command line arguments.
125e5c89e4eSSatish Balay 
126e5c89e4eSSatish Balay    Collective
127e5c89e4eSSatish Balay 
128e5c89e4eSSatish Balay    Level: advanced
129e5c89e4eSSatish Balay 
130db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`
131e5c89e4eSSatish Balay @*/
132d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoArguments(void)
133d71ae5a4SJacob Faibussowitsch {
134e5c89e4eSSatish Balay   int    argc = 0;
13502c9f0b5SLisandro Dalcin   char **args = NULL;
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay   PetscFunctionBegin;
1389566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, NULL, NULL));
1393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
140e5c89e4eSSatish Balay }
141e5c89e4eSSatish Balay 
142e5c89e4eSSatish Balay /*@
143e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
144e5c89e4eSSatish Balay 
14593b6d2d1SJed Brown    Level: beginner
146e5c89e4eSSatish Balay 
147db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
148e5c89e4eSSatish Balay @*/
149d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialized(PetscBool *isInitialized)
150d71ae5a4SJacob Faibussowitsch {
15127104ee2SJacob Faibussowitsch   PetscFunctionBegin;
15227104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized, 1);
153e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
1543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155e5c89e4eSSatish Balay }
156e5c89e4eSSatish Balay 
157e5c89e4eSSatish Balay /*@
158811af0c4SBarry Smith       PetscFinalized - Determine whether `PetscFinalize()` has been called yet
159e5c89e4eSSatish Balay 
160e5c89e4eSSatish Balay    Level: developer
161e5c89e4eSSatish Balay 
162db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()`
163e5c89e4eSSatish Balay @*/
164d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalized(PetscBool *isFinalized)
165d71ae5a4SJacob Faibussowitsch {
16627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
16727104ee2SJacob Faibussowitsch   if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized, 1);
168e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170e5c89e4eSSatish Balay }
171e5c89e4eSSatish Balay 
17257171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]);
173e5c89e4eSSatish Balay 
174e5c89e4eSSatish Balay /*
175e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
176e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
177e5c89e4eSSatish Balay */
178367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP               = 0;
17962e5d2d2SJDBetteridge MPI_Op Petsc_Garbage_SetIntersectOp = 0;
180e5c89e4eSSatish Balay 
181d71ae5a4SJacob Faibussowitsch PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, int *cnt, MPI_Datatype *datatype)
182d71ae5a4SJacob Faibussowitsch {
183e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out, i, count = *cnt;
184e5c89e4eSSatish Balay 
185e5c89e4eSSatish Balay   PetscFunctionBegin;
186e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
1873ba16761SJacob Faibussowitsch     PetscErrorCode ierr = (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
1883ba16761SJacob Faibussowitsch     (void)ierr;
18941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
190e5c89e4eSSatish Balay   }
191e5c89e4eSSatish Balay 
192e5c89e4eSSatish Balay   for (i = 0; i < count; i++) {
193e5c89e4eSSatish Balay     xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]);
194e5c89e4eSSatish Balay     xout[2 * i + 1] += xin[2 * i + 1];
195e5c89e4eSSatish Balay   }
196812af9f3SBarry Smith   PetscFunctionReturnVoid();
197e5c89e4eSSatish Balay }
198e5c89e4eSSatish Balay 
199e5c89e4eSSatish Balay /*
200e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
201e5c89e4eSSatish Balay sum of the second entry.
202b693b147SBarry Smith 
20376f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
204367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
205b693b147SBarry Smith there would be no place to store the both needed results.
206e5c89e4eSSatish Balay */
207d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt sizes[], PetscInt *max, PetscInt *sum)
208d71ae5a4SJacob Faibussowitsch {
209e5c89e4eSSatish Balay   PetscFunctionBegin;
210d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
211d6e4c47cSJed Brown   {
2129371c9d4SSatish Balay     struct {
2139371c9d4SSatish Balay       PetscInt max, sum;
2149371c9d4SSatish Balay     } work;
2159566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce_scatter_block((void *)sizes, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm));
216d6e4c47cSJed Brown     *max = work.max;
217d6e4c47cSJed Brown     *sum = work.sum;
218d6e4c47cSJed Brown   }
219d6e4c47cSJed Brown #else
220d6e4c47cSJed Brown   {
221d6e4c47cSJed Brown     PetscMPIInt size, rank;
2229371c9d4SSatish Balay     struct {
2239371c9d4SSatish Balay       PetscInt max, sum;
2249371c9d4SSatish Balay     } *work;
2259566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
2269566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
2279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &work));
2281c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce((void *)sizes, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm));
2296ac3741eSJed Brown     *max = work[rank].max;
2306ac3741eSJed Brown     *sum = work[rank].sum;
2319566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
232d6e4c47cSJed Brown   }
233d6e4c47cSJed Brown #endif
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
235e5c89e4eSSatish Balay }
236e5c89e4eSSatish Balay 
237e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
238e5c89e4eSSatish Balay 
2399e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
240613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
241613bf2b2SPierre Jolivet     #include <quadmath.h>
242613bf2b2SPierre Jolivet   #endif
2439e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0;
244de272c7aSSatish Balay   #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
24506a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
246613bf2b2SPierre Jolivet   #endif
247e5c89e4eSSatish Balay 
248d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
249d71ae5a4SJacob Faibussowitsch {
250e5c89e4eSSatish Balay   PetscInt i, count = *cnt;
251e5c89e4eSSatish Balay 
252e5c89e4eSSatish Balay   PetscFunctionBegin;
2537c2de775SJed Brown   if (*datatype == MPIU_REAL) {
254e2e03761SBarry Smith     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
255a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2567c2de775SJed Brown   }
2577c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
258c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2597c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
260a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] += xin[i];
2617c2de775SJed Brown   }
2627c2de775SJed Brown   #endif
263613bf2b2SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FLOAT128)
264613bf2b2SPierre Jolivet   else if (*datatype == MPIU___FLOAT128) {
265613bf2b2SPierre Jolivet     __float128 *xin = (__float128 *)in, *xout = (__float128 *)out;
266613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2679371c9d4SSatish Balay   } else if (*datatype == MPIU___COMPLEX128) {
268613bf2b2SPierre Jolivet     __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out;
269613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
270613bf2b2SPierre Jolivet   }
271613bf2b2SPierre Jolivet   #endif
2729e517322SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FP16)
2739e517322SPierre Jolivet   else if (*datatype == MPIU___FP16) {
2749e517322SPierre Jolivet     __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out;
2759e517322SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2769e517322SPierre Jolivet   }
2779e517322SPierre Jolivet   #endif
2787c2de775SJed Brown   else {
2799e517322SPierre Jolivet   #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16)
2803ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SElF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
2819e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FP16)
2823ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types"));
2839e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FLOAT128)
2843ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types"));
2859e517322SPierre Jolivet   #else
2863ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types"));
287613bf2b2SPierre Jolivet   #endif
28841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
289e2e03761SBarry Smith   }
290812af9f3SBarry Smith   PetscFunctionReturnVoid();
291e5c89e4eSSatish Balay }
292e5c89e4eSSatish Balay #endif
293e5c89e4eSSatish Balay 
294570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
295d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
296d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
297d9822059SBarry Smith 
298d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
299d71ae5a4SJacob Faibussowitsch {
300d9822059SBarry Smith   PetscInt i, count = *cnt;
301d9822059SBarry Smith 
302d9822059SBarry Smith   PetscFunctionBegin;
3037c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3048c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
305a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]);
3067c2de775SJed Brown   }
3077c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3087c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3097c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
310ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3117c2de775SJed Brown   }
3127c2de775SJed Brown   #endif
3137c2de775SJed Brown   else {
3143ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
31541e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3168c764dc5SJose Roman   }
317d9822059SBarry Smith   PetscFunctionReturnVoid();
318d9822059SBarry Smith }
319d9822059SBarry Smith 
320d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
321d71ae5a4SJacob Faibussowitsch {
322d9822059SBarry Smith   PetscInt i, count = *cnt;
323d9822059SBarry Smith 
324d9822059SBarry Smith   PetscFunctionBegin;
3257c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3268c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
327a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]);
3287c2de775SJed Brown   }
3297c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3307c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3317c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
332ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3337c2de775SJed Brown   }
3347c2de775SJed Brown   #endif
3357c2de775SJed Brown   else {
3363ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"));
33741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3388c764dc5SJose Roman   }
339d9822059SBarry Smith   PetscFunctionReturnVoid();
340d9822059SBarry Smith }
341d9822059SBarry Smith #endif
342d9822059SBarry Smith 
343480cf27aSJed Brown /*
344480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
345480cf27aSJed Brown 
346ff0e51ddSBarry 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.
347480cf27aSJed Brown 
34812801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
349480cf27aSJed Brown 
350480cf27aSJed Brown */
351d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state)
352d71ae5a4SJacob Faibussowitsch {
35305342407SJunchao Zhang   PetscCommCounter      *counter = (PetscCommCounter *)count_val;
35457f21012SBarry Smith   struct PetscCommStash *comms   = counter->comms, *pcomm;
355480cf27aSJed Brown 
356480cf27aSJed Brown   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm));
3589566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter->iflags));
35957f21012SBarry Smith   while (comms) {
3609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&comms->comm));
36157f21012SBarry Smith     pcomm = comms;
36257f21012SBarry Smith     comms = comms->next;
3639566063dSJacob Faibussowitsch     PetscCall(PetscFree(pcomm));
36457f21012SBarry Smith   }
3659566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter));
366480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
367480cf27aSJed Brown }
368480cf27aSJed Brown 
369480cf27aSJed Brown /*
37047435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
371da3039f7SJed Brown   calls MPI_Comm_free().
372da3039f7SJed Brown 
373da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
374480cf27aSJed Brown 
375ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
376480cf27aSJed Brown 
37712801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
378480cf27aSJed Brown 
379480cf27aSJed Brown */
380d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
381d71ae5a4SJacob Faibussowitsch {
3829371c9d4SSatish Balay   union
383480cf27aSJed Brown   {
3849371c9d4SSatish Balay     MPI_Comm comm;
3859371c9d4SSatish Balay     void    *ptr;
3869371c9d4SSatish Balay   } icomm;
387480cf27aSJed Brown 
388480cf27aSJed Brown   PetscFunctionBegin;
38912801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval");
390ec4fadc2SJed Brown   icomm.ptr = attr_val;
39176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
39233779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
39333779a13SJunchao Zhang     PetscMPIInt flg;
3949371c9d4SSatish Balay     union
3959371c9d4SSatish Balay     {
3969371c9d4SSatish Balay       MPI_Comm comm;
3979371c9d4SSatish Balay       void    *ptr;
3989371c9d4SSatish Balay     } ocomm;
3999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg));
40033779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute");
40133779a13SJunchao 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");
40233779a13SJunchao Zhang   }
4039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval));
4049566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm));
405da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
406b89831e5SBarry Smith }
407da3039f7SJed Brown 
408da3039f7SJed Brown /*
40933779a13SJunchao 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.
410da3039f7SJed Brown  */
411d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
412d71ae5a4SJacob Faibussowitsch {
413da3039f7SJed Brown   PetscFunctionBegin;
4149566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm));
415480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
416480cf27aSJed Brown }
417480cf27aSJed Brown 
41833779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *);
41942218b76SBarry Smith 
420951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
4218cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *);
4228cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
4238cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
424e39fd77fSBarry Smith #endif
425e39fd77fSBarry Smith 
4260f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE;
42712801b39SBarry Smith 
428eb27c7c8SSatish Balay PETSC_INTERN int    PetscGlobalArgc;
429eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4306ae9a8a6SBarry Smith int                 PetscGlobalArgc = 0;
43102c9f0b5SLisandro Dalcin char              **PetscGlobalArgs = NULL;
432dff31646SBarry Smith PetscSegBuffer      PetscCitationsList;
433e5c89e4eSSatish Balay 
434d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void)
435d71ae5a4SJacob Faibussowitsch {
436051e4cf2SJed Brown   PetscFunctionBegin;
4379566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList));
43830c35bf2SSatish Balay 
43930c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\
44030c35bf2SSatish Balay   Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\
44130c35bf2SSatish Balay     and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\
4423f81df79SBarry Smith     and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\
44330c35bf2SSatish Balay     and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\
44430c35bf2SSatish Balay     and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\
44530c35bf2SSatish Balay     and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\
44630c35bf2SSatish Balay     and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\
44730c35bf2SSatish Balay   Title = {{PETSc/TAO} Users Manual},\n\
4488f23a7d0SSatish Balay   Number = {ANL-21/39 - Revision 3.18},\n\
44930c35bf2SSatish Balay   Institution = {Argonne National Laboratory},\n\
4509371c9d4SSatish Balay   Year = {2022}\n}\n",
4519371c9d4SSatish Balay                                    NULL));
45230c35bf2SSatish Balay 
45330c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\
45430c35bf2SSatish Balay   Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\
45530c35bf2SSatish Balay   Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\
45630c35bf2SSatish Balay   Booktitle = {Modern Software Tools in Scientific Computing},\n\
45730c35bf2SSatish Balay   Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\
45830c35bf2SSatish Balay   Pages = {163--202},\n\
45930c35bf2SSatish Balay   Publisher = {Birkh{\\\"{a}}user Press},\n\
4609371c9d4SSatish Balay   Year = {1997}\n}\n",
4619371c9d4SSatish Balay                                    NULL));
46230c35bf2SSatish Balay 
4633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
464051e4cf2SJed Brown }
465e5c89e4eSSatish Balay 
4662d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4672d747510SLisandro Dalcin 
468d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[])
469d71ae5a4SJacob Faibussowitsch {
4702d747510SLisandro Dalcin   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(programname, name, sizeof(programname)));
4723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4732d747510SLisandro Dalcin }
4742d747510SLisandro Dalcin 
4752d747510SLisandro Dalcin /*@C
4762d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4772d747510SLisandro Dalcin 
4782d747510SLisandro Dalcin     Not Collective
4792d747510SLisandro Dalcin 
4802d747510SLisandro Dalcin     Input Parameter:
4812d747510SLisandro Dalcin .   len - length of the string name
4822d747510SLisandro Dalcin 
4832d747510SLisandro Dalcin     Output Parameter:
484811af0c4SBarry Smith .   name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN`
4852d747510SLisandro Dalcin 
4862d747510SLisandro Dalcin    Level: advanced
4872d747510SLisandro Dalcin 
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 
514db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`
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 
531793721a6SBarry Smith    Output Parameters:
532793721a6SBarry Smith .  args - the command line arguments
533793721a6SBarry Smith 
534793721a6SBarry Smith    Level: intermediate
535793721a6SBarry Smith 
536793721a6SBarry Smith    Notes:
537793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
538793721a6SBarry Smith 
539db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`
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 
562793721a6SBarry Smith    Output Parameters:
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));
9099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128));
9109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128));
91185649d77SJunchao Zhang #endif
9129e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
9139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16));
9149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FP16));
91585649d77SJunchao Zhang #endif
91685649d77SJunchao Zhang 
91785649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
9189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM));
919613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX));
920613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN));
9219e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
9229e517322SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128));
92385649d77SJunchao Zhang #endif
92485649d77SJunchao Zhang 
9259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR));
92662e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp));
9279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR));
92885649d77SJunchao Zhang 
92985649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
93085649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
93185649d77SJunchao Zhang   {
93285649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
93393d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)};
93485649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_REAL, MPIU_INT}, tmpStruct;
93585649d77SJunchao Zhang 
9369566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
93793d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT));
9389566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9399566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT));
94085649d77SJunchao Zhang   }
94185649d77SJunchao Zhang   {
94285649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
94393d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)};
94485649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_SCALAR, MPIU_INT}, tmpStruct;
94585649d77SJunchao Zhang 
9469566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
94793d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT));
9489566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9499566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
95085649d77SJunchao Zhang   }
95185649d77SJunchao Zhang #endif
95285649d77SJunchao Zhang 
95385649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT));
9559566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2INT));
95685649d77SJunchao Zhang #endif
9579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT));
9589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPI_4INT));
9599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT));
9609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_4INT));
96185649d77SJunchao Zhang 
96285649d77SJunchao Zhang   /*
96385649d77SJunchao Zhang      Attributes to be set on PETSc communicators
96485649d77SJunchao Zhang   */
9659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0));
9669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0));
9679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0));
9689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0));
96962e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, (void *)0));
97062e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, (void *)0));
97185649d77SJunchao Zhang 
97285649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
9739566063dSJacob Faibussowitsch   if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len));
97485649d77SJunchao Zhang   else
97585649d77SJunchao Zhang #endif
9769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file));
97785649d77SJunchao Zhang 
97885649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9799566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
98085649d77SJunchao Zhang 
98185649d77SJunchao Zhang   /*
98285649d77SJunchao Zhang      Check system options and print help
98385649d77SJunchao Zhang   */
9849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCheckInitial_Private(help));
98585649d77SJunchao Zhang 
986a4af0ceeSJacob Faibussowitsch   /*
9870e6b6b59SJacob Faibussowitsch     Creates the logging data structures; this is enabled even if logging is not turned on
9880e6b6b59SJacob Faibussowitsch     This is the last thing we do before returning to the user code to prevent having the
9890e6b6b59SJacob Faibussowitsch     logging numbers contaminated by any startup time associated with MPI
9900e6b6b59SJacob Faibussowitsch   */
9910e6b6b59SJacob Faibussowitsch #if defined(PETSC_USE_LOG)
9920e6b6b59SJacob Faibussowitsch   PetscCall(PetscLogInitialize());
9930e6b6b59SJacob Faibussowitsch #endif
9940e6b6b59SJacob Faibussowitsch 
9950e6b6b59SJacob Faibussowitsch   /*
996a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
997a4af0ceeSJacob Faibussowitsch 
998a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
999a4af0ceeSJacob Faibussowitsch    1. MPI initialized
1000a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
10010e6b6b59SJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up
10020e6b6b59SJacob Faibussowitsch       its own SIGSEV handler via the push/pop interface.
10030e6b6b59SJacob Faibussowitsch    4. Logging initialized
1004a4af0ceeSJacob Faibussowitsch   */
10059566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
1006a4af0ceeSJacob Faibussowitsch 
1007a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
1008a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
10099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_summary", &flg));
10109566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg));
10119566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL));
1012a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
10139566063dSJacob Faibussowitsch   PetscCall(PetscViennaCLInit());
1014a4af0ceeSJacob Faibussowitsch #endif
1015a4af0ceeSJacob Faibussowitsch 
10169566063dSJacob Faibussowitsch   PetscCall(PetscCitationsInitialize());
101785649d77SJunchao Zhang 
101885649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
10199566063dSJacob Faibussowitsch   PetscCall(PetscInitializeSAWs(ftn ? NULL : help));
102027104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
10219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg));
10229566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscStackViewSAWs());
102385649d77SJunchao Zhang #endif
102485649d77SJunchao Zhang 
102585649d77SJunchao Zhang   /*
102685649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
102785649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
102885649d77SJunchao Zhang   */
10299566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_DynamicLibraries());
103085649d77SJunchao Zhang 
10319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10329566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size));
103396dcfd6cSBarry Smith   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
10349566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname));
103585649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
103685649d77SJunchao Zhang   {
103785649d77SJunchao Zhang     PetscBool omp_view_flag;
103885649d77SJunchao Zhang     char     *threads = getenv("OMP_NUM_THREADS");
103985649d77SJunchao Zhang 
104085649d77SJunchao Zhang     if (threads) {
10419566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads));
104285649d77SJunchao Zhang       (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads);
104385649d77SJunchao Zhang     } else {
10442f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt)omp_get_max_threads();
10459566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads));
104685649d77SJunchao Zhang     }
1047d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys");
10489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg));
10499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag));
1050d0609cedSBarry Smith     PetscOptionsEnd();
105185649d77SJunchao Zhang     if (flg) {
10529566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads));
105385649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
105485649d77SJunchao Zhang     }
105548a46eb9SPierre Jolivet     if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads));
105685649d77SJunchao Zhang   }
105785649d77SJunchao Zhang #endif
105885649d77SJunchao Zhang 
105985649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
106085649d77SJunchao Zhang   /*
106185649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
106285649d77SJunchao Zhang 
106385649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
106485649d77SJunchao Zhang   */
10659566063dSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL));
106685649d77SJunchao Zhang #endif
106785649d77SJunchao Zhang 
106885649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10699566063dSJacob Faibussowitsch   PetscCall(PetscFPTCreate(10000));
107085649d77SJunchao Zhang #endif
107185649d77SJunchao Zhang 
107285649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
107385649d77SJunchao Zhang   {
107485649d77SJunchao Zhang     PetscViewer viewer;
10759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg));
107685649d77SJunchao Zhang     if (flg) {
10779566063dSJacob Faibussowitsch       PetscCall(PetscProcessPlacementView(viewer));
10789566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
107985649d77SJunchao Zhang     }
108085649d77SJunchao Zhang   }
108185649d77SJunchao Zhang #endif
108285649d77SJunchao Zhang 
108385649d77SJunchao Zhang   flg = PETSC_TRUE;
10849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL));
10859566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE));
108685649d77SJunchao Zhang 
108785649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10883ba16761SJacob Faibussowitsch   PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD);
10893ba16761SJacob Faibussowitsch   PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default);
10903ba16761SJacob Faibussowitsch   PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", "");
10913ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, "");
109285649d77SJunchao Zhang #endif
109385649d77SJunchao Zhang 
109485649d77SJunchao Zhang #if defined(__VALGRIND_H)
109585649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE;
109685649d77SJunchao Zhang   #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
10979566063dSJacob Faibussowitsch   if (PETSC_RUNNING_ON_VALGRIND) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING: Running valgrind with the MacOS native BLAS and LAPACK can fail. If it fails suggest configuring with --download-fblaslapack or --download-f2cblaslapack"));
109885649d77SJunchao Zhang   #endif
109985649d77SJunchao Zhang #endif
110085649d77SJunchao Zhang   /*
110185649d77SJunchao Zhang       Set flag that we are completely initialized
110285649d77SJunchao Zhang   */
110385649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
110485649d77SJunchao Zhang 
11059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg));
11069566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPythonInitialize(NULL, NULL));
1107f1f2ae84SBarry Smith 
1108f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1109f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin());
1110f1f2ae84SBarry Smith   else PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_SUP, "PETSc configured using -with-single-library=0; -mpi_linear_solver_server not supported in that case");
11113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111285649d77SJunchao Zhang }
111385649d77SJunchao Zhang 
1114e5c89e4eSSatish Balay /*@C
1115e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1116811af0c4SBarry Smith    `PetscInitialize()` calls MPI_Init() if that has yet to be called,
1117e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1118e5c89e4eSSatish Balay    your program -- usually the very first line!
1119e5c89e4eSSatish Balay 
1120811af0c4SBarry Smith    Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set
1121e5c89e4eSSatish Balay 
1122e5c89e4eSSatish Balay    Input Parameters:
1123e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1124e5c89e4eSSatish Balay .  args - the command line arguments
1125be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1126be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1127be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1128c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
11290298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1130e5c89e4eSSatish Balay 
1131811af0c4SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that
1132811af0c4SBarry Smith    communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a
1133811af0c4SBarry Smith    four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not,
1134811af0c4SBarry Smith    then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even
113505827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1136e5c89e4eSSatish Balay 
1137e5c89e4eSSatish Balay    Options Database Keys:
11387ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
11397ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1140e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11418a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
1142811af0c4SBarry Smith .  -on_error_abort - calls `abort()` when error detected (no traceback)
1143811af0c4SBarry Smith .  -on_error_mpiabort - calls `MPI_abort()` when error detected
11441af3601dSBarry Smith .  -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr
11458a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1146bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1147e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1148e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1149e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
115079dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
115179dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
1152811af0c4SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see `PetscMallocSetDebug()`
1153aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
115492f119d6SBarry 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
1155811af0c4SBarry Smith .  -malloc_view - show a list of all allocated memory during `PetscFinalize()`
115692f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1157608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
115892f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1159e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1160e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1161e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1162e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1163e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
11640841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1165e5c89e4eSSatish Balay 
1166c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1167c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1168c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1169811af0c4SBarry Smith -  -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()`
1170c5b5d8d5SVaclav Hapla 
1171c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1172c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1173c5b5d8d5SVaclav Hapla    They can be used also in option files.
1174c5b5d8d5SVaclav Hapla 
1175811af0c4SBarry Smith    See `PetscOptionsMonitorSet()` to do monitoring programmatically.
1176c5b5d8d5SVaclav Hapla 
1177e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1178a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1179811af0c4SBarry Smith +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`.
1180217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1181217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1182495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1183811af0c4SBarry Smith         hangs without running in the debugger).  See `PetscLogTraceBegin()`.
1184811af0c4SBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see `PetscLogView()`.
1185811af0c4SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`.
1186811af0c4SBarry Smith .  -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView().
11879a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1188495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1189f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1190811af0c4SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See `PetscLogDump()`.
1191811af0c4SBarry Smith .  -log [filename] - Logs basic profiline information  See `PetscLogDump()`.
1192c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
1193811af0c4SBarry Smith .  -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off
1194c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1195495fc317SBarry Smith 
11966a77f485SPierre Jolivet     Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time
1197e5c89e4eSSatish Balay 
1198ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1199ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1200ffbd1cfbSBarry 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
1201ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1202ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1203ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1204ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1205ffbd1cfbSBarry Smith 
1206e5c89e4eSSatish Balay    Environmental Variables:
1207811af0c4SBarry Smith +   `PETSC_TMP` - alternative tmp directory
1208811af0c4SBarry Smith .   `PETSC_SHARED_TMP` - tmp is shared by all processes
1209811af0c4SBarry Smith .   `PETSC_NOT_SHARED_TMP` - each process has its own private tmp
1210811af0c4SBarry Smith .   `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs
1211811af0c4SBarry Smith .   `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1212811af0c4SBarry Smith .   `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer
1213811af0c4SBarry Smith -   `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to
1214e5c89e4eSSatish Balay 
1215e5c89e4eSSatish Balay    Level: beginner
1216e5c89e4eSSatish Balay 
1217811af0c4SBarry Smith    Note:
1218811af0c4SBarry Smith    If for some reason you must call `MPI_Init()` separately, call
1219811af0c4SBarry Smith    it before `PetscInitialize()`.
1220e5c89e4eSSatish Balay 
1221811af0c4SBarry Smith    Fortran Notes:
1222811af0c4SBarry Smith    In Fortran this routine can be called with
1223811af0c4SBarry Smith .vb
1224811af0c4SBarry Smith        call PetscInitialize(ierr)
1225811af0c4SBarry Smith        call PetscInitialize(file,ierr) or
1226811af0c4SBarry Smith        call PetscInitialize(file,help,ierr)
1227811af0c4SBarry Smith .ve
1228e5c89e4eSSatish Balay 
1229811af0c4SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after
1230811af0c4SBarry Smith    calling `PetscInitialize()`.
1231e5c89e4eSSatish Balay 
12325dedd997SBarry Smith    Options Database Key for Developers:
12335dedd997SBarry Smith .  -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form:
12345dedd997SBarry Smith     "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py)
12355dedd997SBarry Smith 
1236db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()`
1237e5c89e4eSSatish Balay @*/
1238d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[])
1239d71ae5a4SJacob Faibussowitsch {
124085649d77SJunchao Zhang   PetscMPIInt flag;
124121ef0414SBarry Smith   const char *prog = "Unknown Name", *mpienv;
1242e5c89e4eSSatish Balay 
124327104ee2SJacob Faibussowitsch   PetscFunctionBegin;
12443ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
12459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Initialized(&flag));
1246e5c89e4eSSatish Balay   if (!flag) {
124739a651e2SJacob 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");
12489566063dSJacob Faibussowitsch     PetscCall(PetscPreMPIInit_Private());
12495e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12505e765c61SJed Brown     {
125139a651e2SJacob Faibussowitsch       PetscMPIInt PETSC_UNUSED provided;
12529566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED, &provided));
12535e765c61SJed Brown     }
12545e765c61SJed Brown #else
12559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Init(argc, args));
12565e765c61SJed Brown #endif
125721ef0414SBarry Smith     if (PetscDefined(HAVE_MPIUNI)) {
125821ef0414SBarry Smith       mpienv = getenv("PMI_SIZE");
125921ef0414SBarry Smith       if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE");
126021ef0414SBarry Smith       if (mpienv) {
126121ef0414SBarry Smith         PetscInt isize;
126221ef0414SBarry Smith         PetscCall(PetscOptionsStringToInt(mpienv, &isize));
126321ef0414SBarry 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");
126421ef0414SBarry 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");
126521ef0414SBarry Smith       }
126621ef0414SBarry Smith     }
1267e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1268e5c89e4eSSatish Balay   }
12694dfee713SSatish Balay 
127085649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1271e5c89e4eSSatish Balay   if (argc && args) {
1272e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1273e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1274e5c89e4eSSatish Balay   }
1275811af0c4SBarry Smith   PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0));
12763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1277e5c89e4eSSatish Balay }
1278e5c89e4eSSatish Balay 
1279a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
128095c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1281ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsCounts;
1282ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsMaxCounts;
128395c0884eSLisandro Dalcin PETSC_INTERN PetscBool    PetscObjectsLog;
12844097062eSBarry Smith #endif
1285e5c89e4eSSatish Balay 
1286008a6e76SBarry Smith /*
1287008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1288008a6e76SBarry Smith */
1289d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void)
1290d71ae5a4SJacob Faibussowitsch {
1291008a6e76SBarry Smith   PetscFunctionBegin;
1292613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
12939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128));
12949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128));
1295008a6e76SBarry Smith #endif
12969e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
12979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FP16));
1298008a6e76SBarry Smith #endif
1299008a6e76SBarry Smith 
1300de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
13019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_SUM));
1302613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MAX));
1303613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MIN));
13049e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
13059e517322SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128));
1306008a6e76SBarry Smith #endif
1307008a6e76SBarry Smith 
13089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR));
13099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT));
13109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT));
131140df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
13129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2INT));
1313008a6e76SBarry Smith #endif
13149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPI_4INT));
13159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_4INT));
13169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
131762e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp));
13183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1319008a6e76SBarry Smith }
1320008a6e76SBarry Smith 
1321a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
1322a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1323a4af0ceeSJacob Faibussowitsch #endif
1324a4af0ceeSJacob Faibussowitsch 
1325e5c89e4eSSatish Balay /*@C
1326e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1327811af0c4SBarry Smith    of the program. `MPI_Finalize()` is called only if the user had not
1328811af0c4SBarry Smith    called `MPI_Init()` before calling `PetscInitialize()`.
1329e5c89e4eSSatish Balay 
1330811af0c4SBarry Smith    Collective on `PETSC_COMM_WORLD`
1331e5c89e4eSSatish Balay 
1332e5c89e4eSSatish Balay    Options Database Keys:
1333811af0c4SBarry Smith +  -options_view - Calls `PetscOptionsView()`
1334e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
13357eb1d149SBarry 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
1336e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1337811af0c4SBarry Smith .  -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed
1338e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
133992f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1340e5c89e4eSSatish Balay 
1341e5c89e4eSSatish Balay    Level: beginner
1342e5c89e4eSSatish Balay 
1343e5c89e4eSSatish Balay    Note:
1344811af0c4SBarry Smith    See `PetscInitialize()` for other runtime options.
1345e5c89e4eSSatish Balay 
1346db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()`
1347e5c89e4eSSatish Balay @*/
1348d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void)
1349d71ae5a4SJacob Faibussowitsch {
13504bb5149bSJed Brown   PetscMPIInt rank;
1351a8d2bbe5SBarry Smith   PetscInt    nopt;
13522bf49c77SBarry Smith   PetscBool   flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE;
1353dff31646SBarry Smith   PetscBool   flg;
135410463e74SBarry Smith #if defined(PETSC_USE_LOG)
135510463e74SBarry Smith   char mname[PETSC_MAX_PATH_LEN];
135610463e74SBarry Smith #endif
1357e5c89e4eSSatish Balay 
13583cbbc5ffSBarry Smith   PetscFunctionBegin;
135939a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()");
13609566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PetscFinalize() called\n"));
1361b022a5c1SBarry Smith 
1362f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1363f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd());
1364f1f2ae84SBarry Smith 
136562e5d2d2SJDBetteridge   /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */
136662e5d2d2SJDBetteridge   {
136762e5d2d2SJDBetteridge     union
136862e5d2d2SJDBetteridge     {
136962e5d2d2SJDBetteridge       MPI_Comm comm;
137062e5d2d2SJDBetteridge       void    *ptr;
137162e5d2d2SJDBetteridge     } ucomm;
137262e5d2d2SJDBetteridge     PetscMPIInt flg;
137362e5d2d2SJDBetteridge     void       *tmp;
137462e5d2d2SJDBetteridge 
137562e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
137662e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
137762e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF));
137862e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
137962e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
138062e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD));
138162e5d2d2SJDBetteridge   }
138262e5d2d2SJDBetteridge 
13839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1384a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13853ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE);
13863ba16761SJacob Faibussowitsch   PetscCallExternal(adios_finalize, rank);
1387a56f64adSBarry Smith #endif
13889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg));
1389dff31646SBarry Smith   if (flg) {
13901f817a21SBarry Smith     char *cits, filename[PETSC_MAX_PATH_LEN];
13911f817a21SBarry Smith     FILE *fd = PETSC_STDOUT;
13921f817a21SBarry Smith 
13939566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL));
139448a46eb9SPierre Jolivet     if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd));
13959566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits));
1396dff31646SBarry Smith     cits[0] = 0;
13979566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits));
13989566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n"));
13999566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14009566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits));
14019566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14029566063dSJacob Faibussowitsch     PetscCall(PetscFClose(PETSC_COMM_WORLD, fd));
14039566063dSJacob Faibussowitsch     PetscCall(PetscFree(cits));
1404dff31646SBarry Smith   }
14059566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&PetscCitationsList));
1406dff31646SBarry Smith 
1407c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
140804102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
140904102261SBarry Smith   {
141004102261SBarry Smith     PetscInt nmax = 2;
141104102261SBarry Smith     char   **buffs;
14129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
14139566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-textbelt", buffs, &nmax, &flg1));
141404102261SBarry Smith     if (flg1) {
141528b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-textbelt requires either the phone number or number,\"message\"");
141604102261SBarry Smith       if (nmax == 1) {
1417c6a7a370SJeremy L Thompson         size_t len = 128;
1418c6a7a370SJeremy L Thompson         PetscCall(PetscMalloc1(len, &buffs[1]));
14199566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
1420c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(buffs[1], " has completed", len));
142104102261SBarry Smith       }
14229566063dSJacob Faibussowitsch       PetscCall(PetscTextBelt(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
14239566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
14249566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
142504102261SBarry Smith     }
14269566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
142704102261SBarry Smith   }
1428763ec7b1SBarry Smith   {
1429763ec7b1SBarry Smith     PetscInt nmax = 2;
1430763ec7b1SBarry Smith     char   **buffs;
14319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
14329566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-tellmycell", buffs, &nmax, &flg1));
1433763ec7b1SBarry Smith     if (flg1) {
143428b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-tellmycell requires either the phone number or number,\"message\"");
1435763ec7b1SBarry Smith       if (nmax == 1) {
1436c6a7a370SJeremy L Thompson         size_t len = 128;
1437c6a7a370SJeremy L Thompson         PetscCall(PetscMalloc1(len, &buffs[1]));
14389566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
1439c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(buffs[1], " has completed", len));
1440763ec7b1SBarry Smith       }
14419566063dSJacob Faibussowitsch       PetscCall(PetscTellMyCell(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
14429566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
14439566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
1444763ec7b1SBarry Smith     }
14459566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
1446763ec7b1SBarry Smith   }
144704102261SBarry Smith #endif
144804102261SBarry Smith 
14492d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14509566063dSJacob Faibussowitsch   PetscCall(PetscFPTDestroy());
14512d53ad75SBarry Smith #endif
14522d53ad75SBarry Smith 
1453e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1454dff31646SBarry Smith   flg = PETSC_FALSE;
14559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL));
14561baa6e33SBarry Smith   if (flg) PetscCall(PetscOptionsSAWsDestroy());
1457d5649816SBarry Smith #endif
1458d5649816SBarry Smith 
1459681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1460681455b2SBarry Smith   flg1 = PETSC_FALSE;
14619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL));
1462681455b2SBarry Smith   if (flg1) {
1463681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
14649566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL));
1465681455b2SBarry Smith   }
1466681455b2SBarry Smith #endif
1467681455b2SBarry Smith 
146867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_info", &flg2, NULL));
1470e5c89e4eSSatish Balay   if (!flg2) {
147190d69ab7SBarry Smith     flg2 = PETSC_FALSE;
14729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL));
1473e5c89e4eSSatish Balay   }
147448a46eb9SPierre Jolivet   if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n"));
147567584ceeSBarry Smith #endif
1476e5c89e4eSSatish Balay 
1477e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
147890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL));
1480e5c89e4eSSatish Balay   if (flg1) {
1481e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
14829566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD));
14839566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops));
1484e5c89e4eSSatish Balay   }
1485e5c89e4eSSatish Balay #endif
1486e5c89e4eSSatish Balay 
1487e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1488e5c89e4eSSatish Balay   #if defined(PETSC_HAVE_MPE)
1489e5c89e4eSSatish Balay   mname[0] = 0;
14909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1));
1491e5c89e4eSSatish Balay   if (flg1) {
14929566063dSJacob Faibussowitsch     if (mname[0]) PetscCall(PetscLogMPEDump(mname));
14939566063dSJacob Faibussowitsch     else PetscCall(PetscLogMPEDump(0));
1494e5c89e4eSSatish Balay   }
1495e5c89e4eSSatish Balay   #endif
1496356e5837SBarry Smith #endif
1497a297a907SKarl Rupp 
1498dd710f27SBarry Smith   /*
1499dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1500dd710f27SBarry Smith   */
15019566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1502dd710f27SBarry Smith 
1503356e5837SBarry Smith #if defined(PETSC_USE_LOG)
15049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE));
15059566063dSJacob Faibussowitsch   PetscCall(PetscLogViewFromOptions());
15069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPopGetViewerOff());
150787c3beb6SLisandro Dalcin 
1508356e5837SBarry Smith   mname[0] = 0;
15099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_summary", mname, sizeof(mname), &flg1));
1510e5c89e4eSSatish Balay   if (flg1) {
151191eabc43SBarry Smith     PetscViewer viewer;
15129566063dSJacob Faibussowitsch     PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD, "\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n"));
151391eabc43SBarry Smith     if (mname[0]) {
15149566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, mname, &viewer));
15159566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
15169566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
151733f85c2fSBarry Smith     } else {
151833f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
15199566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
15209566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
15219566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
152233f85c2fSBarry Smith     }
1523e5c89e4eSSatish Balay   }
1524a297a907SKarl Rupp 
1525dd710f27SBarry Smith   /*
1526dd710f27SBarry Smith      Free any objects created by the last block of code.
1527dd710f27SBarry Smith   */
15289566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1529dd710f27SBarry Smith 
1530dd710f27SBarry Smith   mname[0] = 0;
15319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1));
15329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2));
15339566063dSJacob Faibussowitsch   if (flg1 || flg2) PetscCall(PetscLogDump(mname));
1534e5c89e4eSSatish Balay #endif
153510463e74SBarry Smith 
153690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL));
15389566063dSJacob Faibussowitsch   if (!flg1) PetscCall(PetscPopSignalHandler());
153990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL));
15411baa6e33SBarry Smith   if (flg1) PetscCall(PetscMPIDump(stdout));
154290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
154390d69ab7SBarry Smith   flg2 = PETSC_FALSE;
15448bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
15459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1));
15469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
15479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL));
1548e4c476e2SSatish Balay 
1549e5c89e4eSSatish Balay   if (flg2) {
1550be56827dSJed Brown     PetscViewer viewer;
15519566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15529566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15539566063dSJacob Faibussowitsch     PetscCall(PetscOptionsView(NULL, viewer));
15549566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1555e5c89e4eSSatish Balay   }
1556e5c89e4eSSatish Balay 
1557e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
15589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1));
15599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1));
1560e5c89e4eSSatish Balay 
156133fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
15629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1));
15639245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1564e5c89e4eSSatish Balay   if (flg3) {
15653de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1566be56827dSJed Brown       PetscViewer viewer;
15679566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15689566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15699566063dSJacob Faibussowitsch       PetscCall(PetscOptionsView(NULL, viewer));
15709566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
1571e5c89e4eSSatish Balay     }
15729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsAllUsed(NULL, &nopt));
15733de2bfdfSBarry Smith     if (nopt) {
15749566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n"));
15759566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n"));
15763de2bfdfSBarry Smith       if (nopt == 1) {
15779566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n"));
1578e5c89e4eSSatish Balay       } else {
15799566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt));
1580e5c89e4eSSatish Balay       }
15813de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15829566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n"));
1583df12ba86SBarry Smith     }
15849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsLeft(NULL));
1585e5c89e4eSSatish Balay   }
1586e5c89e4eSSatish Balay 
1587e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1588d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15899566063dSJacob Faibussowitsch     PetscCall(PetscStackSAWsViewOff());
1590792fecdfSBarry Smith     PetscCallSAWs(SAWs_Finalize, ());
1591d45a07a7SBarry Smith   }
1592ec957eceSBarry Smith #endif
1593ec957eceSBarry Smith 
15944097062eSBarry Smith #if defined(PETSC_USE_LOG)
159510463e74SBarry Smith   /*
1596dbc8283eSBarry Smith        List all objects the user may have forgot to free
15972eff7a51SBarry Smith   */
159805df10baSBarry Smith   if (PetscObjectsLog) {
15999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
1600a64a8e02SBarry Smith     if (flg1) {
1601a64a8e02SBarry Smith       MPI_Comm local_comm;
16027eb1d149SBarry Smith       char     string[64];
1603a64a8e02SBarry Smith 
16049566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL));
1605afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16069566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16079566063dSJacob Faibussowitsch       PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
16089566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16099566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
16100a1571b3SBarry Smith     }
161105df10baSBarry Smith   }
16124097062eSBarry Smith #endif
16134097062eSBarry Smith 
16144097062eSBarry Smith #if defined(PETSC_USE_LOG)
1615dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1616dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
16179566063dSJacob Faibussowitsch   PetscCall(PetscFree(PetscObjects));
16184097062eSBarry Smith #endif
16192eff7a51SBarry Smith 
162033f85c2fSBarry Smith   /*
162133f85c2fSBarry Smith      Destroy any packages that registered a finalize
162233f85c2fSBarry Smith   */
16239566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalizeAll());
162433f85c2fSBarry Smith 
1625101409b8SToby Isaac #if defined(PETSC_USE_LOG)
16269566063dSJacob Faibussowitsch   PetscCall(PetscLogFinalize());
1627101409b8SToby Isaac #endif
1628101409b8SToby Isaac 
162933f85c2fSBarry Smith   /*
163048dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
163148dd1dffSBarry Smith   */
16322e956fe4SStefano Zampini   if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll());
163337e93019SBarry Smith 
16344028d114SSatish Balay   if (petsc_history) {
16359566063dSJacob Faibussowitsch     PetscCall(PetscCloseHistoryFile(&petsc_history));
163602c9f0b5SLisandro Dalcin     petsc_history = NULL;
1637e5c89e4eSSatish Balay   }
16389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
16399566063dSJacob Faibussowitsch   PetscCall(PetscInfoDestroy());
1640e5c89e4eSSatish Balay 
164167584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
164292f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1643e5c89e4eSSatish Balay     char  fname[PETSC_MAX_PATH_LEN];
164492f119d6SBarry Smith     char  sname[PETSC_MAX_PATH_LEN];
1645e5c89e4eSSatish Balay     FILE *fd;
1646ed9cf6e9SBarry Smith     int   err;
1647e5c89e4eSSatish Balay 
1648dc92acbaSJed Brown     flg2 = PETSC_FALSE;
164992f119d6SBarry Smith     flg3 = PETSC_FALSE;
16509566063dSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL));
16519566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL));
165292f119d6SBarry Smith     fname[0] = 0;
16539566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1));
1654e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
16553ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
16569371c9d4SSatish Balay       fd = fopen(sname, "w");
16579371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16589566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(fd));
1659ed9cf6e9SBarry Smith       err = fclose(fd);
166028b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
166192f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1662e5c89e4eSSatish Balay       MPI_Comm local_comm;
1663e5c89e4eSSatish Balay 
1664afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16659566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16669566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(stdout));
16679566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16689566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1669e5c89e4eSSatish Balay     }
1670e5c89e4eSSatish Balay     fname[0] = 0;
16719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1));
1672e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
16733ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
16749371c9d4SSatish Balay       fd = fopen(sname, "w");
16759371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16769566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(fd));
1677ed9cf6e9SBarry Smith       err = fclose(fd);
167828b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
167992f119d6SBarry Smith     } else if (flg1) {
168092f119d6SBarry Smith       MPI_Comm local_comm;
168192f119d6SBarry Smith 
1682afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16839566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16849566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(stdout));
16859566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16869566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1687e5c89e4eSSatish Balay     }
1688e5c89e4eSSatish Balay   }
168967584ceeSBarry Smith #endif
169020e2c332SMatthew G. Knepley 
16915486ca60SMatthew G. Knepley   /*
16925486ca60SMatthew G. Knepley      Close any open dynamic libraries
16935486ca60SMatthew G. Knepley   */
16949566063dSJacob Faibussowitsch   PetscCall(PetscFinalize_DynamicLibraries());
16955486ca60SMatthew G. Knepley 
1696e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroyDefault());
1698e5c89e4eSSatish Balay 
1699e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
170002c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1701e5c89e4eSSatish Balay 
1702c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1703c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
17049566063dSJacob Faibussowitsch     PetscCall(PetscKokkosFinalize_Private());
1705c2b86a48SJunchao Zhang     PetscBeganKokkos       = PETSC_FALSE;
170645639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1707c2b86a48SJunchao Zhang   }
1708c2b86a48SJunchao Zhang #endif
1709c2b86a48SJunchao Zhang 
171071438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
171171438e86SJunchao Zhang   if (PetscBeganNvshmem) {
17129566063dSJacob Faibussowitsch     PetscCall(PetscNvshmemFinalize());
171371438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
171471438e86SJunchao Zhang   }
171571438e86SJunchao Zhang #endif
1716a0e72f99SJunchao Zhang 
17179566063dSJacob Faibussowitsch   PetscCall(PetscFreeMPIResources());
1718e5c89e4eSSatish Balay 
1719dbc8283eSBarry Smith   /*
1720efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1721efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1722efb80d3cSBarry Smith 
1723efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1724efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1725dbc8283eSBarry Smith  */
1726b770b1f6SSatish Balay   {
1727dbc8283eSBarry Smith     PetscCommCounter *counter;
1728dbc8283eSBarry Smith     PetscMPIInt       flg;
1729dbc8283eSBarry Smith     MPI_Comm          icomm;
17309371c9d4SSatish Balay     union
17319371c9d4SSatish Balay     {
17329371c9d4SSatish Balay       MPI_Comm comm;
17339371c9d4SSatish Balay       void    *ptr;
17349371c9d4SSatish Balay     } ucomm;
17359566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
1736dbc8283eSBarry Smith     if (flg) {
1737265f3f35SJed Brown       icomm = ucomm.comm;
17389566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
173928b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1740dbc8283eSBarry Smith 
17419566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval));
17429566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17439566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1744dbc8283eSBarry Smith     }
17459566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
1746dbc8283eSBarry Smith     if (flg) {
1747265f3f35SJed Brown       icomm = ucomm.comm;
17489566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
174928b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1750dbc8283eSBarry Smith 
17519566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval));
17529566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17539566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1754dbc8283eSBarry Smith     }
1755b770b1f6SSatish Balay   }
1756dbc8283eSBarry Smith 
17579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
17589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
17599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
17609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
176162e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval));
176262e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval));
1763480cf27aSJed Brown 
1764*5ea2b939SDuncan Campbell   // Free keyvals which may be silently created by some routines
1765*5ea2b939SDuncan Campbell   if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval));
1766*5ea2b939SDuncan Campbell   if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval));
1767*5ea2b939SDuncan Campbell 
17689566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
17699566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
17709566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
17719566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock));
1772ef19f930SBarry Smith 
1773e5c89e4eSSatish Balay   if (PetscBeganMPI) {
177499b1327fSBarry Smith     PetscMPIInt flag;
17759566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalized(&flag));
177639a651e2SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
177739a651e2SJacob Faibussowitsch     /* wait until the very last moment to disable error handling */
177839a651e2SJacob Faibussowitsch     PetscErrorHandlingInitialized = PETSC_FALSE;
17799566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalize());
178039a651e2SJacob Faibussowitsch   } else PetscErrorHandlingInitialized = PETSC_FALSE;
178139a651e2SJacob Faibussowitsch 
1782e5c89e4eSSatish Balay   /*
1783e5c89e4eSSatish Balay 
1784e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1785e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1786e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1787e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1788e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17890e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1790e5c89e4eSSatish Balay    memory was not freed.
1791e5c89e4eSSatish Balay 
1792e5c89e4eSSatish Balay */
17939566063dSJacob Faibussowitsch   PetscCall(PetscMallocClear());
17949566063dSJacob Faibussowitsch   PetscCall(PetscStackReset());
1795a297a907SKarl Rupp 
1796e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1797e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
17987ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE)
179956883f60SBarry Smith   /*
180056883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
180156883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
180256883f60SBarry Smith    */
180356883f60SBarry Smith   __gcov_flush();
180456883f60SBarry Smith #endif
18051724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
18061724198aSStefano Zampini   PetscStackClearTop;
18073ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1808e5c89e4eSSatish Balay }
1809e5c89e4eSSatish Balay 
181043db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
1811d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b)
1812d71ae5a4SJacob Faibussowitsch {
181343db4dbbSBarry Smith   if (*a == *b) return 1;
181443db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
181543db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
181643db4dbbSBarry Smith   return 0;
181743db4dbbSBarry Smith }
1818a70650f6SBarry Smith #endif
181943db4dbbSBarry Smith 
182043db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
1821d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b)
1822d71ae5a4SJacob Faibussowitsch {
182343db4dbbSBarry Smith   if (*a == *b) return 1;
182443db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
182543db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
182643db4dbbSBarry Smith   return 0;
182743db4dbbSBarry Smith }
182843db4dbbSBarry Smith #endif
1829