xref: /petsc/src/sys/objects/pinit.c (revision 74c01ffa5774e67e0fa00e95c4737f6527ecf615)
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 
605ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedWD_keyval  = MPI_KEYVAL_INVALID;
615ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedTmp_keyval = MPI_KEYVAL_INVALID;
625ea2b939SDuncan 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];
267*74c01ffaSSatish Balay     #if defined(PETSC_HAVE_COMPLEX)
2689371c9d4SSatish Balay   } else if (*datatype == MPIU___COMPLEX128) {
269613bf2b2SPierre Jolivet     __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out;
270613bf2b2SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
271*74c01ffaSSatish Balay     #endif
272613bf2b2SPierre Jolivet   }
273613bf2b2SPierre Jolivet   #endif
2749e517322SPierre Jolivet   #if defined(PETSC_HAVE_REAL___FP16)
2759e517322SPierre Jolivet   else if (*datatype == MPIU___FP16) {
2769e517322SPierre Jolivet     __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out;
2779e517322SPierre Jolivet     for (i = 0; i < count; i++) xout[i] += xin[i];
2789e517322SPierre Jolivet   }
2799e517322SPierre Jolivet   #endif
2807c2de775SJed Brown   else {
2819e517322SPierre Jolivet   #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16)
2823ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SElF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
2839e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FP16)
2843ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types"));
2859e517322SPierre Jolivet   #elif !defined(PETSC_HAVE_REAL___FLOAT128)
2863ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types"));
2879e517322SPierre Jolivet   #else
2883ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types"));
289613bf2b2SPierre Jolivet   #endif
29041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
291e2e03761SBarry Smith   }
292812af9f3SBarry Smith   PetscFunctionReturnVoid();
293e5c89e4eSSatish Balay }
294e5c89e4eSSatish Balay #endif
295e5c89e4eSSatish Balay 
296570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
297d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
298d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
299d9822059SBarry Smith 
300d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
301d71ae5a4SJacob Faibussowitsch {
302d9822059SBarry Smith   PetscInt i, count = *cnt;
303d9822059SBarry Smith 
304d9822059SBarry Smith   PetscFunctionBegin;
3057c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3068c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
307a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]);
3087c2de775SJed Brown   }
3097c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3107c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3117c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
312ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3137c2de775SJed Brown   }
3147c2de775SJed Brown   #endif
3157c2de775SJed Brown   else {
3163ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"));
31741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3188c764dc5SJose Roman   }
319d9822059SBarry Smith   PetscFunctionReturnVoid();
320d9822059SBarry Smith }
321d9822059SBarry Smith 
322d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype)
323d71ae5a4SJacob Faibussowitsch {
324d9822059SBarry Smith   PetscInt i, count = *cnt;
325d9822059SBarry Smith 
326d9822059SBarry Smith   PetscFunctionBegin;
3277c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3288c764dc5SJose Roman     PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out;
329a297a907SKarl Rupp     for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]);
3307c2de775SJed Brown   }
3317c2de775SJed Brown   #if defined(PETSC_HAVE_COMPLEX)
3327c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3337c2de775SJed Brown     PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out;
334ad540459SPierre Jolivet     for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3357c2de775SJed Brown   }
3367c2de775SJed Brown   #endif
3377c2de775SJed Brown   else {
3383ba16761SJacob Faibussowitsch     PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"));
33941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG);
3408c764dc5SJose Roman   }
341d9822059SBarry Smith   PetscFunctionReturnVoid();
342d9822059SBarry Smith }
343d9822059SBarry Smith #endif
344d9822059SBarry Smith 
345480cf27aSJed Brown /*
346480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
347480cf27aSJed Brown 
348ff0e51ddSBarry 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.
349480cf27aSJed Brown 
35012801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
351480cf27aSJed Brown 
352480cf27aSJed Brown */
353d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state)
354d71ae5a4SJacob Faibussowitsch {
35505342407SJunchao Zhang   PetscCommCounter      *counter = (PetscCommCounter *)count_val;
35657f21012SBarry Smith   struct PetscCommStash *comms   = counter->comms, *pcomm;
357480cf27aSJed Brown 
358480cf27aSJed Brown   PetscFunctionBegin;
3599566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm));
3609566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter->iflags));
36157f21012SBarry Smith   while (comms) {
3629566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&comms->comm));
36357f21012SBarry Smith     pcomm = comms;
36457f21012SBarry Smith     comms = comms->next;
3659566063dSJacob Faibussowitsch     PetscCall(PetscFree(pcomm));
36657f21012SBarry Smith   }
3679566063dSJacob Faibussowitsch   PetscCallMPI(PetscFree(counter));
368480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
369480cf27aSJed Brown }
370480cf27aSJed Brown 
371480cf27aSJed Brown /*
37247435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
373da3039f7SJed Brown   calls MPI_Comm_free().
374da3039f7SJed Brown 
375da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
376480cf27aSJed Brown 
377ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
378480cf27aSJed Brown 
37912801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
380480cf27aSJed Brown 
381480cf27aSJed Brown */
382d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
383d71ae5a4SJacob Faibussowitsch {
3849371c9d4SSatish Balay   union
385480cf27aSJed Brown   {
3869371c9d4SSatish Balay     MPI_Comm comm;
3879371c9d4SSatish Balay     void    *ptr;
3889371c9d4SSatish Balay   } icomm;
389480cf27aSJed Brown 
390480cf27aSJed Brown   PetscFunctionBegin;
39112801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval");
392ec4fadc2SJed Brown   icomm.ptr = attr_val;
39376bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
39433779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
39533779a13SJunchao Zhang     PetscMPIInt flg;
3969371c9d4SSatish Balay     union
3979371c9d4SSatish Balay     {
3989371c9d4SSatish Balay       MPI_Comm comm;
3999371c9d4SSatish Balay       void    *ptr;
4009371c9d4SSatish Balay     } ocomm;
4019566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg));
40233779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute");
40333779a13SJunchao 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");
40433779a13SJunchao Zhang   }
4059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval));
4069566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm));
407da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
408b89831e5SBarry Smith }
409da3039f7SJed Brown 
410da3039f7SJed Brown /*
41133779a13SJunchao 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.
412da3039f7SJed Brown  */
413d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state)
414d71ae5a4SJacob Faibussowitsch {
415da3039f7SJed Brown   PetscFunctionBegin;
4169566063dSJacob Faibussowitsch   PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm));
417480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
418480cf27aSJed Brown }
419480cf27aSJed Brown 
42033779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *);
42142218b76SBarry Smith 
422951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
4238cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *);
4248cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
4258cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *);
426e39fd77fSBarry Smith #endif
427e39fd77fSBarry Smith 
4280f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE;
42912801b39SBarry Smith 
430eb27c7c8SSatish Balay PETSC_INTERN int    PetscGlobalArgc;
431eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4326ae9a8a6SBarry Smith int                 PetscGlobalArgc = 0;
43302c9f0b5SLisandro Dalcin char              **PetscGlobalArgs = NULL;
434dff31646SBarry Smith PetscSegBuffer      PetscCitationsList;
435e5c89e4eSSatish Balay 
436d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void)
437d71ae5a4SJacob Faibussowitsch {
438051e4cf2SJed Brown   PetscFunctionBegin;
4399566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList));
44030c35bf2SSatish Balay 
44130c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\
44230c35bf2SSatish Balay   Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\
44330c35bf2SSatish Balay     and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\
4443f81df79SBarry Smith     and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\
44530c35bf2SSatish Balay     and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\
44630c35bf2SSatish Balay     and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\
44730c35bf2SSatish Balay     and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\
44830c35bf2SSatish Balay     and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\
44930c35bf2SSatish Balay   Title = {{PETSc/TAO} Users Manual},\n\
450c4dc7257SSatish Balay   Number = {ANL-21/39 - Revision 3.19},\n\
45130c35bf2SSatish Balay   Institution = {Argonne National Laboratory},\n\
452c4dc7257SSatish Balay   Year = {2023}\n}\n",
4539371c9d4SSatish Balay                                    NULL));
45430c35bf2SSatish Balay 
45530c35bf2SSatish Balay   PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\
45630c35bf2SSatish Balay   Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\
45730c35bf2SSatish Balay   Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\
45830c35bf2SSatish Balay   Booktitle = {Modern Software Tools in Scientific Computing},\n\
45930c35bf2SSatish Balay   Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\
46030c35bf2SSatish Balay   Pages = {163--202},\n\
46130c35bf2SSatish Balay   Publisher = {Birkh{\\\"{a}}user Press},\n\
4629371c9d4SSatish Balay   Year = {1997}\n}\n",
4639371c9d4SSatish Balay                                    NULL));
46430c35bf2SSatish Balay 
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
466051e4cf2SJed Brown }
467e5c89e4eSSatish Balay 
4682d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4692d747510SLisandro Dalcin 
470d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[])
471d71ae5a4SJacob Faibussowitsch {
4722d747510SLisandro Dalcin   PetscFunctionBegin;
4739566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(programname, name, sizeof(programname)));
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4752d747510SLisandro Dalcin }
4762d747510SLisandro Dalcin 
4772d747510SLisandro Dalcin /*@C
4782d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4792d747510SLisandro Dalcin 
4802d747510SLisandro Dalcin     Not Collective
4812d747510SLisandro Dalcin 
4822d747510SLisandro Dalcin     Input Parameter:
4832d747510SLisandro Dalcin .   len - length of the string name
4842d747510SLisandro Dalcin 
4852d747510SLisandro Dalcin     Output Parameter:
486811af0c4SBarry Smith .   name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN`
4872d747510SLisandro Dalcin 
4882d747510SLisandro Dalcin    Level: advanced
4892d747510SLisandro Dalcin 
4902d747510SLisandro Dalcin @*/
491d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProgramName(char name[], size_t len)
492d71ae5a4SJacob Faibussowitsch {
4932d747510SLisandro Dalcin   PetscFunctionBegin;
4949566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(name, programname, len));
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4962d747510SLisandro Dalcin }
4972d747510SLisandro Dalcin 
498e5c89e4eSSatish Balay /*@C
499e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
500811af0c4SBarry Smith      after PetscInitialize() is called but before `PetscFinalize()`.
501e5c89e4eSSatish Balay 
502e5c89e4eSSatish Balay    Not Collective
503e5c89e4eSSatish Balay 
504e5c89e4eSSatish Balay    Output Parameters:
505e5c89e4eSSatish Balay +  argc - count of number of command line arguments
506e5c89e4eSSatish Balay -  args - the command line arguments
507e5c89e4eSSatish Balay 
508e5c89e4eSSatish Balay    Level: intermediate
509e5c89e4eSSatish Balay 
510e5c89e4eSSatish Balay    Notes:
511e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
512e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
513e5c89e4eSSatish Balay 
514f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
515f177e3b1SBarry Smith 
516db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`
517e5c89e4eSSatish Balay @*/
518d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args)
519d71ae5a4SJacob Faibussowitsch {
520e5c89e4eSSatish Balay   PetscFunctionBegin;
52139a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
522e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
523e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
5243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
525e5c89e4eSSatish Balay }
526e5c89e4eSSatish Balay 
527793721a6SBarry Smith /*@C
528793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
529811af0c4SBarry Smith      after `PetscInitialize()` is called but before `PetscFinalize()`.
530793721a6SBarry Smith 
531793721a6SBarry Smith    Not Collective
532793721a6SBarry Smith 
5332fe279fdSBarry Smith    Output Parameter:
534793721a6SBarry Smith .  args - the command line arguments
535793721a6SBarry Smith 
536793721a6SBarry Smith    Level: intermediate
537793721a6SBarry Smith 
538793721a6SBarry Smith    Notes:
539793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
540793721a6SBarry Smith 
541db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`
542793721a6SBarry Smith @*/
543d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args)
544d71ae5a4SJacob Faibussowitsch {
545793721a6SBarry Smith   PetscInt i, argc = PetscGlobalArgc;
546793721a6SBarry Smith 
547793721a6SBarry Smith   PetscFunctionBegin;
54839a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()");
5499371c9d4SSatish Balay   if (!argc) {
5509371c9d4SSatish Balay     *args = NULL;
5513ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5529371c9d4SSatish Balay   }
5539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(argc, args));
5549566063dSJacob Faibussowitsch   for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i]));
5552d747510SLisandro Dalcin   (*args)[argc - 1] = NULL;
5563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557793721a6SBarry Smith }
558793721a6SBarry Smith 
559793721a6SBarry Smith /*@C
560811af0c4SBarry Smith    PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()`
561793721a6SBarry Smith 
562793721a6SBarry Smith    Not Collective
563793721a6SBarry Smith 
5642fe279fdSBarry Smith    Output Parameter:
565793721a6SBarry Smith .  args - the command line arguments
566793721a6SBarry Smith 
567793721a6SBarry Smith    Level: intermediate
568793721a6SBarry Smith 
569db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()`
570793721a6SBarry Smith @*/
571d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args)
572d71ae5a4SJacob Faibussowitsch {
57339a651e2SJacob Faibussowitsch   PetscFunctionBegin;
57439a651e2SJacob Faibussowitsch   if (args) {
575793721a6SBarry Smith     PetscInt i = 0;
576793721a6SBarry Smith 
5779566063dSJacob Faibussowitsch     while (args[i]) PetscCall(PetscFree(args[i++]));
5789566063dSJacob Faibussowitsch     PetscCall(PetscFree(args));
57939a651e2SJacob Faibussowitsch   }
5803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
581793721a6SBarry Smith }
582793721a6SBarry Smith 
58327104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
58430befbd2SBarry Smith   #include <petscconfiginfo.h>
58530befbd2SBarry Smith 
586d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
587d71ae5a4SJacob Faibussowitsch {
58827104ee2SJacob Faibussowitsch   PetscFunctionBegin;
58911525c0dSBarry Smith   if (!PetscGlobalRank) {
59030befbd2SBarry Smith     char      cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64];
59111525c0dSBarry Smith     int       port;
592ffbd1cfbSBarry Smith     PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE;
59311525c0dSBarry Smith     size_t    applinelen, introlen;
594ffbd1cfbSBarry Smith     char      sawsurl[256];
59511525c0dSBarry Smith 
5969566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg));
59711525c0dSBarry Smith     if (flg) {
59811525c0dSBarry Smith       char sawslog[PETSC_MAX_PATH_LEN];
59911525c0dSBarry Smith 
6009566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL));
60111525c0dSBarry Smith       if (sawslog[0]) {
602792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog));
60311525c0dSBarry Smith       } else {
604792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL));
60511525c0dSBarry Smith       }
60611525c0dSBarry Smith     }
6079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg));
60848a46eb9SPierre Jolivet     if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert));
6099566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL));
610ffbd1cfbSBarry Smith     if (selectport) {
611792fecdfSBarry Smith       PetscCallSAWs(SAWs_Get_Available_Port, (&port));
612792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Port, (port));
613ffbd1cfbSBarry Smith     } else {
6149566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg));
61548a46eb9SPierre Jolivet       if (flg) PetscCallSAWs(SAWs_Set_Port, (port));
616ffbd1cfbSBarry Smith     }
6179566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg));
61811525c0dSBarry Smith     if (flg) {
619792fecdfSBarry Smith       PetscCallSAWs(SAWs_Set_Document_Root, (root));
6209566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(root, ".", &rootlocal));
6219c1e0ce8SBarry Smith     } else {
6229566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg));
6239c1e0ce8SBarry Smith       if (flg) {
6249566063dSJacob Faibussowitsch         PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root)));
625792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Document_Root, (root));
6269c1e0ce8SBarry Smith       }
62711525c0dSBarry Smith     }
6289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2));
62911525c0dSBarry Smith     if (flg2) {
63011525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
63128b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option");
6329566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root));
6339566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(jsdir, 'r', &flg));
63428b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory");
635792fecdfSBarry Smith       PetscCallSAWs(SAWs_Push_Local_Header, ());
63611525c0dSBarry Smith     }
6379566063dSJacob Faibussowitsch     PetscCall(PetscGetProgramName(programname, sizeof(programname)));
6389566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(help, &applinelen));
63911525c0dSBarry Smith     introlen = 4096 + applinelen;
64030a8c9c0SSurtai Han     applinelen += 1024;
6419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(applinelen, &appline));
6429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(introlen, &intro));
64311525c0dSBarry Smith 
64411525c0dSBarry Smith     if (rootlocal) {
6459566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname));
6469566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(appline, 'r', &rootlocal));
64711525c0dSBarry Smith     }
6489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetAll(NULL, &options));
64911525c0dSBarry Smith     if (rootlocal && help) {
6509566063dSJacob 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));
65111525c0dSBarry Smith     } else if (help) {
6529566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help));
65311525c0dSBarry Smith     } else {
6549566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options));
65511525c0dSBarry Smith     }
6569566063dSJacob Faibussowitsch     PetscCall(PetscFree(options));
6579566063dSJacob Faibussowitsch     PetscCall(PetscGetVersion(version, sizeof(version)));
6589371c9d4SSatish Balay     PetscCall(PetscSNPrintf(intro, introlen,
6599371c9d4SSatish Balay                             "<body>\n"
660a17b96a8SKyle 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"
661df62c222SBarry 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"
6629371c9d4SSatish Balay                             "%s",
6639371c9d4SSatish Balay                             version, petscconfigureoptions, appline));
664792fecdfSBarry Smith     PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro));
6659566063dSJacob Faibussowitsch     PetscCall(PetscFree(intro));
6669566063dSJacob Faibussowitsch     PetscCall(PetscFree(appline));
667ffbd1cfbSBarry Smith     if (selectport) {
668aa573868SBarry Smith       PetscBool silent;
6697d812c46SBarry Smith 
6707d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
67139a651e2SJacob Faibussowitsch       while (SAWs_Initialize()) {
672792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_Available_Port, (&port));
673792fecdfSBarry Smith         PetscCallSAWs(SAWs_Set_Port, (port));
6747d812c46SBarry Smith       }
6757d812c46SBarry Smith 
6769566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL));
677aa573868SBarry Smith       if (!silent) {
678792fecdfSBarry Smith         PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl));
6799566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl));
680ffbd1cfbSBarry Smith       }
6817d812c46SBarry Smith     } else {
682792fecdfSBarry Smith       PetscCallSAWs(SAWs_Initialize, ());
683aa573868SBarry Smith     }
6849566063dSJacob Faibussowitsch     PetscCall(PetscCitationsRegister("@TechReport{ saws,\n"
6850af79b04SBarry Smith                                      "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6860af79b04SBarry Smith                                      "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6870af79b04SBarry Smith                                      "  Institution = {Argonne National Laboratory},\n"
6889371c9d4SSatish Balay                                      "  Year   = 2013\n}\n",
6899371c9d4SSatish Balay                                      NULL));
69011525c0dSBarry Smith   }
6913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69211525c0dSBarry Smith }
69311525c0dSBarry Smith #endif
69411525c0dSBarry Smith 
6954dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
696d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
697d71ae5a4SJacob Faibussowitsch {
6984dfee713SSatish Balay   PetscFunctionBegin;
6994dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
7004dfee713SSatish Balay   /* see MPI.py for details on this bug */
7014dfee713SSatish Balay   (void)setenv("HWLOC_COMPONENTS", "-x86", 1);
7024dfee713SSatish Balay #endif
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7044dfee713SSatish Balay }
7054dfee713SSatish Balay 
706a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
707a56f64adSBarry Smith   #include <adios.h>
70822580e64SBarry Smith   #include <adios_read.h>
7097b56e58cSSatish Balay int64_t Petsc_adios_group;
710a56f64adSBarry Smith #endif
711a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
712cd1b99a6SBarry Smith   #include <omp.h>
713f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
714cd1b99a6SBarry Smith #endif
715a56f64adSBarry Smith 
716a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
717a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA)
7180e6b6b59SJacob Faibussowitsch   #include <petscdevice_cuda.h>
719a4af0ceeSJacob Faibussowitsch // REMOVE ME
720a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
721a4af0ceeSJacob Faibussowitsch #endif
722a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP)
7230e6b6b59SJacob Faibussowitsch   #include <petscdevice_hip.h>
724a4af0ceeSJacob Faibussowitsch // REMOVE ME
725a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
726a4af0ceeSJacob Faibussowitsch #endif
727a4af0ceeSJacob Faibussowitsch 
72827104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
729bc8a24ecSBarry Smith   #include <dlfcn.h>
730bc8a24ecSBarry Smith #endif
731a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
732a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
733a4af0ceeSJacob Faibussowitsch #endif
734a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
7353274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscViennaCLInit(void);
736a4af0ceeSJacob Faibussowitsch PetscBool                   PetscViennaCLSynchronize = PETSC_FALSE;
737a4af0ceeSJacob Faibussowitsch #endif
738bc8a24ecSBarry Smith 
739660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE;
740660278c0SBarry Smith 
74185649d77SJunchao Zhang /*
74285649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
74385649d77SJunchao Zhang 
74485649d77SJunchao Zhang   prog:     program name
74502101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
74685649d77SJunchao Zhang   help:     program help message
747da81f932SPierre Jolivet   ftn:      is it called from Fortran initialization (petscinitializef_)?
74885649d77SJunchao Zhang   readarguments,len: used when fortran is true
74985649d77SJunchao Zhang */
750d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len)
751d71ae5a4SJacob Faibussowitsch {
75285649d77SJunchao Zhang   PetscMPIInt size;
75385649d77SJunchao Zhang   PetscBool   flg = PETSC_TRUE;
75485649d77SJunchao Zhang   char        hostname[256];
75585649d77SJunchao Zhang 
75627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
7573ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
75839a651e2SJacob Faibussowitsch   /* these must be initialized in a routine, not as a constant declaration */
75939a651e2SJacob Faibussowitsch   PETSC_STDOUT = stdout;
76039a651e2SJacob Faibussowitsch   PETSC_STDERR = stderr;
76139a651e2SJacob Faibussowitsch 
7629566063dSJacob Faibussowitsch   /* PetscCall can be used from now */
76339a651e2SJacob Faibussowitsch   PetscErrorHandlingInitialized = PETSC_TRUE;
76439a651e2SJacob Faibussowitsch 
76585649d77SJunchao Zhang   /*
76685649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
76785649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
76885649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
76985649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
77085649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
77185649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
77285649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
77385649d77SJunchao Zhang       listed above and since that time are compatible.
77485649d77SJunchao Zhang 
77585649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
77685649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
77785649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
77885649d77SJunchao Zhang       to perform the checking.
77985649d77SJunchao Zhang 
78085649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
78185649d77SJunchao Zhang 
78285649d77SJunchao Zhang       Questions:
78385649d77SJunchao Zhang 
78485649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
78585649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
78685649d77SJunchao Zhang   */
78785649d77SJunchao Zhang 
78885649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
78985649d77SJunchao Zhang   {
79085649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
79185649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
79239a651e2SJacob Faibussowitsch 
7939566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength));
79485649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
79585649d77SJunchao Zhang   #if defined(MPICH_VERSION)
79685649d77SJunchao Zhang     #if MPICH_NUMVERSION < 30100000
79785649d77SJunchao Zhang     {
79885649d77SJunchao Zhang       char     *ver, *lf;
79985649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
80039a651e2SJacob Faibussowitsch 
8019566063dSJacob Faibussowitsch       PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver));
80239a651e2SJacob Faibussowitsch       if (ver) {
8039566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(ver, '\n', &lf));
80439a651e2SJacob Faibussowitsch         if (lf) {
80585649d77SJunchao Zhang           *lf = 0;
8069566063dSJacob Faibussowitsch           PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg));
80785649d77SJunchao Zhang         }
80885649d77SJunchao Zhang       }
80985649d77SJunchao Zhang       if (!flg) {
810d5b396e8SSatish Balay         PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION));
81185649d77SJunchao Zhang         flg = PETSC_TRUE;
81285649d77SJunchao Zhang       }
81385649d77SJunchao Zhang     }
81485649d77SJunchao Zhang     #endif
81585649d77SJunchao 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?) */
81685649d77SJunchao Zhang   #elif defined(OMPI_MAJOR_VERSION)
81785649d77SJunchao Zhang     {
81885649d77SJunchao Zhang       char     *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf;
81985649d77SJunchao Zhang       PetscBool flg                                              = PETSC_FALSE;
82085649d77SJunchao Zhang     #define PSTRSZ 2
82185649d77SJunchao Zhang       char      ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"};
82285649d77SJunchao Zhang       char      ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "};
82385649d77SJunchao Zhang       int       i;
82485649d77SJunchao Zhang       for (i = 0; i < PSTRSZ; i++) {
8259566063dSJacob Faibussowitsch         PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver));
82639a651e2SJacob Faibussowitsch         if (ver) {
8279566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION));
8289566063dSJacob Faibussowitsch           PetscCall(PetscStrstr(ver, bs, &bsf));
82939a651e2SJacob Faibussowitsch           if (bsf) flg = PETSC_TRUE;
83085649d77SJunchao Zhang           break;
83185649d77SJunchao Zhang         }
83285649d77SJunchao Zhang       }
83385649d77SJunchao Zhang       if (!flg) {
8343ba16761SJacob 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));
83585649d77SJunchao Zhang         flg = PETSC_TRUE;
83685649d77SJunchao Zhang       }
83785649d77SJunchao Zhang     }
83885649d77SJunchao Zhang   #endif
83985649d77SJunchao Zhang   }
84085649d77SJunchao Zhang #endif
84185649d77SJunchao Zhang 
842e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__))
84385649d77SJunchao 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 */
84439a651e2SJacob 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");
84585649d77SJunchao Zhang #endif
84685649d77SJunchao Zhang 
84785649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
84885649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
84985649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
85085649d77SJunchao Zhang #endif
85185649d77SJunchao Zhang 
8529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
85385649d77SJunchao Zhang 
85485649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
85585649d77SJunchao Zhang 
8569566063dSJacob Faibussowitsch   PetscCall(PetscSetProgramName(prog));
8579566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen));
8589566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout));
8599566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr));
8609566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockCreate(&PetscCommSpinLock));
86185649d77SJunchao Zhang 
86285649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
8639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN));
86485649d77SJunchao Zhang 
86585649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
8669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS));
8679566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE));
86885649d77SJunchao Zhang   }
86985649d77SJunchao Zhang 
87085649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
8719566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
87285649d77SJunchao Zhang 
8739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank));
8749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize));
87585649d77SJunchao Zhang 
87685649d77SJunchao Zhang   MPIU_BOOL        = MPI_INT;
87785649d77SJunchao Zhang   MPIU_ENUM        = MPI_INT;
87885649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64;
87985649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
88085649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
88185649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
88285649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
88385649d77SJunchao Zhang #endif
88439a651e2SJacob Faibussowitsch   else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t");
88585649d77SJunchao Zhang 
88685649d77SJunchao Zhang     /*
88785649d77SJunchao Zhang      Initialized the global complex variable; this is because with
88885649d77SJunchao Zhang      shared libraries the constructors for global variables
88985649d77SJunchao Zhang      are not called; at least on IRIX.
89085649d77SJunchao Zhang   */
89185649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
89285649d77SJunchao Zhang   {
89385649d77SJunchao Zhang   #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
89485649d77SJunchao Zhang     PetscComplex ic(0.0, 1.0);
89585649d77SJunchao Zhang     PETSC_i = ic;
89685649d77SJunchao Zhang   #else
89785649d77SJunchao Zhang     PETSC_i = _Complex_I;
89885649d77SJunchao Zhang   #endif
89985649d77SJunchao Zhang   }
90085649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
90185649d77SJunchao Zhang 
90285649d77SJunchao Zhang   /*
90385649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
90485649d77SJunchao Zhang      half of the entries and maxes the second half.
90585649d77SJunchao Zhang   */
9069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP));
90785649d77SJunchao Zhang 
908613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
9099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128));
9109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128));
911*74c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
9129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128));
9139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128));
91485649d77SJunchao Zhang   #endif
915*74c01ffaSSatish Balay #endif
9169e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
9179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16));
9189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU___FP16));
91985649d77SJunchao Zhang #endif
92085649d77SJunchao Zhang 
92185649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
9229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM));
923613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX));
924613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN));
9259e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
9269e517322SPierre Jolivet   PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128));
92785649d77SJunchao Zhang #endif
92885649d77SJunchao Zhang 
9299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR));
93062e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp));
9319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR));
93285649d77SJunchao Zhang 
93385649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
93485649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
93585649d77SJunchao Zhang   {
93685649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
93793d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)};
93885649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_REAL, MPIU_INT}, tmpStruct;
93985649d77SJunchao Zhang 
9409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
94193d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT));
9429566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9439566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT));
94485649d77SJunchao Zhang   }
94585649d77SJunchao Zhang   {
94685649d77SJunchao Zhang     PetscMPIInt  blockSizes[2]   = {1, 1};
94793d501b3SJacob Faibussowitsch     MPI_Aint     blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)};
94885649d77SJunchao Zhang     MPI_Datatype blockTypes[2]   = {MPIU_SCALAR, MPIU_INT}, tmpStruct;
94985649d77SJunchao Zhang 
9509566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct));
95193d501b3SJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT));
9529566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&tmpStruct));
9539566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT));
95485649d77SJunchao Zhang   }
95585649d77SJunchao Zhang #endif
95685649d77SJunchao Zhang 
95785649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
9589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT));
9599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_2INT));
96085649d77SJunchao Zhang #endif
9619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT));
9629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPI_4INT));
9639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT));
9649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&MPIU_4INT));
96585649d77SJunchao Zhang 
96685649d77SJunchao Zhang   /*
96785649d77SJunchao Zhang      Attributes to be set on PETSc communicators
96885649d77SJunchao Zhang   */
9699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0));
9709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0));
9719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0));
9729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0));
97362e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, (void *)0));
97462e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, (void *)0));
97585649d77SJunchao Zhang 
97685649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
9779566063dSJacob Faibussowitsch   if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len));
97885649d77SJunchao Zhang   else
97985649d77SJunchao Zhang #endif
9809566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file));
98185649d77SJunchao Zhang 
98285649d77SJunchao Zhang   /* call a second time so it can look in the options database */
9839566063dSJacob Faibussowitsch   PetscCall(PetscErrorPrintfInitialize());
98485649d77SJunchao Zhang 
98585649d77SJunchao Zhang   /*
98685649d77SJunchao Zhang      Check system options and print help
98785649d77SJunchao Zhang   */
9889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCheckInitial_Private(help));
98985649d77SJunchao Zhang 
990a4af0ceeSJacob Faibussowitsch   /*
9910e6b6b59SJacob Faibussowitsch     Creates the logging data structures; this is enabled even if logging is not turned on
9920e6b6b59SJacob Faibussowitsch     This is the last thing we do before returning to the user code to prevent having the
9930e6b6b59SJacob Faibussowitsch     logging numbers contaminated by any startup time associated with MPI
9940e6b6b59SJacob Faibussowitsch   */
9950e6b6b59SJacob Faibussowitsch #if defined(PETSC_USE_LOG)
9960e6b6b59SJacob Faibussowitsch   PetscCall(PetscLogInitialize());
9970e6b6b59SJacob Faibussowitsch #endif
9980e6b6b59SJacob Faibussowitsch 
9990e6b6b59SJacob Faibussowitsch   /*
1000a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
1001a4af0ceeSJacob Faibussowitsch 
1002a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
1003a4af0ceeSJacob Faibussowitsch    1. MPI initialized
1004a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
10050e6b6b59SJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up
10060e6b6b59SJacob Faibussowitsch       its own SIGSEV handler via the push/pop interface.
10070e6b6b59SJacob Faibussowitsch    4. Logging initialized
1008a4af0ceeSJacob Faibussowitsch   */
10099566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD));
1010a4af0ceeSJacob Faibussowitsch 
1011a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
1012a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
10139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_summary", &flg));
10149566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg));
10159566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL));
1016a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
10179566063dSJacob Faibussowitsch   PetscCall(PetscViennaCLInit());
1018a4af0ceeSJacob Faibussowitsch #endif
1019a4af0ceeSJacob Faibussowitsch 
10209566063dSJacob Faibussowitsch   PetscCall(PetscCitationsInitialize());
102185649d77SJunchao Zhang 
102285649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
10239566063dSJacob Faibussowitsch   PetscCall(PetscInitializeSAWs(ftn ? NULL : help));
102427104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
10259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg));
10269566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscStackViewSAWs());
102785649d77SJunchao Zhang #endif
102885649d77SJunchao Zhang 
102985649d77SJunchao Zhang   /*
103085649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
103185649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
103285649d77SJunchao Zhang   */
10339566063dSJacob Faibussowitsch   PetscCall(PetscInitialize_DynamicLibraries());
103485649d77SJunchao Zhang 
10359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10369566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size));
103796dcfd6cSBarry Smith   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
10389566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname));
103985649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
104085649d77SJunchao Zhang   {
104185649d77SJunchao Zhang     PetscBool omp_view_flag;
104285649d77SJunchao Zhang     char     *threads = getenv("OMP_NUM_THREADS");
104385649d77SJunchao Zhang 
104485649d77SJunchao Zhang     if (threads) {
10459566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads));
104685649d77SJunchao Zhang       (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads);
104785649d77SJunchao Zhang     } else {
10482f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt)omp_get_max_threads();
10499566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads));
105085649d77SJunchao Zhang     }
1051d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys");
10529566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg));
10539566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag));
1054d0609cedSBarry Smith     PetscOptionsEnd();
105585649d77SJunchao Zhang     if (flg) {
10569566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads));
105785649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
105885649d77SJunchao Zhang     }
105948a46eb9SPierre Jolivet     if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads));
106085649d77SJunchao Zhang   }
106185649d77SJunchao Zhang #endif
106285649d77SJunchao Zhang 
106385649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
106485649d77SJunchao Zhang   /*
106585649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
106685649d77SJunchao Zhang 
106785649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
106885649d77SJunchao Zhang   */
10699566063dSJacob Faibussowitsch   if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL));
107085649d77SJunchao Zhang #endif
107185649d77SJunchao Zhang 
107285649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
10739566063dSJacob Faibussowitsch   PetscCall(PetscFPTCreate(10000));
107485649d77SJunchao Zhang #endif
107585649d77SJunchao Zhang 
107685649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
107785649d77SJunchao Zhang   {
107885649d77SJunchao Zhang     PetscViewer viewer;
10799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg));
108085649d77SJunchao Zhang     if (flg) {
10819566063dSJacob Faibussowitsch       PetscCall(PetscProcessPlacementView(viewer));
10829566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
108385649d77SJunchao Zhang     }
108485649d77SJunchao Zhang   }
108585649d77SJunchao Zhang #endif
108685649d77SJunchao Zhang 
108785649d77SJunchao Zhang   flg = PETSC_TRUE;
10889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL));
10899566063dSJacob Faibussowitsch   if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE));
109085649d77SJunchao Zhang 
109185649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
10923ba16761SJacob Faibussowitsch   PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD);
10933ba16761SJacob Faibussowitsch   PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default);
10943ba16761SJacob Faibussowitsch   PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", "");
10953ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, "");
109685649d77SJunchao Zhang #endif
109785649d77SJunchao Zhang 
109885649d77SJunchao Zhang #if defined(__VALGRIND_H)
109985649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE;
110085649d77SJunchao Zhang   #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
11019566063dSJacob 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"));
110285649d77SJunchao Zhang   #endif
110385649d77SJunchao Zhang #endif
110485649d77SJunchao Zhang   /*
110585649d77SJunchao Zhang       Set flag that we are completely initialized
110685649d77SJunchao Zhang   */
110785649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
110885649d77SJunchao Zhang 
11099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg));
11109566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscPythonInitialize(NULL, NULL));
1111f1f2ae84SBarry Smith 
1112f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1113f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin());
1114f1f2ae84SBarry 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");
11153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111685649d77SJunchao Zhang }
111785649d77SJunchao Zhang 
1118e5c89e4eSSatish Balay /*@C
1119e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1120811af0c4SBarry Smith    `PetscInitialize()` calls MPI_Init() if that has yet to be called,
1121e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1122e5c89e4eSSatish Balay    your program -- usually the very first line!
1123e5c89e4eSSatish Balay 
1124811af0c4SBarry Smith    Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set
1125e5c89e4eSSatish Balay 
1126e5c89e4eSSatish Balay    Input Parameters:
1127e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1128e5c89e4eSSatish Balay .  args - the command line arguments
1129be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1130be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1131be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1132c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
11330298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1134e5c89e4eSSatish Balay 
1135811af0c4SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that
1136811af0c4SBarry Smith    communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a
1137811af0c4SBarry Smith    four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not,
1138811af0c4SBarry Smith    then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even
113905827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1140e5c89e4eSSatish Balay 
1141e5c89e4eSSatish Balay    Options Database Keys:
11427ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
11437ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1144e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11458a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
1146811af0c4SBarry Smith .  -on_error_abort - calls `abort()` when error detected (no traceback)
1147811af0c4SBarry Smith .  -on_error_mpiabort - calls `MPI_abort()` when error detected
11481af3601dSBarry Smith .  -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr
11498a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1150bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1151e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1152e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1153e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
115479dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
115579dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
1156811af0c4SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see `PetscMallocSetDebug()`
1157aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
115892f119d6SBarry 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
1159811af0c4SBarry Smith .  -malloc_view - show a list of all allocated memory during `PetscFinalize()`
116092f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1161608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
116292f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1163e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1164e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1165e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1166e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1167e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
11680841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1169e5c89e4eSSatish Balay 
1170c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1171c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1172c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1173811af0c4SBarry Smith -  -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()`
1174c5b5d8d5SVaclav Hapla 
1175c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1176c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1177c5b5d8d5SVaclav Hapla    They can be used also in option files.
1178c5b5d8d5SVaclav Hapla 
1179811af0c4SBarry Smith    See `PetscOptionsMonitorSet()` to do monitoring programmatically.
1180c5b5d8d5SVaclav Hapla 
1181e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1182a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1183811af0c4SBarry Smith +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`.
1184217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1185217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1186495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1187811af0c4SBarry Smith         hangs without running in the debugger).  See `PetscLogTraceBegin()`.
1188811af0c4SBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see `PetscLogView()`.
1189811af0c4SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`.
1190811af0c4SBarry Smith .  -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView().
11919a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1192495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1193f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1194811af0c4SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See `PetscLogDump()`.
1195811af0c4SBarry Smith .  -log [filename] - Logs basic profiline information  See `PetscLogDump()`.
1196c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
1197811af0c4SBarry Smith .  -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off
1198c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1199495fc317SBarry Smith 
12006a77f485SPierre Jolivet     Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time
1201e5c89e4eSSatish Balay 
1202ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1203ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1204ffbd1cfbSBarry 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
1205ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1206ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1207ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1208ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1209ffbd1cfbSBarry Smith 
1210e5c89e4eSSatish Balay    Environmental Variables:
1211811af0c4SBarry Smith +   `PETSC_TMP` - alternative tmp directory
1212811af0c4SBarry Smith .   `PETSC_SHARED_TMP` - tmp is shared by all processes
1213811af0c4SBarry Smith .   `PETSC_NOT_SHARED_TMP` - each process has its own private tmp
1214811af0c4SBarry Smith .   `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs
1215811af0c4SBarry Smith .   `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1216811af0c4SBarry Smith .   `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer
1217811af0c4SBarry Smith -   `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to
1218e5c89e4eSSatish Balay 
1219e5c89e4eSSatish Balay    Level: beginner
1220e5c89e4eSSatish Balay 
1221811af0c4SBarry Smith    Note:
1222811af0c4SBarry Smith    If for some reason you must call `MPI_Init()` separately, call
1223811af0c4SBarry Smith    it before `PetscInitialize()`.
1224e5c89e4eSSatish Balay 
1225811af0c4SBarry Smith    Fortran Notes:
1226811af0c4SBarry Smith    In Fortran this routine can be called with
1227811af0c4SBarry Smith .vb
1228811af0c4SBarry Smith        call PetscInitialize(ierr)
1229811af0c4SBarry Smith        call PetscInitialize(file,ierr) or
1230811af0c4SBarry Smith        call PetscInitialize(file,help,ierr)
1231811af0c4SBarry Smith .ve
1232e5c89e4eSSatish Balay 
1233811af0c4SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after
1234811af0c4SBarry Smith    calling `PetscInitialize()`.
1235e5c89e4eSSatish Balay 
12365dedd997SBarry Smith    Options Database Key for Developers:
12375dedd997SBarry Smith .  -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form:
12385dedd997SBarry Smith     "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py)
12395dedd997SBarry Smith 
1240db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()`
1241e5c89e4eSSatish Balay @*/
1242d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[])
1243d71ae5a4SJacob Faibussowitsch {
124485649d77SJunchao Zhang   PetscMPIInt flag;
124521ef0414SBarry Smith   const char *prog = "Unknown Name", *mpienv;
1246e5c89e4eSSatish Balay 
124727104ee2SJacob Faibussowitsch   PetscFunctionBegin;
12483ba16761SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
12499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Initialized(&flag));
1250e5c89e4eSSatish Balay   if (!flag) {
125139a651e2SJacob 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");
12529566063dSJacob Faibussowitsch     PetscCall(PetscPreMPIInit_Private());
12535e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12545e765c61SJed Brown     {
125539a651e2SJacob Faibussowitsch       PetscMPIInt PETSC_UNUSED provided;
12569566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED, &provided));
12575e765c61SJed Brown     }
12585e765c61SJed Brown #else
12599566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Init(argc, args));
12605e765c61SJed Brown #endif
126121ef0414SBarry Smith     if (PetscDefined(HAVE_MPIUNI)) {
126221ef0414SBarry Smith       mpienv = getenv("PMI_SIZE");
126321ef0414SBarry Smith       if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE");
126421ef0414SBarry Smith       if (mpienv) {
126521ef0414SBarry Smith         PetscInt isize;
126621ef0414SBarry Smith         PetscCall(PetscOptionsStringToInt(mpienv, &isize));
126721ef0414SBarry 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");
126821ef0414SBarry 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");
126921ef0414SBarry Smith       }
127021ef0414SBarry Smith     }
1271e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1272e5c89e4eSSatish Balay   }
12734dfee713SSatish Balay 
127485649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1275e5c89e4eSSatish Balay   if (argc && args) {
1276e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1277e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1278e5c89e4eSSatish Balay   }
1279811af0c4SBarry Smith   PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0));
12803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1281e5c89e4eSSatish Balay }
1282e5c89e4eSSatish Balay 
1283a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
128495c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1285ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsCounts;
1286ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt     PetscObjectsMaxCounts;
128795c0884eSLisandro Dalcin PETSC_INTERN PetscBool    PetscObjectsLog;
12884097062eSBarry Smith #endif
1289e5c89e4eSSatish Balay 
1290008a6e76SBarry Smith /*
1291008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1292008a6e76SBarry Smith */
1293d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void)
1294d71ae5a4SJacob Faibussowitsch {
1295008a6e76SBarry Smith   PetscFunctionBegin;
1296613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128)
12979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128));
1298*74c01ffaSSatish Balay   #if defined(PETSC_HAVE_COMPLEX)
12999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128));
1300008a6e76SBarry Smith   #endif
1301*74c01ffaSSatish Balay #endif
13029e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16)
13039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU___FP16));
1304008a6e76SBarry Smith #endif
1305008a6e76SBarry Smith 
1306de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
13079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_SUM));
1308613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MAX));
1309613bf2b2SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_MIN));
13109e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
13119e517322SPierre Jolivet   PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128));
1312008a6e76SBarry Smith #endif
1313008a6e76SBarry Smith 
13149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR));
13159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT));
13169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT));
131740df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
13189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_2INT));
1319008a6e76SBarry Smith #endif
13209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPI_4INT));
13219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&MPIU_4INT));
13229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP));
132362e5d2d2SJDBetteridge   PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp));
13243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1325008a6e76SBarry Smith }
1326008a6e76SBarry Smith 
1327a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
1328a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1329a4af0ceeSJacob Faibussowitsch #endif
1330a4af0ceeSJacob Faibussowitsch 
1331e5c89e4eSSatish Balay /*@C
1332e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1333811af0c4SBarry Smith    of the program. `MPI_Finalize()` is called only if the user had not
1334811af0c4SBarry Smith    called `MPI_Init()` before calling `PetscInitialize()`.
1335e5c89e4eSSatish Balay 
1336811af0c4SBarry Smith    Collective on `PETSC_COMM_WORLD`
1337e5c89e4eSSatish Balay 
1338e5c89e4eSSatish Balay    Options Database Keys:
1339811af0c4SBarry Smith +  -options_view - Calls `PetscOptionsView()`
1340e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
13417eb1d149SBarry 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
1342e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1343811af0c4SBarry Smith .  -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed
1344e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
134592f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1346e5c89e4eSSatish Balay 
1347e5c89e4eSSatish Balay    Level: beginner
1348e5c89e4eSSatish Balay 
1349e5c89e4eSSatish Balay    Note:
1350811af0c4SBarry Smith    See `PetscInitialize()` for other runtime options.
1351e5c89e4eSSatish Balay 
1352db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()`
1353e5c89e4eSSatish Balay @*/
1354d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void)
1355d71ae5a4SJacob Faibussowitsch {
13564bb5149bSJed Brown   PetscMPIInt rank;
1357a8d2bbe5SBarry Smith   PetscInt    nopt;
13582bf49c77SBarry Smith   PetscBool   flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE;
1359dff31646SBarry Smith   PetscBool   flg;
136010463e74SBarry Smith #if defined(PETSC_USE_LOG)
136110463e74SBarry Smith   char mname[PETSC_MAX_PATH_LEN];
136210463e74SBarry Smith #endif
1363e5c89e4eSSatish Balay 
13643cbbc5ffSBarry Smith   PetscFunctionBegin;
136539a651e2SJacob Faibussowitsch   PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()");
13669566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "PetscFinalize() called\n"));
1367b022a5c1SBarry Smith 
1368f1f2ae84SBarry Smith   PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg));
1369f1f2ae84SBarry Smith   if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd());
1370f1f2ae84SBarry Smith 
137162e5d2d2SJDBetteridge   /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */
137262e5d2d2SJDBetteridge   {
137362e5d2d2SJDBetteridge     union
137462e5d2d2SJDBetteridge     {
137562e5d2d2SJDBetteridge       MPI_Comm comm;
137662e5d2d2SJDBetteridge       void    *ptr;
137762e5d2d2SJDBetteridge     } ucomm;
137862e5d2d2SJDBetteridge     PetscMPIInt flg;
137962e5d2d2SJDBetteridge     void       *tmp;
138062e5d2d2SJDBetteridge 
138162e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
138262e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
138362e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF));
138462e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
138562e5d2d2SJDBetteridge     if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg));
138662e5d2d2SJDBetteridge     if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD));
138762e5d2d2SJDBetteridge   }
138862e5d2d2SJDBetteridge 
13899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1390a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
13913ba16761SJacob Faibussowitsch   PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE);
13923ba16761SJacob Faibussowitsch   PetscCallExternal(adios_finalize, rank);
1393a56f64adSBarry Smith #endif
13949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg));
1395dff31646SBarry Smith   if (flg) {
13961f817a21SBarry Smith     char *cits, filename[PETSC_MAX_PATH_LEN];
13971f817a21SBarry Smith     FILE *fd = PETSC_STDOUT;
13981f817a21SBarry Smith 
13999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL));
140048a46eb9SPierre Jolivet     if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd));
14019566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits));
1402dff31646SBarry Smith     cits[0] = 0;
14039566063dSJacob Faibussowitsch     PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits));
14049566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n"));
14059566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14069566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits));
14079566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n"));
14089566063dSJacob Faibussowitsch     PetscCall(PetscFClose(PETSC_COMM_WORLD, fd));
14099566063dSJacob Faibussowitsch     PetscCall(PetscFree(cits));
1410dff31646SBarry Smith   }
14119566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&PetscCitationsList));
1412dff31646SBarry Smith 
1413c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
141404102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
141504102261SBarry Smith   {
141604102261SBarry Smith     PetscInt nmax = 2;
141704102261SBarry Smith     char   **buffs;
14189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
14199566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-textbelt", buffs, &nmax, &flg1));
142004102261SBarry Smith     if (flg1) {
142128b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-textbelt requires either the phone number or number,\"message\"");
142204102261SBarry Smith       if (nmax == 1) {
1423c6a7a370SJeremy L Thompson         size_t len = 128;
1424c6a7a370SJeremy L Thompson         PetscCall(PetscMalloc1(len, &buffs[1]));
14259566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
1426c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(buffs[1], " has completed", len));
142704102261SBarry Smith       }
14289566063dSJacob Faibussowitsch       PetscCall(PetscTextBelt(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
14299566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
14309566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
143104102261SBarry Smith     }
14329566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
143304102261SBarry Smith   }
1434763ec7b1SBarry Smith   {
1435763ec7b1SBarry Smith     PetscInt nmax = 2;
1436763ec7b1SBarry Smith     char   **buffs;
14379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2, &buffs));
14389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-tellmycell", buffs, &nmax, &flg1));
1439763ec7b1SBarry Smith     if (flg1) {
144028b400f6SJacob Faibussowitsch       PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-tellmycell requires either the phone number or number,\"message\"");
1441763ec7b1SBarry Smith       if (nmax == 1) {
1442c6a7a370SJeremy L Thompson         size_t len = 128;
1443c6a7a370SJeremy L Thompson         PetscCall(PetscMalloc1(len, &buffs[1]));
14449566063dSJacob Faibussowitsch         PetscCall(PetscGetProgramName(buffs[1], 32));
1445c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(buffs[1], " has completed", len));
1446763ec7b1SBarry Smith       }
14479566063dSJacob Faibussowitsch       PetscCall(PetscTellMyCell(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL));
14489566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[0]));
14499566063dSJacob Faibussowitsch       PetscCall(PetscFree(buffs[1]));
1450763ec7b1SBarry Smith     }
14519566063dSJacob Faibussowitsch     PetscCall(PetscFree(buffs));
1452763ec7b1SBarry Smith   }
145304102261SBarry Smith #endif
145404102261SBarry Smith 
14552d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14569566063dSJacob Faibussowitsch   PetscCall(PetscFPTDestroy());
14572d53ad75SBarry Smith #endif
14582d53ad75SBarry Smith 
1459e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1460dff31646SBarry Smith   flg = PETSC_FALSE;
14619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL));
14621baa6e33SBarry Smith   if (flg) PetscCall(PetscOptionsSAWsDestroy());
1463d5649816SBarry Smith #endif
1464d5649816SBarry Smith 
1465681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1466681455b2SBarry Smith   flg1 = PETSC_FALSE;
14679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL));
1468681455b2SBarry Smith   if (flg1) {
1469681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
14709566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL));
1471681455b2SBarry Smith   }
1472681455b2SBarry Smith #endif
1473681455b2SBarry Smith 
147467584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_info", &flg2, NULL));
1476e5c89e4eSSatish Balay   if (!flg2) {
147790d69ab7SBarry Smith     flg2 = PETSC_FALSE;
14789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL));
1479e5c89e4eSSatish Balay   }
148048a46eb9SPierre Jolivet   if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n"));
148167584ceeSBarry Smith #endif
1482e5c89e4eSSatish Balay 
1483e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
148490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
14859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL));
1486e5c89e4eSSatish Balay   if (flg1) {
1487e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
14889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD));
14899566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops));
1490e5c89e4eSSatish Balay   }
1491e5c89e4eSSatish Balay #endif
1492e5c89e4eSSatish Balay 
1493e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1494e5c89e4eSSatish Balay   #if defined(PETSC_HAVE_MPE)
1495e5c89e4eSSatish Balay   mname[0] = 0;
14969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1));
1497e5c89e4eSSatish Balay   if (flg1) {
14989566063dSJacob Faibussowitsch     if (mname[0]) PetscCall(PetscLogMPEDump(mname));
14999566063dSJacob Faibussowitsch     else PetscCall(PetscLogMPEDump(0));
1500e5c89e4eSSatish Balay   }
1501e5c89e4eSSatish Balay   #endif
1502356e5837SBarry Smith #endif
1503a297a907SKarl Rupp 
1504dd710f27SBarry Smith   /*
1505dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1506dd710f27SBarry Smith   */
15079566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1508dd710f27SBarry Smith 
1509356e5837SBarry Smith #if defined(PETSC_USE_LOG)
15109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE));
15119566063dSJacob Faibussowitsch   PetscCall(PetscLogViewFromOptions());
15129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsPopGetViewerOff());
151387c3beb6SLisandro Dalcin 
1514356e5837SBarry Smith   mname[0] = 0;
15159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_summary", mname, sizeof(mname), &flg1));
1516e5c89e4eSSatish Balay   if (flg1) {
151791eabc43SBarry Smith     PetscViewer viewer;
15189566063dSJacob Faibussowitsch     PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD, "\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n"));
151991eabc43SBarry Smith     if (mname[0]) {
15209566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, mname, &viewer));
15219566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
15229566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
152333f85c2fSBarry Smith     } else {
152433f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
15259566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
15269566063dSJacob Faibussowitsch       PetscCall(PetscLogView(viewer));
15279566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
152833f85c2fSBarry Smith     }
1529e5c89e4eSSatish Balay   }
1530a297a907SKarl Rupp 
1531dd710f27SBarry Smith   /*
1532dd710f27SBarry Smith      Free any objects created by the last block of code.
1533dd710f27SBarry Smith   */
15349566063dSJacob Faibussowitsch   PetscCall(PetscObjectRegisterDestroyAll());
1535dd710f27SBarry Smith 
1536dd710f27SBarry Smith   mname[0] = 0;
15379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1));
15389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2));
15399566063dSJacob Faibussowitsch   if (flg1 || flg2) PetscCall(PetscLogDump(mname));
1540e5c89e4eSSatish Balay #endif
154110463e74SBarry Smith 
154290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL));
15449566063dSJacob Faibussowitsch   if (!flg1) PetscCall(PetscPopSignalHandler());
154590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
15469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL));
15471baa6e33SBarry Smith   if (flg1) PetscCall(PetscMPIDump(stdout));
154890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
154990d69ab7SBarry Smith   flg2 = PETSC_FALSE;
15508bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
15519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1));
15529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
15539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL));
1554e4c476e2SSatish Balay 
1555e5c89e4eSSatish Balay   if (flg2) {
1556be56827dSJed Brown     PetscViewer viewer;
15579566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15589566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsView(NULL, viewer));
15609566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1561e5c89e4eSSatish Balay   }
1562e5c89e4eSSatish Balay 
1563e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
15649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1));
15659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1));
1566e5c89e4eSSatish Balay 
156733fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
15689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1));
15699245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1570e5c89e4eSSatish Balay   if (flg3) {
15713de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1572be56827dSJed Brown       PetscViewer viewer;
15739566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
15749566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
15759566063dSJacob Faibussowitsch       PetscCall(PetscOptionsView(NULL, viewer));
15769566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
1577e5c89e4eSSatish Balay     }
15789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsAllUsed(NULL, &nopt));
15793de2bfdfSBarry Smith     if (nopt) {
15809566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n"));
15819566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n"));
15823de2bfdfSBarry Smith       if (nopt == 1) {
15839566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n"));
1584e5c89e4eSSatish Balay       } else {
15859566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt));
1586e5c89e4eSSatish Balay       }
15873de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15889566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n"));
1589df12ba86SBarry Smith     }
15909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsLeft(NULL));
1591e5c89e4eSSatish Balay   }
1592e5c89e4eSSatish Balay 
1593e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1594d45a07a7SBarry Smith   if (!PetscGlobalRank) {
15959566063dSJacob Faibussowitsch     PetscCall(PetscStackSAWsViewOff());
1596792fecdfSBarry Smith     PetscCallSAWs(SAWs_Finalize, ());
1597d45a07a7SBarry Smith   }
1598ec957eceSBarry Smith #endif
1599ec957eceSBarry Smith 
16004097062eSBarry Smith #if defined(PETSC_USE_LOG)
160110463e74SBarry Smith   /*
1602dbc8283eSBarry Smith        List all objects the user may have forgot to free
16032eff7a51SBarry Smith   */
160405df10baSBarry Smith   if (PetscObjectsLog) {
16059566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1));
1606a64a8e02SBarry Smith     if (flg1) {
1607a64a8e02SBarry Smith       MPI_Comm local_comm;
16087eb1d149SBarry Smith       char     string[64];
1609a64a8e02SBarry Smith 
16109566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL));
1611afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16129566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16139566063dSJacob Faibussowitsch       PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE));
16149566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16159566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
16160a1571b3SBarry Smith     }
161705df10baSBarry Smith   }
16184097062eSBarry Smith #endif
16194097062eSBarry Smith 
16204097062eSBarry Smith #if defined(PETSC_USE_LOG)
1621dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1622dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
16239566063dSJacob Faibussowitsch   PetscCall(PetscFree(PetscObjects));
16244097062eSBarry Smith #endif
16252eff7a51SBarry Smith 
162633f85c2fSBarry Smith   /*
162733f85c2fSBarry Smith      Destroy any packages that registered a finalize
162833f85c2fSBarry Smith   */
16299566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalizeAll());
163033f85c2fSBarry Smith 
1631101409b8SToby Isaac #if defined(PETSC_USE_LOG)
16329566063dSJacob Faibussowitsch   PetscCall(PetscLogFinalize());
1633101409b8SToby Isaac #endif
1634101409b8SToby Isaac 
163533f85c2fSBarry Smith   /*
163648dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
163748dd1dffSBarry Smith   */
16382e956fe4SStefano Zampini   if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll());
163937e93019SBarry Smith 
16404028d114SSatish Balay   if (petsc_history) {
16419566063dSJacob Faibussowitsch     PetscCall(PetscCloseHistoryFile(&petsc_history));
164202c9f0b5SLisandro Dalcin     petsc_history = NULL;
1643e5c89e4eSSatish Balay   }
16449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton));
16459566063dSJacob Faibussowitsch   PetscCall(PetscInfoDestroy());
1646e5c89e4eSSatish Balay 
164767584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
164892f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1649e5c89e4eSSatish Balay     char  fname[PETSC_MAX_PATH_LEN];
165092f119d6SBarry Smith     char  sname[PETSC_MAX_PATH_LEN];
1651e5c89e4eSSatish Balay     FILE *fd;
1652ed9cf6e9SBarry Smith     int   err;
1653e5c89e4eSSatish Balay 
1654dc92acbaSJed Brown     flg2 = PETSC_FALSE;
165592f119d6SBarry Smith     flg3 = PETSC_FALSE;
16569566063dSJacob Faibussowitsch     if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL));
16579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL));
165892f119d6SBarry Smith     fname[0] = 0;
16599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1));
1660e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
16613ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
16629371c9d4SSatish Balay       fd = fopen(sname, "w");
16639371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16649566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(fd));
1665ed9cf6e9SBarry Smith       err = fclose(fd);
166628b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
166792f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1668e5c89e4eSSatish Balay       MPI_Comm local_comm;
1669e5c89e4eSSatish Balay 
1670afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16719566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16729566063dSJacob Faibussowitsch       PetscCall(PetscMallocDump(stdout));
16739566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16749566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1675e5c89e4eSSatish Balay     }
1676e5c89e4eSSatish Balay     fname[0] = 0;
16779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1));
1678e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
16793ba16761SJacob Faibussowitsch       PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank));
16809371c9d4SSatish Balay       fd = fopen(sname, "w");
16819371c9d4SSatish Balay       PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname);
16829566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(fd));
1683ed9cf6e9SBarry Smith       err = fclose(fd);
168428b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
168592f119d6SBarry Smith     } else if (flg1) {
168692f119d6SBarry Smith       MPI_Comm local_comm;
168792f119d6SBarry Smith 
1688afd7ed4bSMatthew G. Knepley       PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm));
16899566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1));
16909566063dSJacob Faibussowitsch       PetscCall(PetscMallocView(stdout));
16919566063dSJacob Faibussowitsch       PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1));
16929566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&local_comm));
1693e5c89e4eSSatish Balay     }
1694e5c89e4eSSatish Balay   }
169567584ceeSBarry Smith #endif
169620e2c332SMatthew G. Knepley 
16975486ca60SMatthew G. Knepley   /*
16985486ca60SMatthew G. Knepley      Close any open dynamic libraries
16995486ca60SMatthew G. Knepley   */
17009566063dSJacob Faibussowitsch   PetscCall(PetscFinalize_DynamicLibraries());
17015486ca60SMatthew G. Knepley 
1702e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
17039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroyDefault());
1704e5c89e4eSSatish Balay 
1705e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
170602c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1707e5c89e4eSSatish Balay 
1708c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1709c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
17109566063dSJacob Faibussowitsch     PetscCall(PetscKokkosFinalize_Private());
1711c2b86a48SJunchao Zhang     PetscBeganKokkos       = PETSC_FALSE;
171245639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1713c2b86a48SJunchao Zhang   }
1714c2b86a48SJunchao Zhang #endif
1715c2b86a48SJunchao Zhang 
171671438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
171771438e86SJunchao Zhang   if (PetscBeganNvshmem) {
17189566063dSJacob Faibussowitsch     PetscCall(PetscNvshmemFinalize());
171971438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
172071438e86SJunchao Zhang   }
172171438e86SJunchao Zhang #endif
1722a0e72f99SJunchao Zhang 
17239566063dSJacob Faibussowitsch   PetscCall(PetscFreeMPIResources());
1724e5c89e4eSSatish Balay 
1725dbc8283eSBarry Smith   /*
1726efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1727efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1728efb80d3cSBarry Smith 
1729efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1730efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1731dbc8283eSBarry Smith  */
1732b770b1f6SSatish Balay   {
1733dbc8283eSBarry Smith     PetscCommCounter *counter;
1734dbc8283eSBarry Smith     PetscMPIInt       flg;
1735dbc8283eSBarry Smith     MPI_Comm          icomm;
17369371c9d4SSatish Balay     union
17379371c9d4SSatish Balay     {
17389371c9d4SSatish Balay       MPI_Comm comm;
17399371c9d4SSatish Balay       void    *ptr;
17409371c9d4SSatish Balay     } ucomm;
17419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg));
1742dbc8283eSBarry Smith     if (flg) {
1743265f3f35SJed Brown       icomm = ucomm.comm;
17449566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
174528b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1746dbc8283eSBarry Smith 
17479566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval));
17489566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17499566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1750dbc8283eSBarry Smith     }
17519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg));
1752dbc8283eSBarry Smith     if (flg) {
1753265f3f35SJed Brown       icomm = ucomm.comm;
17549566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg));
175528b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1756dbc8283eSBarry Smith 
17579566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval));
17589566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval));
17599566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_free(&icomm));
1760dbc8283eSBarry Smith     }
1761b770b1f6SSatish Balay   }
1762dbc8283eSBarry Smith 
17639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval));
17649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval));
17659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval));
17669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval));
176762e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval));
176862e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval));
1769480cf27aSJed Brown 
17705ea2b939SDuncan Campbell   // Free keyvals which may be silently created by some routines
17715ea2b939SDuncan Campbell   if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval));
17725ea2b939SDuncan Campbell   if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval));
17735ea2b939SDuncan Campbell 
17749566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen));
17759566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout));
17769566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr));
17779566063dSJacob Faibussowitsch   PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock));
1778ef19f930SBarry Smith 
1779e5c89e4eSSatish Balay   if (PetscBeganMPI) {
178099b1327fSBarry Smith     PetscMPIInt flag;
17819566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalized(&flag));
178239a651e2SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
178339a651e2SJacob Faibussowitsch     /* wait until the very last moment to disable error handling */
178439a651e2SJacob Faibussowitsch     PetscErrorHandlingInitialized = PETSC_FALSE;
17859566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Finalize());
178639a651e2SJacob Faibussowitsch   } else PetscErrorHandlingInitialized = PETSC_FALSE;
178739a651e2SJacob Faibussowitsch 
1788e5c89e4eSSatish Balay   /*
1789e5c89e4eSSatish Balay 
1790e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1791e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1792e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1793e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1794e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17950e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1796e5c89e4eSSatish Balay    memory was not freed.
1797e5c89e4eSSatish Balay 
1798e5c89e4eSSatish Balay */
17999566063dSJacob Faibussowitsch   PetscCall(PetscMallocClear());
18009566063dSJacob Faibussowitsch   PetscCall(PetscStackReset());
1801a297a907SKarl Rupp 
1802e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1803e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
18047ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE)
180556883f60SBarry Smith   /*
180656883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
180756883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
180856883f60SBarry Smith    */
180956883f60SBarry Smith   __gcov_flush();
181056883f60SBarry Smith #endif
18111724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
18121724198aSStefano Zampini   PetscStackClearTop;
18133ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1814e5c89e4eSSatish Balay }
1815e5c89e4eSSatish Balay 
181643db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
1817d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b)
1818d71ae5a4SJacob Faibussowitsch {
181943db4dbbSBarry Smith   if (*a == *b) return 1;
182043db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
182143db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
182243db4dbbSBarry Smith   return 0;
182343db4dbbSBarry Smith }
1824a70650f6SBarry Smith #endif
182543db4dbbSBarry Smith 
182643db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
1827d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b)
1828d71ae5a4SJacob Faibussowitsch {
182943db4dbbSBarry Smith   if (*a == *b) return 1;
183043db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
183143db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
183243db4dbbSBarry Smith   return 0;
183343db4dbbSBarry Smith }
183443db4dbbSBarry Smith #endif
1835