xref: /petsc/src/sys/objects/pinit.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>        /*I  "petscsys.h"   I*/
6665c2dedSJed Brown #include <petscviewer.h>
7a0e72f99SJunchao Zhang 
86edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
96edef35eSSatish Balay #include <petsc/private/valgrind/valgrind.h>
106edef35eSSatish Balay #endif
116edef35eSSatish Balay 
1285649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
1385649d77SJunchao Zhang #include <petsc/private/fortranimpl.h>
1485649d77SJunchao Zhang #endif
1585649d77SJunchao Zhang 
1656883f60SBarry Smith #if defined(PETSC_USE_GCOV)
1756883f60SBarry Smith EXTERN_C_BEGIN
1856883f60SBarry Smith void  __gcov_flush(void);
1956883f60SBarry Smith EXTERN_C_END
2056883f60SBarry Smith #endif
218101f56cSMatthew Knepley 
222d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
2395c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
242d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
252d53ad75SBarry Smith #endif
262d53ad75SBarry Smith 
2727104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
2816ad0300SBarry Smith #include <petscviewersaws.h>
29a6790183SBarry Smith #endif
30eb02082bSJunchao Zhang 
31e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
32e5c89e4eSSatish Balay 
3395c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
34e5c89e4eSSatish Balay 
3595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
3795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
410069ddf5SShri Abhyankar 
426de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
43e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
4427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
456de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
466de5d289SStefano Zampini #else
476de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
486de5d289SStefano Zampini #endif
49e5c89e4eSSatish Balay 
50480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
51480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
52480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
535f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
54480cf27aSJed Brown 
55e5c89e4eSSatish Balay /*
56e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
57e5c89e4eSSatish Balay */
5802c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
5902c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
60e5c89e4eSSatish Balay 
61ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
62ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
630f8e0872SSatish Balay 
64a2f94806SJed Brown PetscInt PetscHotRegionDepth;
65a2f94806SJed Brown 
666edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
676edef35eSSatish Balay 
68b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
69b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
70b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
71b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
72b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
73b22622e2STadeu Manoel #endif
74b22622e2STadeu Manoel 
75e5c89e4eSSatish Balay /*
76945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
7772a42c3cSBarry Smith 
7872a42c3cSBarry Smith    Collective
7972a42c3cSBarry Smith 
8072a42c3cSBarry Smith    Level: advanced
8172a42c3cSBarry Smith 
8295452b02SPatrick Sanan     Notes:
83a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
840f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
85a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
860f11a792SBarry Smith 
87a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
881ea5a559SBarry Smith 
8972a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
900f11a792SBarry Smith */
91945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
9272a42c3cSBarry Smith {
9372a42c3cSBarry Smith   PetscErrorCode ierr;
9472a42c3cSBarry Smith   int            myargc   = argc;
9572a42c3cSBarry Smith   char           **myargs = args;
9672a42c3cSBarry Smith 
9772a42c3cSBarry Smith   PetscFunctionBegin;
9827104ee2SJacob Faibussowitsch   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) PetscFunctionReturn(ierr);
991ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
100df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
10127104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
10272a42c3cSBarry Smith }
10372a42c3cSBarry Smith 
104f0865b08SBarry Smith /*
105a56f64adSBarry Smith       Used by Julia interface to get communicator
106f0865b08SBarry Smith */
107945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
108f0865b08SBarry Smith {
109f0865b08SBarry Smith   PetscFunctionBegin;
11027104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidPointer(comm,1);
111f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
112f0865b08SBarry Smith   PetscFunctionReturn(0);
113f0865b08SBarry Smith }
114f0865b08SBarry Smith 
115e5c89e4eSSatish Balay /*@C
116e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
117e5c89e4eSSatish Balay         the command line arguments.
118e5c89e4eSSatish Balay 
119e5c89e4eSSatish Balay    Collective
120e5c89e4eSSatish Balay 
121e5c89e4eSSatish Balay    Level: advanced
122e5c89e4eSSatish Balay 
123e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
124e5c89e4eSSatish Balay @*/
1257087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
126e5c89e4eSSatish Balay {
127e5c89e4eSSatish Balay   PetscErrorCode ierr;
128e5c89e4eSSatish Balay   int            argc   = 0;
12902c9f0b5SLisandro Dalcin   char           **args = NULL;
130e5c89e4eSSatish Balay 
131e5c89e4eSSatish Balay   PetscFunctionBegin;
1320298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
133e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
134e5c89e4eSSatish Balay }
135e5c89e4eSSatish Balay 
136e5c89e4eSSatish Balay /*@
137e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
138e5c89e4eSSatish Balay 
13993b6d2d1SJed Brown    Level: beginner
140e5c89e4eSSatish Balay 
141e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
142e5c89e4eSSatish Balay @*/
1437087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized)
144e5c89e4eSSatish Balay {
14527104ee2SJacob Faibussowitsch   PetscFunctionBegin;
14627104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized,1);
147e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
14827104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
149e5c89e4eSSatish Balay }
150e5c89e4eSSatish Balay 
151e5c89e4eSSatish Balay /*@
152e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
153e5c89e4eSSatish Balay 
154e5c89e4eSSatish Balay    Level: developer
155e5c89e4eSSatish Balay 
156e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
157e5c89e4eSSatish Balay @*/
1587087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
159e5c89e4eSSatish Balay {
16027104ee2SJacob Faibussowitsch   PetscFunctionBegin;
16127104ee2SJacob Faibussowitsch   if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized,1);
162e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
16327104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
164e5c89e4eSSatish Balay }
165e5c89e4eSSatish Balay 
16657171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []);
167e5c89e4eSSatish Balay 
168e5c89e4eSSatish Balay /*
169e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
170e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
171e5c89e4eSSatish Balay */
172367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
173e5c89e4eSSatish Balay 
174367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
175e5c89e4eSSatish Balay {
176e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
177e5c89e4eSSatish Balay 
178e5c89e4eSSatish Balay   PetscFunctionBegin;
179e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
180e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
18141e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
182e5c89e4eSSatish Balay   }
183e5c89e4eSSatish Balay 
184e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
185e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
186e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
187e5c89e4eSSatish Balay   }
188812af9f3SBarry Smith   PetscFunctionReturnVoid();
189e5c89e4eSSatish Balay }
190e5c89e4eSSatish Balay 
191e5c89e4eSSatish Balay /*
192e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
193e5c89e4eSSatish Balay sum of the second entry.
194b693b147SBarry Smith 
19576f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
196367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
197b693b147SBarry Smith there would be no place to store the both needed results.
198e5c89e4eSSatish Balay */
19976ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
200e5c89e4eSSatish Balay {
201e5c89e4eSSatish Balay   PetscErrorCode ierr;
202e5c89e4eSSatish Balay 
203e5c89e4eSSatish Balay   PetscFunctionBegin;
204d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
205d6e4c47cSJed Brown   {
206d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
207ffc4695bSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRMPI(ierr);
208d6e4c47cSJed Brown     *max = work.max;
209d6e4c47cSJed Brown     *sum = work.sum;
210d6e4c47cSJed Brown   }
211d6e4c47cSJed Brown #else
212d6e4c47cSJed Brown   {
213d6e4c47cSJed Brown     PetscMPIInt    size,rank;
214d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
215ffc4695bSBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
216ffc4695bSBarry Smith     ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
217785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
218820f2d46SBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRMPI(ierr);
2196ac3741eSJed Brown     *max = work[rank].max;
2206ac3741eSJed Brown     *sum = work[rank].sum;
221e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
222d6e4c47cSJed Brown   }
223d6e4c47cSJed Brown #endif
224e5c89e4eSSatish Balay   PetscFunctionReturn(0);
225e5c89e4eSSatish Balay }
226e5c89e4eSSatish Balay 
227e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
228e5c89e4eSSatish Balay 
229de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
23006a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
231e5c89e4eSSatish Balay 
232092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
233e5c89e4eSSatish Balay {
234e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
235e5c89e4eSSatish Balay 
236e5c89e4eSSatish Balay   PetscFunctionBegin;
2377c2de775SJed Brown   if (*datatype == MPIU_REAL) {
238e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
239a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2407c2de775SJed Brown   }
2417c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
242c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2437c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
244a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2457c2de775SJed Brown   }
2467c2de775SJed Brown #endif
2477c2de775SJed Brown   else {
2487c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
24941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
250e2e03761SBarry Smith   }
251812af9f3SBarry Smith   PetscFunctionReturnVoid();
252e5c89e4eSSatish Balay }
253e5c89e4eSSatish Balay #endif
254e5c89e4eSSatish Balay 
255570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
256d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
257d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
258d9822059SBarry Smith 
259092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
260d9822059SBarry Smith {
261d9822059SBarry Smith   PetscInt i,count = *cnt;
262d9822059SBarry Smith 
263d9822059SBarry Smith   PetscFunctionBegin;
2647c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2658c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
266a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2677c2de775SJed Brown   }
2687c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2697c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2707c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2717c2de775SJed Brown     for (i=0; i<count; i++) {
2727c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2737c2de775SJed Brown     }
2747c2de775SJed Brown   }
2757c2de775SJed Brown #endif
2767c2de775SJed Brown   else {
2777c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
27841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2798c764dc5SJose Roman   }
280d9822059SBarry Smith   PetscFunctionReturnVoid();
281d9822059SBarry Smith }
282d9822059SBarry Smith 
283092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
284d9822059SBarry Smith {
285d9822059SBarry Smith   PetscInt    i,count = *cnt;
286d9822059SBarry Smith 
287d9822059SBarry Smith   PetscFunctionBegin;
2887c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2898c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
290a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2917c2de775SJed Brown   }
2927c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2937c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2947c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2957c2de775SJed Brown     for (i=0; i<count; i++) {
2967c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2977c2de775SJed Brown     }
2987c2de775SJed Brown   }
2997c2de775SJed Brown #endif
3007c2de775SJed Brown   else {
3018c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
30241e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
3038c764dc5SJose Roman   }
304d9822059SBarry Smith   PetscFunctionReturnVoid();
305d9822059SBarry Smith }
306d9822059SBarry Smith #endif
307d9822059SBarry Smith 
308480cf27aSJed Brown /*
309480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
310480cf27aSJed Brown 
311ff0e51ddSBarry 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.
312480cf27aSJed Brown 
31312801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
314480cf27aSJed Brown 
315480cf27aSJed Brown */
31633779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
317480cf27aSJed Brown {
318480cf27aSJed Brown   PetscErrorCode        ierr;
31905342407SJunchao Zhang   PetscCommCounter      *counter=(PetscCommCounter*)count_val;
32057f21012SBarry Smith   struct PetscCommStash *comms = counter->comms, *pcomm;
321480cf27aSJed Brown 
322480cf27aSJed Brown   PetscFunctionBegin;
3237d3de750SJacob Faibussowitsch   ierr = PetscInfo(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
32405342407SJunchao Zhang   ierr = PetscFree(counter->iflags);CHKERRMPI(ierr);
32557f21012SBarry Smith   while (comms) {
32657f21012SBarry Smith     ierr  = MPI_Comm_free(&comms->comm);CHKERRMPI(ierr);
32757f21012SBarry Smith     pcomm = comms;
32857f21012SBarry Smith     comms = comms->next;
32957f21012SBarry Smith     ierr  = PetscFree(pcomm);CHKERRQ(ierr);
33057f21012SBarry Smith   }
33105342407SJunchao Zhang   ierr = PetscFree(counter);CHKERRMPI(ierr);
332480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
333480cf27aSJed Brown }
334480cf27aSJed Brown 
335480cf27aSJed Brown /*
33647435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
337da3039f7SJed Brown   calls MPI_Comm_free().
338da3039f7SJed Brown 
339da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
340480cf27aSJed Brown 
341ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
342480cf27aSJed Brown 
34312801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
344480cf27aSJed Brown 
345480cf27aSJed Brown */
34633779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
347480cf27aSJed Brown {
348480cf27aSJed Brown   PetscErrorCode                    ierr;
34933779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
350480cf27aSJed Brown 
351480cf27aSJed Brown   PetscFunctionBegin;
35212801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
353ec4fadc2SJed Brown   icomm.ptr = attr_val;
35476bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
35533779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
35633779a13SJunchao Zhang     PetscMPIInt flg;
35733779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
35847435625SJed Brown     ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
35933779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
36033779a13SJunchao 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");
36133779a13SJunchao Zhang   }
36247435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
3637d3de750SJacob Faibussowitsch   ierr = PetscInfo(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);CHKERRMPI(ierr);
364da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
365b89831e5SBarry Smith }
366da3039f7SJed Brown 
367da3039f7SJed Brown /*
36833779a13SJunchao 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.
369da3039f7SJed Brown  */
37033779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
371da3039f7SJed Brown {
372da3039f7SJed Brown   PetscErrorCode ierr;
373da3039f7SJed Brown 
374da3039f7SJed Brown   PetscFunctionBegin;
3757d3de750SJacob Faibussowitsch   ierr = PetscInfo(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
376480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
377480cf27aSJed Brown }
378480cf27aSJed Brown 
37933779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
38042218b76SBarry Smith 
381951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3828cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3838cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3848cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
385e39fd77fSBarry Smith #endif
386e39fd77fSBarry Smith 
3870f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
38812801b39SBarry Smith 
389eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
390eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3916ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
39202c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
393dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
394e5c89e4eSSatish Balay 
395dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
396051e4cf2SJed Brown {
397051e4cf2SJed Brown   PetscErrorCode ierr;
398051e4cf2SJed Brown 
399051e4cf2SJed Brown   PetscFunctionBegin;
400051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
401a1601671SBarry Smith   ierr = PetscCitationsRegister("@TechReport{petsc-user-ref,\n  Author = {Satish Balay and Shrirang Abhyankar and Mark F. Adams and Jed Brown \n            and Peter Brune and Kris Buschelman and Lisandro Dalcin and\n            Victor Eijkhout and William D. Gropp and Dmitry Karpeyev and\n            Dinesh Kaushik and Matthew G. Knepley and Dave A. May and Lois Curfman McInnes\n            and Richard Tran Mills and Todd Munson and Karl Rupp and Patrick Sanan\n            and Barry F. Smith and Stefano Zampini and Hong Zhang and Hong Zhang},\n  Title = {{PETS}c Users Manual},\n  Number = {ANL-95/11 - Revision 3.11},\n  Institution = {Argonne National Laboratory},\n  Year = {2019}\n}\n",NULL);CHKERRQ(ierr);
402051e4cf2SJed Brown   ierr = PetscCitationsRegister("@InProceedings{petsc-efficient,\n  Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n  Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n  Booktitle = {Modern Software Tools in Scientific Computing},\n  Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n  Pages = {163--202},\n  Publisher = {Birkh{\\\"{a}}user Press},\n  Year = {1997}\n}\n",NULL);CHKERRQ(ierr);
403051e4cf2SJed Brown   PetscFunctionReturn(0);
404051e4cf2SJed Brown }
405e5c89e4eSSatish Balay 
4062d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4072d747510SLisandro Dalcin 
4082d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4092d747510SLisandro Dalcin {
4102d747510SLisandro Dalcin   PetscErrorCode ierr;
4112d747510SLisandro Dalcin 
4122d747510SLisandro Dalcin   PetscFunctionBegin;
413589a23caSBarry Smith   ierr  = PetscStrncpy(programname,name,sizeof(programname));CHKERRQ(ierr);
4142d747510SLisandro Dalcin   PetscFunctionReturn(0);
4152d747510SLisandro Dalcin }
4162d747510SLisandro Dalcin 
4172d747510SLisandro Dalcin /*@C
4182d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4192d747510SLisandro Dalcin 
4202d747510SLisandro Dalcin     Not Collective
4212d747510SLisandro Dalcin 
4222d747510SLisandro Dalcin     Input Parameter:
4232d747510SLisandro Dalcin .   len - length of the string name
4242d747510SLisandro Dalcin 
4252d747510SLisandro Dalcin     Output Parameter:
4262d747510SLisandro Dalcin .   name - the name of the running program
4272d747510SLisandro Dalcin 
4282d747510SLisandro Dalcin    Level: advanced
4292d747510SLisandro Dalcin 
4302d747510SLisandro Dalcin     Notes:
4312d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4322d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4332d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4342d747510SLisandro Dalcin @*/
4352d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4362d747510SLisandro Dalcin {
4372d747510SLisandro Dalcin   PetscErrorCode ierr;
4382d747510SLisandro Dalcin 
4392d747510SLisandro Dalcin   PetscFunctionBegin;
4402d747510SLisandro Dalcin   ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4412d747510SLisandro Dalcin   PetscFunctionReturn(0);
4422d747510SLisandro Dalcin }
4432d747510SLisandro Dalcin 
444e5c89e4eSSatish Balay /*@C
445e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
446e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
447e5c89e4eSSatish Balay 
448e5c89e4eSSatish Balay    Not Collective
449e5c89e4eSSatish Balay 
450e5c89e4eSSatish Balay    Output Parameters:
451e5c89e4eSSatish Balay +  argc - count of number of command line arguments
452e5c89e4eSSatish Balay -  args - the command line arguments
453e5c89e4eSSatish Balay 
454e5c89e4eSSatish Balay    Level: intermediate
455e5c89e4eSSatish Balay 
456e5c89e4eSSatish Balay    Notes:
457e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
458e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
459e5c89e4eSSatish Balay 
460f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
461f177e3b1SBarry Smith 
462793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
463e5c89e4eSSatish Balay 
464e5c89e4eSSatish Balay @*/
4657087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
466e5c89e4eSSatish Balay {
467e5c89e4eSSatish Balay   PetscFunctionBegin;
468*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!PetscInitializeCalled && PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
469e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
470e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
471e5c89e4eSSatish Balay   PetscFunctionReturn(0);
472e5c89e4eSSatish Balay }
473e5c89e4eSSatish Balay 
474793721a6SBarry Smith /*@C
475793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
476793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
477793721a6SBarry Smith 
478793721a6SBarry Smith    Not Collective
479793721a6SBarry Smith 
480793721a6SBarry Smith    Output Parameters:
481793721a6SBarry Smith .  args - the command line arguments
482793721a6SBarry Smith 
483793721a6SBarry Smith    Level: intermediate
484793721a6SBarry Smith 
485793721a6SBarry Smith    Notes:
486793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
487793721a6SBarry Smith 
488793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
489793721a6SBarry Smith 
490793721a6SBarry Smith @*/
4917087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
492793721a6SBarry Smith {
493793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
494793721a6SBarry Smith   PetscErrorCode ierr;
495793721a6SBarry Smith 
496793721a6SBarry Smith   PetscFunctionBegin;
497*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!PetscInitializeCalled && PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4982d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
499785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
500793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
501793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
502793721a6SBarry Smith   }
5032d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
504793721a6SBarry Smith   PetscFunctionReturn(0);
505793721a6SBarry Smith }
506793721a6SBarry Smith 
507793721a6SBarry Smith /*@C
508793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
509793721a6SBarry Smith 
510793721a6SBarry Smith    Not Collective
511793721a6SBarry Smith 
512793721a6SBarry Smith    Output Parameters:
513793721a6SBarry Smith .  args - the command line arguments
514793721a6SBarry Smith 
515793721a6SBarry Smith    Level: intermediate
516793721a6SBarry Smith 
517793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
518793721a6SBarry Smith 
519793721a6SBarry Smith @*/
5207087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
521793721a6SBarry Smith {
522793721a6SBarry Smith   PetscInt       i = 0;
523793721a6SBarry Smith   PetscErrorCode ierr;
524793721a6SBarry Smith 
525793721a6SBarry Smith   PetscFunctionBegin;
526a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
527793721a6SBarry Smith   while (args[i]) {
528793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
529793721a6SBarry Smith     i++;
530793721a6SBarry Smith   }
531793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
532793721a6SBarry Smith   PetscFunctionReturn(0);
533793721a6SBarry Smith }
534793721a6SBarry Smith 
53527104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
53630befbd2SBarry Smith #include <petscconfiginfo.h>
53730befbd2SBarry Smith 
53895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
53911525c0dSBarry Smith {
54027104ee2SJacob Faibussowitsch   PetscFunctionBegin;
54111525c0dSBarry Smith   if (!PetscGlobalRank) {
54230befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
54311525c0dSBarry Smith     int            port;
544ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
54511525c0dSBarry Smith     size_t         applinelen,introlen;
54611525c0dSBarry Smith     PetscErrorCode ierr;
547ffbd1cfbSBarry Smith     char           sawsurl[256];
54811525c0dSBarry Smith 
549c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
55011525c0dSBarry Smith     if (flg) {
55111525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
55211525c0dSBarry Smith 
553589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL);CHKERRQ(ierr);
55411525c0dSBarry Smith       if (sawslog[0]) {
55511525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
55611525c0dSBarry Smith       } else {
55711525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
55811525c0dSBarry Smith       }
55911525c0dSBarry Smith     }
560589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg);CHKERRQ(ierr);
56111525c0dSBarry Smith     if (flg) {
56211525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
56311525c0dSBarry Smith     }
564c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
565ffbd1cfbSBarry Smith     if (selectport) {
566ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
567ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
568ffbd1cfbSBarry Smith     } else {
569c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
57011525c0dSBarry Smith       if (flg) {
57111525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
57211525c0dSBarry Smith       }
573ffbd1cfbSBarry Smith     }
574589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg);CHKERRQ(ierr);
57511525c0dSBarry Smith     if (flg) {
5762f613bf5SBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
57711525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5789c1e0ce8SBarry Smith     } else {
579c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5809c1e0ce8SBarry Smith       if (flg) {
581589a23caSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root));CHKERRQ(ierr);
5822f613bf5SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
5839c1e0ce8SBarry Smith       }
58411525c0dSBarry Smith     }
585c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
58611525c0dSBarry Smith     if (flg2) {
58711525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
588*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
589589a23caSBarry Smith       ierr = PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root);CHKERRQ(ierr);
59011525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
591*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
5922f613bf5SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());
59311525c0dSBarry Smith     }
594589a23caSBarry Smith     ierr = PetscGetProgramName(programname,sizeof(programname));CHKERRQ(ierr);
59511525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
59611525c0dSBarry Smith     introlen   = 4096 + applinelen;
59730a8c9c0SSurtai Han     applinelen += 1024;
59811525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
59911525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
60011525c0dSBarry Smith 
60111525c0dSBarry Smith     if (rootlocal) {
60211525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
60311525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
60411525c0dSBarry Smith     }
60576a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
60611525c0dSBarry Smith     if (rootlocal && help) {
607928bb9adSStefano Zampini       ierr = 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);CHKERRQ(ierr);
60811525c0dSBarry Smith     } else if (help) {
609928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
61011525c0dSBarry Smith     } else {
611928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
61211525c0dSBarry Smith     }
613b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
61430befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
61511525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
616a17b96a8SKyle 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"
617df62c222SBarry 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"
618928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
61943da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
62011525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
62111525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
622ffbd1cfbSBarry Smith     if (selectport) {
623aa573868SBarry Smith       PetscBool silent;
6247d812c46SBarry Smith 
6257d812c46SBarry Smith       ierr = SAWs_Initialize();
6267d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6277d812c46SBarry Smith       while (ierr) {
6287d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6297d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6307d812c46SBarry Smith         ierr = SAWs_Initialize();
6317d812c46SBarry Smith       }
6327d812c46SBarry Smith 
633aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
634aa573868SBarry Smith       if (!silent) {
635ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
636ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
637ffbd1cfbSBarry Smith       }
6387d812c46SBarry Smith     } else {
6397d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
640aa573868SBarry Smith     }
6410af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6420af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6430af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6440af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6450af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
64611525c0dSBarry Smith   }
64711525c0dSBarry Smith   PetscFunctionReturn(0);
64811525c0dSBarry Smith }
64911525c0dSBarry Smith #endif
65011525c0dSBarry Smith 
6514dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6524dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6534dfee713SSatish Balay {
6544dfee713SSatish Balay   PetscFunctionBegin;
6554dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6564dfee713SSatish Balay     /* see MPI.py for details on this bug */
6574dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6584dfee713SSatish Balay #endif
6594dfee713SSatish Balay   PetscFunctionReturn(0);
6604dfee713SSatish Balay }
6614dfee713SSatish Balay 
662a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS)
663a56f64adSBarry Smith #include <adios.h>
66422580e64SBarry Smith #include <adios_read.h>
6657b56e58cSSatish Balay int64_t Petsc_adios_group;
666a56f64adSBarry Smith #endif
667a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP)
668cd1b99a6SBarry Smith #include <omp.h>
669f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
670cd1b99a6SBarry Smith #endif
671a56f64adSBarry Smith 
672a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_DEVICE)
673a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
674a4af0ceeSJacob Faibussowitsch #  if PetscDefined(HAVE_CUDA)
675a4af0ceeSJacob Faibussowitsch // REMOVE ME
676a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL;
677a4af0ceeSJacob Faibussowitsch #  endif
678a4af0ceeSJacob Faibussowitsch #  if PetscDefined(HAVE_HIP)
679a4af0ceeSJacob Faibussowitsch // REMOVE ME
680a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL;
681a4af0ceeSJacob Faibussowitsch #  endif
682a4af0ceeSJacob Faibussowitsch #endif
683a4af0ceeSJacob Faibussowitsch 
68427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
685bc8a24ecSBarry Smith #include <dlfcn.h>
686bc8a24ecSBarry Smith #endif
687a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
688a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
689a4af0ceeSJacob Faibussowitsch #endif
690a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
691a4af0ceeSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViennaCLInit();
692a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE;
693a4af0ceeSJacob Faibussowitsch #endif
694bc8a24ecSBarry Smith 
69585649d77SJunchao Zhang /*
69685649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
69785649d77SJunchao Zhang 
69885649d77SJunchao Zhang   prog:     program name
69902101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
70085649d77SJunchao Zhang   help:     program help message
70102101c96SSatish Balay   ftn:      is it called from Fortran initilization (petscinitializef_)?
70285649d77SJunchao Zhang   readarguments,len: used when fortran is true
70385649d77SJunchao Zhang */
70402101c96SSatish Balay PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char* prog,const char* file,const char *help,PetscBool ftn,PetscBool readarguments,PetscInt len)
70585649d77SJunchao Zhang {
70685649d77SJunchao Zhang   PetscErrorCode ierr;
70785649d77SJunchao Zhang   PetscMPIInt    size;
70885649d77SJunchao Zhang   PetscBool      flg = PETSC_TRUE;
70985649d77SJunchao Zhang   char           hostname[256];
71085649d77SJunchao Zhang 
71127104ee2SJacob Faibussowitsch   PetscFunctionBegin;
71227104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
71385649d77SJunchao Zhang   /*
71485649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
71585649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
71685649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
71785649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
71885649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
71985649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
72085649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
72185649d77SJunchao Zhang       listed above and since that time are compatible.
72285649d77SJunchao Zhang 
72385649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
72485649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
72585649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
72685649d77SJunchao Zhang       to perform the checking.
72785649d77SJunchao Zhang 
72885649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
72985649d77SJunchao Zhang 
73085649d77SJunchao Zhang       Questions:
73185649d77SJunchao Zhang 
73285649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
73385649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
73485649d77SJunchao Zhang   */
73585649d77SJunchao Zhang 
73685649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
73785649d77SJunchao Zhang   {
73885649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
73985649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
74027104ee2SJacob Faibussowitsch     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);
74127104ee2SJacob Faibussowitsch     if (ierr) PetscFunctionReturn(ierr);
74285649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
74385649d77SJunchao Zhang #if defined(MPICH_VERSION)
74485649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000
74585649d77SJunchao Zhang     {
74685649d77SJunchao Zhang       char      *ver,*lf;
74785649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
74827104ee2SJacob Faibussowitsch       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);
74927104ee2SJacob Faibussowitsch       if (ierr) PetscFunctionReturn(ierr);
75027104ee2SJacob Faibussowitsch       else if (ver) {
75127104ee2SJacob Faibussowitsch         ierr = PetscStrchr(ver,'\n',&lf);
75227104ee2SJacob Faibussowitsch         if (ierr) PetscFunctionReturn(ierr);
75327104ee2SJacob Faibussowitsch         else if (lf) {
75485649d77SJunchao Zhang           *lf = 0;
75527104ee2SJacob Faibussowitsch           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) PetscFunctionReturn(ierr);
75685649d77SJunchao Zhang         }
75785649d77SJunchao Zhang       }
75885649d77SJunchao Zhang       if (!flg) {
7597d3de750SJacob Faibussowitsch         PetscInfo(NULL,"PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n",mpilibraryversion,MPICH_VESION);
76085649d77SJunchao Zhang         flg = PETSC_TRUE;
76185649d77SJunchao Zhang       }
76285649d77SJunchao Zhang     }
76385649d77SJunchao Zhang #endif
76485649d77SJunchao 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?) */
76585649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION)
76685649d77SJunchao Zhang     {
76785649d77SJunchao Zhang       char *ver,bs[MPI_MAX_LIBRARY_VERSION_STRING],*bsf;
76885649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
76985649d77SJunchao Zhang #define PSTRSZ 2
77085649d77SJunchao Zhang       char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI","FUJITSU MPI"};
77185649d77SJunchao Zhang       char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v","Library "};
77285649d77SJunchao Zhang       int i;
77385649d77SJunchao Zhang       for (i=0; i<PSTRSZ; i++) {
77427104ee2SJacob Faibussowitsch         ierr = PetscStrstr(mpilibraryversion,ompistr1[i],&ver);
77527104ee2SJacob Faibussowitsch         if (ierr) PetscFunctionReturn(ierr);
77627104ee2SJacob Faibussowitsch         else if (ver) {
77785649d77SJunchao Zhang           PetscSNPrintf(bs,MPI_MAX_LIBRARY_VERSION_STRING,"%s%d.%d",ompistr2[i],OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
77827104ee2SJacob Faibussowitsch           ierr = PetscStrstr(ver,bs,&bsf);
77927104ee2SJacob Faibussowitsch           if (ierr) PetscFunctionReturn(ierr);
78027104ee2SJacob Faibussowitsch           else if (bsf) flg = PETSC_TRUE;
78185649d77SJunchao Zhang           break;
78285649d77SJunchao Zhang         }
78385649d77SJunchao Zhang       }
78485649d77SJunchao Zhang       if (!flg) {
7857d3de750SJacob Faibussowitsch         PetscInfo(NULL,"PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n",mpilibraryversion,OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
78685649d77SJunchao Zhang         flg = PETSC_TRUE;
78785649d77SJunchao Zhang       }
78885649d77SJunchao Zhang     }
78985649d77SJunchao Zhang #endif
79085649d77SJunchao Zhang   }
79185649d77SJunchao Zhang #endif
79285649d77SJunchao Zhang 
79385649d77SJunchao Zhang #if defined(PETSC_HAVE_DLSYM)
79485649d77SJunchao 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 */
79527104ee2SJacob Faibussowitsch   if (PetscUnlikely(dlsym(RTLD_DEFAULT,"ompi_mpi_init") && dlsym(RTLD_DEFAULT,"MPID_Abort"))) {
79685649d77SJunchao Zhang     fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
79727104ee2SJacob Faibussowitsch     ierr = PetscStackView(stderr);CHKERRQ(ierr);
79827104ee2SJacob Faibussowitsch     PetscFunctionReturn(PETSC_ERR_MPI_LIB_INCOMP);
79985649d77SJunchao Zhang   }
80085649d77SJunchao Zhang #endif
80185649d77SJunchao Zhang 
80285649d77SJunchao Zhang   /* these must be initialized in a routine, not as a constant declaration*/
80385649d77SJunchao Zhang   PETSC_STDOUT = stdout;
80485649d77SJunchao Zhang   PETSC_STDERR = stderr;
80585649d77SJunchao Zhang 
80685649d77SJunchao Zhang   /*CHKERRQ can be used from now */
80785649d77SJunchao Zhang   PetscErrorHandlingInitialized = PETSC_TRUE;
80885649d77SJunchao Zhang 
80985649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
81085649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
81185649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
81285649d77SJunchao Zhang #endif
81385649d77SJunchao Zhang 
81485649d77SJunchao Zhang   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
81585649d77SJunchao Zhang 
81685649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
81785649d77SJunchao Zhang 
81885649d77SJunchao Zhang   ierr = PetscSetProgramName(prog);CHKERRQ(ierr);
81985649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
82085649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
82185649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
82285649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
82385649d77SJunchao Zhang 
82485649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
82585649d77SJunchao Zhang   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRMPI(ierr);
82685649d77SJunchao Zhang 
82785649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
82885649d77SJunchao Zhang     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRMPI(ierr);
82985649d77SJunchao Zhang     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRMPI(ierr);
83085649d77SJunchao Zhang   }
83185649d77SJunchao Zhang 
83285649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
83385649d77SJunchao Zhang   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
83485649d77SJunchao Zhang 
83585649d77SJunchao Zhang   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRMPI(ierr);
83685649d77SJunchao Zhang   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRMPI(ierr);
83785649d77SJunchao Zhang 
83885649d77SJunchao Zhang   MPIU_BOOL = MPI_INT;
83985649d77SJunchao Zhang   MPIU_ENUM = MPI_INT;
84085649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
84185649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
84285649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
84385649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
84485649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
84585649d77SJunchao Zhang #endif
84627104ee2SJacob Faibussowitsch   else {
84727104ee2SJacob Faibussowitsch     (*PetscErrorPrintf)("PetscInitialize_Common: Could not find MPI type for size_t\n");
84827104ee2SJacob Faibussowitsch     PetscFunctionReturn(PETSC_ERR_SUP_SYS);
84927104ee2SJacob Faibussowitsch   }
85085649d77SJunchao Zhang 
85185649d77SJunchao Zhang   /*
85285649d77SJunchao Zhang      Initialized the global complex variable; this is because with
85385649d77SJunchao Zhang      shared libraries the constructors for global variables
85485649d77SJunchao Zhang      are not called; at least on IRIX.
85585649d77SJunchao Zhang   */
85685649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
85785649d77SJunchao Zhang   {
85885649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
85985649d77SJunchao Zhang     PetscComplex ic(0.0,1.0);
86085649d77SJunchao Zhang     PETSC_i = ic;
86185649d77SJunchao Zhang #else
86285649d77SJunchao Zhang     PETSC_i = _Complex_I;
86385649d77SJunchao Zhang #endif
86485649d77SJunchao Zhang   }
86585649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
86685649d77SJunchao Zhang 
86785649d77SJunchao Zhang   /*
86885649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
86985649d77SJunchao Zhang      half of the entries and maxes the second half.
87085649d77SJunchao Zhang   */
87185649d77SJunchao Zhang   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRMPI(ierr);
87285649d77SJunchao Zhang 
87385649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128)
87485649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRMPI(ierr);
87585649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRMPI(ierr);
87685649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
87785649d77SJunchao Zhang   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRMPI(ierr);
87885649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRMPI(ierr);
87985649d77SJunchao Zhang #endif
88085649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRMPI(ierr);
88185649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRMPI(ierr);
88285649d77SJunchao Zhang #elif defined(PETSC_USE_REAL___FP16)
88385649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRMPI(ierr);
88485649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRMPI(ierr);
88585649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRMPI(ierr);
88685649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRMPI(ierr);
88785649d77SJunchao Zhang #endif
88885649d77SJunchao Zhang 
88985649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
89085649d77SJunchao Zhang   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRMPI(ierr);
89185649d77SJunchao Zhang #endif
89285649d77SJunchao Zhang 
89385649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRMPI(ierr);
89485649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRMPI(ierr);
89585649d77SJunchao Zhang 
89685649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
89785649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
89885649d77SJunchao Zhang   {
89985649d77SJunchao Zhang     struct PetscRealInt { PetscReal v; PetscInt i; };
90085649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
90185649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscRealInt,v),offsetof(struct PetscRealInt,i)};
90285649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_REAL,MPIU_INT}, tmpStruct;
90385649d77SJunchao Zhang 
90485649d77SJunchao Zhang     ierr = MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct);CHKERRMPI(ierr);
90585649d77SJunchao Zhang     ierr = MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscRealInt),&MPIU_REAL_INT);CHKERRMPI(ierr);
90685649d77SJunchao Zhang     ierr = MPI_Type_free(&tmpStruct);CHKERRMPI(ierr);
90785649d77SJunchao Zhang     ierr = MPI_Type_commit(&MPIU_REAL_INT);CHKERRMPI(ierr);
90885649d77SJunchao Zhang   }
90985649d77SJunchao Zhang   {
91085649d77SJunchao Zhang     struct PetscScalarInt { PetscScalar v; PetscInt i; };
91185649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
91285649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscScalarInt,v),offsetof(struct PetscScalarInt,i)};
91385649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_SCALAR,MPIU_INT}, tmpStruct;
91485649d77SJunchao Zhang 
91585649d77SJunchao Zhang     ierr = MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct);CHKERRMPI(ierr);
91685649d77SJunchao Zhang     ierr = MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscScalarInt),&MPIU_SCALAR_INT);CHKERRMPI(ierr);
91785649d77SJunchao Zhang     ierr = MPI_Type_free(&tmpStruct);CHKERRMPI(ierr);
91885649d77SJunchao Zhang     ierr = MPI_Type_commit(&MPIU_SCALAR_INT);CHKERRMPI(ierr);
91985649d77SJunchao Zhang   }
92085649d77SJunchao Zhang #endif
92185649d77SJunchao Zhang 
92285649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
92385649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRMPI(ierr);
92485649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRMPI(ierr);
92585649d77SJunchao Zhang #endif
92685649d77SJunchao Zhang   ierr = MPI_Type_contiguous(4,MPI_INT,&MPI_4INT);CHKERRMPI(ierr);
92785649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPI_4INT);CHKERRMPI(ierr);
92885649d77SJunchao Zhang   ierr = MPI_Type_contiguous(4,MPIU_INT,&MPIU_4INT);CHKERRMPI(ierr);
92985649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU_4INT);CHKERRMPI(ierr);
93085649d77SJunchao Zhang 
93185649d77SJunchao Zhang   /*
93285649d77SJunchao Zhang      Attributes to be set on PETSc communicators
93385649d77SJunchao Zhang   */
93485649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);CHKERRMPI(ierr);
93585649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);CHKERRMPI(ierr);
93685649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);CHKERRMPI(ierr);
93785649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);CHKERRMPI(ierr);
93885649d77SJunchao Zhang 
93985649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
94002101c96SSatish Balay   if (ftn) {ierr = PetscInitFortran_Private(readarguments,file,len);CHKERRQ(ierr);}
94185649d77SJunchao Zhang   else
94285649d77SJunchao Zhang #endif
94385649d77SJunchao Zhang   {ierr = PetscOptionsInsert(NULL,&PetscGlobalArgc,&PetscGlobalArgs,file);CHKERRQ(ierr);}
94485649d77SJunchao Zhang 
94585649d77SJunchao Zhang   /* call a second time so it can look in the options database */
94685649d77SJunchao Zhang   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
94785649d77SJunchao Zhang 
94885649d77SJunchao Zhang   /*
94985649d77SJunchao Zhang      Check system options and print help
95085649d77SJunchao Zhang   */
95185649d77SJunchao Zhang   ierr = PetscOptionsCheckInitial_Private(help);CHKERRQ(ierr);
95285649d77SJunchao Zhang 
953a4af0ceeSJacob Faibussowitsch   /*
954a4af0ceeSJacob Faibussowitsch    Initialize PetscDevice and PetscDeviceContext
955a4af0ceeSJacob Faibussowitsch 
956a4af0ceeSJacob Faibussowitsch    Note to any future devs thinking of moving this, proper initialization requires:
957a4af0ceeSJacob Faibussowitsch    1. MPI initialized
958a4af0ceeSJacob Faibussowitsch    2. Options DB initialized
959a4af0ceeSJacob Faibussowitsch    3. Petsc error handling initialized, specifically signal handlers. This expects to set up its own SIGSEV handler via
960a4af0ceeSJacob Faibussowitsch       the push/pop interface.
961a4af0ceeSJacob Faibussowitsch   */
962a2158755SJunchao Zhang #if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP) || PetscDefined(HAVE_SYCL))
963a4af0ceeSJacob Faibussowitsch   ierr = PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD);CHKERRQ(ierr);
964a4af0ceeSJacob Faibussowitsch #endif
965a4af0ceeSJacob Faibussowitsch 
966a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL)
967a4af0ceeSJacob Faibussowitsch   flg = PETSC_FALSE;
968a4af0ceeSJacob Faibussowitsch   ierr = PetscOptionsHasName(NULL,NULL,"-log_summary",&flg);CHKERRQ(ierr);
969a4af0ceeSJacob Faibussowitsch   if (!flg) {ierr = PetscOptionsHasName(NULL,NULL,"-log_view",&flg);CHKERRQ(ierr);}
970a4af0ceeSJacob Faibussowitsch   if (!flg) {ierr = PetscOptionsGetBool(NULL,NULL,"-viennacl_synchronize",&flg,NULL);CHKERRQ(ierr);}
971a4af0ceeSJacob Faibussowitsch   PetscViennaCLSynchronize = flg;
972a4af0ceeSJacob Faibussowitsch   ierr = PetscViennaCLInit();CHKERRQ(ierr);
973a4af0ceeSJacob Faibussowitsch #endif
974a4af0ceeSJacob Faibussowitsch 
975a4af0ceeSJacob Faibussowitsch   /*
976a4af0ceeSJacob Faibussowitsch      Creates the logging data structures; this is enabled even if logging is not turned on
977a4af0ceeSJacob Faibussowitsch      This is the last thing we do before returning to the user code to prevent having the
978a4af0ceeSJacob Faibussowitsch      logging numbers contaminated by any startup time associated with MPI
979a4af0ceeSJacob Faibussowitsch   */
980a4af0ceeSJacob Faibussowitsch #if defined(PETSC_USE_LOG)
981a4af0ceeSJacob Faibussowitsch   ierr = PetscLogInitialize();CHKERRQ(ierr);
982a4af0ceeSJacob Faibussowitsch #endif
983a4af0ceeSJacob Faibussowitsch 
98485649d77SJunchao Zhang   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
98585649d77SJunchao Zhang 
98685649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
98702101c96SSatish Balay   ierr = PetscInitializeSAWs(ftn ? NULL : help);CHKERRQ(ierr);
98827104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
98927104ee2SJacob Faibussowitsch   ierr = PetscOptionsHasName(NULL,NULL,"-stack_view",&flg);CHKERRQ(ierr);
99027104ee2SJacob Faibussowitsch   if (flg) PetscStackViewSAWs();
99185649d77SJunchao Zhang #endif
99285649d77SJunchao Zhang 
99385649d77SJunchao Zhang   /*
99485649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
99585649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
99685649d77SJunchao Zhang   */
99785649d77SJunchao Zhang   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
99885649d77SJunchao Zhang 
99985649d77SJunchao Zhang   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
10007d3de750SJacob Faibussowitsch   ierr = PetscInfo(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
100185649d77SJunchao Zhang   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
10027d3de750SJacob Faibussowitsch   ierr = PetscInfo(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
100385649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
100485649d77SJunchao Zhang   {
100585649d77SJunchao Zhang     PetscBool omp_view_flag;
100685649d77SJunchao Zhang     char      *threads = getenv("OMP_NUM_THREADS");
100785649d77SJunchao Zhang 
100885649d77SJunchao Zhang     if (threads) {
10097d3de750SJacob Faibussowitsch       ierr = PetscInfo(NULL,"Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
101085649d77SJunchao Zhang       (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
101185649d77SJunchao Zhang     } else {
10122f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt) omp_get_max_threads();
10137d3de750SJacob Faibussowitsch       ierr = PetscInfo(NULL,"Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n",PetscNumOMPThreads);CHKERRQ(ierr);
101485649d77SJunchao Zhang     }
101585649d77SJunchao Zhang     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
101685649d77SJunchao Zhang     ierr = PetscOptionsInt("-omp_num_threads","Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS","None",PetscNumOMPThreads,&PetscNumOMPThreads,&flg);CHKERRQ(ierr);
101785649d77SJunchao Zhang     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
101885649d77SJunchao Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
101985649d77SJunchao Zhang     if (flg) {
10207d3de750SJacob Faibussowitsch       ierr = PetscInfo(NULL,"Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
102185649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
102285649d77SJunchao Zhang     }
102385649d77SJunchao Zhang     if (omp_view_flag) {
10243ca90d2dSJacob Faibussowitsch       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %" PetscInt_FMT "\n",PetscNumOMPThreads);CHKERRQ(ierr);
102585649d77SJunchao Zhang     }
102685649d77SJunchao Zhang   }
102785649d77SJunchao Zhang #endif
102885649d77SJunchao Zhang 
102985649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
103085649d77SJunchao Zhang   /*
103185649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
103285649d77SJunchao Zhang 
103385649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
103485649d77SJunchao Zhang   */
103585649d77SJunchao Zhang   if (!PetscBinaryBigEndian()) {
103685649d77SJunchao Zhang     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRMPI(ierr);
103785649d77SJunchao Zhang   }
103885649d77SJunchao Zhang #endif
103985649d77SJunchao Zhang 
104085649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
104185649d77SJunchao Zhang   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
104285649d77SJunchao Zhang #endif
104385649d77SJunchao Zhang 
104485649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
104585649d77SJunchao Zhang   {
104685649d77SJunchao Zhang     PetscViewer viewer;
104785649d77SJunchao Zhang     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
104885649d77SJunchao Zhang     if (flg) {
104985649d77SJunchao Zhang       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
105085649d77SJunchao Zhang       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
105185649d77SJunchao Zhang     }
105285649d77SJunchao Zhang   }
105385649d77SJunchao Zhang #endif
105485649d77SJunchao Zhang 
105585649d77SJunchao Zhang   flg  = PETSC_TRUE;
105685649d77SJunchao Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
105785649d77SJunchao Zhang   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);}
105885649d77SJunchao Zhang 
105985649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
106085649d77SJunchao Zhang   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
106185649d77SJunchao Zhang   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
106285649d77SJunchao Zhang   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
106385649d77SJunchao Zhang   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
106485649d77SJunchao Zhang #endif
106585649d77SJunchao Zhang 
106685649d77SJunchao Zhang #if defined(__VALGRIND_H)
106785649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND? PETSC_TRUE: PETSC_FALSE;
106885649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
106985649d77SJunchao Zhang   if (PETSC_RUNNING_ON_VALGRIND) {
107085649d77SJunchao Zhang     ierr = 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");CHKERRQ(ierr);
107185649d77SJunchao Zhang     }
107285649d77SJunchao Zhang #endif
107385649d77SJunchao Zhang #endif
107485649d77SJunchao Zhang   /*
107585649d77SJunchao Zhang       Set flag that we are completely initialized
107685649d77SJunchao Zhang   */
107785649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
107885649d77SJunchao Zhang 
107985649d77SJunchao Zhang   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
108085649d77SJunchao Zhang   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
108127104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
108285649d77SJunchao Zhang }
108385649d77SJunchao Zhang 
1084e5c89e4eSSatish Balay /*@C
1085e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1086e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
1087e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1088e5c89e4eSSatish Balay    your program -- usually the very first line!
1089e5c89e4eSSatish Balay 
1090e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
1091e5c89e4eSSatish Balay 
1092e5c89e4eSSatish Balay    Input Parameters:
1093e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1094e5c89e4eSSatish Balay .  args - the command line arguments
1095be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1096be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1097be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1098c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
10990298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1100e5c89e4eSSatish Balay 
110105827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
110205827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
110305827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
110405827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
110505827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1106e5c89e4eSSatish Balay 
1107e5c89e4eSSatish Balay    Options Database Keys:
11087ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
11097ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1110e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
11118a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
11128a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
11138a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
11148a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
11158a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1116bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1117e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1118e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1119e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
112079dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
112179dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
112279dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
1123aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
112492f119d6SBarry 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
112592f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
112692f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1127608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
112892f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1129e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1130e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1131e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1132e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1133e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
11340841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1135e5c89e4eSSatish Balay 
1136c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1137c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1138c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1139c5b5d8d5SVaclav Hapla -  -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet()
1140c5b5d8d5SVaclav Hapla 
1141c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1142c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1143c5b5d8d5SVaclav Hapla    They can be used also in option files.
1144c5b5d8d5SVaclav Hapla 
1145c5b5d8d5SVaclav Hapla    See PetscOptionsMonitorSet() to do monitoring programmatically.
1146c5b5d8d5SVaclav Hapla 
1147e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1148a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1149fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
1150217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1151217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1152495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1153e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
11549a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
115579dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
11569a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1157495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1158f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1159495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
1160495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
1161c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
116287c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
1163c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1164495fc317SBarry Smith 
1165609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
1166e5c89e4eSSatish Balay 
1167ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1168ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1169ffbd1cfbSBarry 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
1170ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1171ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1172ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1173ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1174ffbd1cfbSBarry Smith 
1175e5c89e4eSSatish Balay    Environmental Variables:
1176e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
1177e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
1178e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
11794a971ea4SToby Isaac .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
11804a971ea4SToby Isaac .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1181e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
1182e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
1183e5c89e4eSSatish Balay 
1184e5c89e4eSSatish Balay    Level: beginner
1185e5c89e4eSSatish Balay 
1186e5c89e4eSSatish Balay    Notes:
1187e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
1188e5c89e4eSSatish Balay    it before PetscInitialize().
1189e5c89e4eSSatish Balay 
1190e5c89e4eSSatish Balay    Fortran Version:
1191e5c89e4eSSatish Balay    In Fortran this routine has the format
1192e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
1193e5c89e4eSSatish Balay 
1194e5c89e4eSSatish Balay +  ierr - error return code
1195c5b5d8d5SVaclav Hapla -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
1196c5b5d8d5SVaclav Hapla           Use PETSC_NULL_CHARACTER to not check for code specific file.
1197c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
1198e5c89e4eSSatish Balay 
1199e5c89e4eSSatish Balay    Important Fortran Note:
12000eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
12010298fd71SBarry Smith    null character string; you CANNOT just use NULL as
1202a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
1203e5c89e4eSSatish Balay 
120401cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
120501cb0274SBarry Smith    calling PetscInitialize().
1206e5c89e4eSSatish Balay 
120701cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
1208e5c89e4eSSatish Balay 
1209e5c89e4eSSatish Balay @*/
12107087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
1211e5c89e4eSSatish Balay {
1212e5c89e4eSSatish Balay   PetscErrorCode ierr;
121385649d77SJunchao Zhang   PetscMPIInt    flag;
121485649d77SJunchao Zhang   const char     *prog = "Unknown Name";
1215e5c89e4eSSatish Balay 
121627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
121727104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
1218ffc4695bSBarry Smith   ierr = MPI_Initialized(&flag);CHKERRMPI(ierr);
1219e5c89e4eSSatish Balay   if (!flag) {
1220*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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");
12214dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
12225e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
12235e765c61SJed Brown     {
12245e765c61SJed Brown       PetscMPIInt provided;
1225ffc4695bSBarry Smith       ierr = MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided);CHKERRMPI(ierr);
12265e765c61SJed Brown     }
12275e765c61SJed Brown #else
1228ffc4695bSBarry Smith     ierr = MPI_Init(argc,args);CHKERRMPI(ierr);
12295e765c61SJed Brown #endif
1230e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1231e5c89e4eSSatish Balay   }
12324dfee713SSatish Balay 
123385649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1234e5c89e4eSSatish Balay   if (argc && args) {
1235e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1236e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1237e5c89e4eSSatish Balay   }
123885649d77SJunchao Zhang   ierr = PetscInitialize_Common(prog,file,help,PETSC_FALSE/*C*/,PETSC_FALSE,0);CHKERRQ(ierr);
123927104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
1240e5c89e4eSSatish Balay }
1241e5c89e4eSSatish Balay 
1242a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
124395c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1244ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1245ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
124695c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
12474097062eSBarry Smith #endif
1248e5c89e4eSSatish Balay 
1249008a6e76SBarry Smith /*
1250008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1251008a6e76SBarry Smith */
1252008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1253008a6e76SBarry Smith {
1254008a6e76SBarry Smith   PetscErrorCode ierr;
1255008a6e76SBarry Smith 
1256008a6e76SBarry Smith   PetscFunctionBegin;
1257008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1258ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRMPI(ierr);
1259008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1260ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRMPI(ierr);
1261008a6e76SBarry Smith #endif
1262ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRMPI(ierr);
1263ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRMPI(ierr);
1264008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1265ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRMPI(ierr);
1266ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRMPI(ierr);
1267ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRMPI(ierr);
1268008a6e76SBarry Smith #endif
1269008a6e76SBarry Smith 
1270de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1271ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRMPI(ierr);
1272008a6e76SBarry Smith #endif
1273008a6e76SBarry Smith 
1274ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRMPI(ierr);
1275092991acSStefano Zampini   ierr = MPI_Type_free(&MPIU_REAL_INT);CHKERRMPI(ierr);
1276092991acSStefano Zampini   ierr = MPI_Type_free(&MPIU_SCALAR_INT);CHKERRMPI(ierr);
127740df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1278ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRMPI(ierr);
1279008a6e76SBarry Smith #endif
1280b5a892a1SMatthew G. Knepley   ierr = MPI_Type_free(&MPI_4INT);CHKERRMPI(ierr);
1281b5a892a1SMatthew G. Knepley   ierr = MPI_Type_free(&MPIU_4INT);CHKERRMPI(ierr);
1282ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRMPI(ierr);
1283008a6e76SBarry Smith   PetscFunctionReturn(0);
1284008a6e76SBarry Smith }
1285008a6e76SBarry Smith 
1286a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG)
1287a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
1288a4af0ceeSJacob Faibussowitsch #endif
1289a4af0ceeSJacob Faibussowitsch 
1290e5c89e4eSSatish Balay /*@C
1291e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1292e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1293e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1294e5c89e4eSSatish Balay 
1295e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1296e5c89e4eSSatish Balay 
1297e5c89e4eSSatish Balay    Options Database Keys:
129826a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1299e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
13007eb1d149SBarry 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
1301e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
130292f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1303e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
130492f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1305e5c89e4eSSatish Balay 
1306e5c89e4eSSatish Balay    Level: beginner
1307e5c89e4eSSatish Balay 
1308e5c89e4eSSatish Balay    Note:
1309e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1310e5c89e4eSSatish Balay 
131188c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1312e5c89e4eSSatish Balay @*/
13137087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1314e5c89e4eSSatish Balay {
1315e5c89e4eSSatish Balay   PetscErrorCode ierr;
13164bb5149bSJed Brown   PetscMPIInt    rank;
1317a8d2bbe5SBarry Smith   PetscInt       nopt;
13182bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1319dff31646SBarry Smith   PetscBool      flg;
132010463e74SBarry Smith #if defined(PETSC_USE_LOG)
132110463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
132210463e74SBarry Smith #endif
1323e5c89e4eSSatish Balay 
13243cbbc5ffSBarry Smith   PetscFunctionBegin;
132527104ee2SJacob Faibussowitsch   if (PetscUnlikely(!PetscInitializeCalled)) {
13261724198aSStefano Zampini     fprintf(PETSC_STDOUT,"PetscInitialize() must be called before PetscFinalize()\n");
13271724198aSStefano Zampini     ierr = PetscStackView(PETSC_STDOUT);CHKERRQ(ierr);
13281724198aSStefano Zampini     PetscStackClearTop;
132927104ee2SJacob Faibussowitsch     return PETSC_ERR_ARG_WRONGSTATE;
133027104ee2SJacob Faibussowitsch   }
13310298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1332b022a5c1SBarry Smith 
1333ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
1334a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
133522580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1336a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1337a56f64adSBarry Smith #endif
1338c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1339dff31646SBarry Smith   if (flg) {
13401f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
13411f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
13421f817a21SBarry Smith 
1343589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL);CHKERRQ(ierr);
13441f817a21SBarry Smith     if (filename[0]) {
13451f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
13461f817a21SBarry Smith     }
1347dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1348dff31646SBarry Smith     cits[0] = 0;
1349dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
13501f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
13511f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
13521f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
13531f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
13541f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1355dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1356dff31646SBarry Smith   }
1357dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1358dff31646SBarry Smith 
1359c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
136004102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
136104102261SBarry Smith   {
136204102261SBarry Smith     PetscInt nmax = 2;
136304102261SBarry Smith     char     **buffs;
136404102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1365c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
136604102261SBarry Smith     if (flg1) {
1367*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
136804102261SBarry Smith       if (nmax == 1) {
136904102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
137004102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
137104102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
137204102261SBarry Smith       }
137304102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
137404102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
137504102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
137604102261SBarry Smith     }
137704102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
137804102261SBarry Smith   }
1379763ec7b1SBarry Smith   {
1380763ec7b1SBarry Smith     PetscInt nmax = 2;
1381763ec7b1SBarry Smith     char     **buffs;
1382763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1383763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1384763ec7b1SBarry Smith     if (flg1) {
1385*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1386763ec7b1SBarry Smith       if (nmax == 1) {
1387763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1388763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1389763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1390763ec7b1SBarry Smith       }
1391763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1392763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1393763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1394763ec7b1SBarry Smith     }
1395763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1396763ec7b1SBarry Smith   }
139704102261SBarry Smith #endif
139804102261SBarry Smith 
13992d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
14002d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
14012d53ad75SBarry Smith #endif
14022d53ad75SBarry Smith 
1403e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1404dff31646SBarry Smith   flg = PETSC_FALSE;
1405c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1406d5649816SBarry Smith   if (flg) {
1407e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1408d5649816SBarry Smith   }
1409d5649816SBarry Smith #endif
1410d5649816SBarry Smith 
1411681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1412681455b2SBarry Smith   flg1 = PETSC_FALSE;
1413c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1414681455b2SBarry Smith   if (flg1) {
1415681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1416681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1417681455b2SBarry Smith   }
1418681455b2SBarry Smith #endif
1419681455b2SBarry Smith 
142067584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1421c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1422e5c89e4eSSatish Balay   if (!flg2) {
142390d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1424c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1425e5c89e4eSSatish Balay   }
1426e5c89e4eSSatish Balay   if (flg2) {
14270841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1428e5c89e4eSSatish Balay   }
142967584ceeSBarry Smith #endif
1430e5c89e4eSSatish Balay 
1431e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
143290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1433c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1434e5c89e4eSSatish Balay   if (flg1) {
1435e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1436ffc4695bSBarry Smith     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRMPI(ierr);
1437e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1438e5c89e4eSSatish Balay   }
1439e5c89e4eSSatish Balay #endif
1440e5c89e4eSSatish Balay 
1441e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1442e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1443e5c89e4eSSatish Balay   mname[0] = 0;
1444589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1445e5c89e4eSSatish Balay   if (flg1) {
1446e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1447e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1448e5c89e4eSSatish Balay   }
1449e5c89e4eSSatish Balay #endif
1450356e5837SBarry Smith #endif
1451a297a907SKarl Rupp 
1452dd710f27SBarry Smith   /*
1453dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1454dd710f27SBarry Smith   */
1455dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1456dd710f27SBarry Smith 
1457356e5837SBarry Smith #if defined(PETSC_USE_LOG)
145887c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1459f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
146087c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
146187c3beb6SLisandro Dalcin 
1462356e5837SBarry Smith   mname[0] = 0;
1463589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1464e5c89e4eSSatish Balay   if (flg1) {
146591eabc43SBarry Smith     PetscViewer viewer;
146620a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
146791eabc43SBarry Smith     if (mname[0]) {
146891eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
146991eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
14706bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
147133f85c2fSBarry Smith     } else {
147233f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
14739a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
147433f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
14759a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
147633f85c2fSBarry Smith     }
1477e5c89e4eSSatish Balay   }
1478a297a907SKarl Rupp 
1479dd710f27SBarry Smith   /*
1480dd710f27SBarry Smith      Free any objects created by the last block of code.
1481dd710f27SBarry Smith   */
1482dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1483dd710f27SBarry Smith 
1484dd710f27SBarry Smith   mname[0] = 0;
1485589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1486589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2);CHKERRQ(ierr);
14877ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1488e5c89e4eSSatish Balay #endif
148910463e74SBarry Smith 
149090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1491c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1492e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
149390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1494c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1495e5c89e4eSSatish Balay   if (flg1) {
1496e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1497e5c89e4eSSatish Balay   }
149890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
149990d69ab7SBarry Smith   flg2 = PETSC_FALSE;
15008bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1501c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1502c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1503c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1504e4c476e2SSatish Balay 
1505e5c89e4eSSatish Balay   if (flg2) {
1506be56827dSJed Brown     PetscViewer viewer;
150702ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
150802ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1509c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1510be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1511e5c89e4eSSatish Balay   }
1512e5c89e4eSSatish Balay 
1513e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1514c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1515c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1516e5c89e4eSSatish Balay 
151733fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1518c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
15199245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1520e5c89e4eSSatish Balay   if (flg3) {
15213de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1522be56827dSJed Brown       PetscViewer viewer;
152302ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
152402ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1525c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1526be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1527e5c89e4eSSatish Balay     }
15283de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
15293de2bfdfSBarry Smith     if (nopt) {
15303de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
15313de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
15323de2bfdfSBarry Smith       if (nopt == 1) {
1533e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1534e5c89e4eSSatish Balay       } else {
15353ca90d2dSJacob Faibussowitsch         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %" PetscInt_FMT " unused database options. They are:\n",nopt);CHKERRQ(ierr);
1536e5c89e4eSSatish Balay       }
15373de2bfdfSBarry Smith     } else if (flg3 && flg1) {
15383de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1539df12ba86SBarry Smith     }
1540c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1541e5c89e4eSSatish Balay   }
1542e5c89e4eSSatish Balay 
1543e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1544d45a07a7SBarry Smith   if (!PetscGlobalRank) {
154587f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
154616ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1547d45a07a7SBarry Smith   }
1548ec957eceSBarry Smith #endif
1549ec957eceSBarry Smith 
15504097062eSBarry Smith #if defined(PETSC_USE_LOG)
155110463e74SBarry Smith   /*
1552dbc8283eSBarry Smith        List all objects the user may have forgot to free
15532eff7a51SBarry Smith   */
155405df10baSBarry Smith   if (PetscObjectsLog) {
1555c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1556a64a8e02SBarry Smith     if (flg1) {
1557a64a8e02SBarry Smith       MPI_Comm local_comm;
15587eb1d149SBarry Smith       char     string[64];
1559a64a8e02SBarry Smith 
1560589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL);CHKERRQ(ierr);
1561ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(ierr);
1562a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
15637eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1564a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1565ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(ierr);
15660a1571b3SBarry Smith     }
156705df10baSBarry Smith   }
15684097062eSBarry Smith #endif
15694097062eSBarry Smith 
15704097062eSBarry Smith #if defined(PETSC_USE_LOG)
1571dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1572dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1573a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
15744097062eSBarry Smith #endif
15752eff7a51SBarry Smith 
157633f85c2fSBarry Smith   /*
157733f85c2fSBarry Smith      Destroy any packages that registered a finalize
157833f85c2fSBarry Smith   */
157933f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
158033f85c2fSBarry Smith 
1581101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1582fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1583101409b8SToby Isaac #endif
1584101409b8SToby Isaac 
158533f85c2fSBarry Smith   /*
158648dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
158748dd1dffSBarry Smith 
158837e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
158948dd1dffSBarry Smith   */
159037e93019SBarry Smith 
15914028d114SSatish Balay   if (petsc_history) {
1592f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
159302c9f0b5SLisandro Dalcin     petsc_history = NULL;
1594e5c89e4eSSatish Balay   }
15959de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1596e94e781bSJacob Faibussowitsch   ierr = PetscInfoDestroy();CHKERRQ(ierr);
1597e5c89e4eSSatish Balay 
159867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
159992f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1600e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
160192f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1602e5c89e4eSSatish Balay     FILE *fd;
1603ed9cf6e9SBarry Smith     int  err;
1604e5c89e4eSSatish Balay 
1605dc92acbaSJed Brown     flg2 = PETSC_FALSE;
160692f119d6SBarry Smith     flg3 = PETSC_FALSE;
1607cf9c20a2SJed Brown     if (PetscDefined(USE_DEBUG)) {ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);}
160892f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
160992f119d6SBarry Smith     fname[0] = 0;
1610589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1611e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1612e5c89e4eSSatish Balay 
1613589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1614*2c71b3e2SJacob Faibussowitsch       fd   = fopen(sname,"w"); PetscCheckFalse(!fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1615e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1616ed9cf6e9SBarry Smith       err  = fclose(fd);
1617*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
161892f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1619e5c89e4eSSatish Balay       MPI_Comm local_comm;
1620e5c89e4eSSatish Balay 
1621ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(ierr);
1622e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1623e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1624e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1625ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(ierr);
1626e5c89e4eSSatish Balay     }
1627e5c89e4eSSatish Balay     fname[0] = 0;
1628589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1629e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1630e5c89e4eSSatish Balay 
1631589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1632*2c71b3e2SJacob Faibussowitsch       fd   = fopen(sname,"w"); PetscCheckFalse(!fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
163392f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1634ed9cf6e9SBarry Smith       err  = fclose(fd);
1635*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
163692f119d6SBarry Smith     } else if (flg1) {
163792f119d6SBarry Smith       MPI_Comm local_comm;
163892f119d6SBarry Smith 
1639ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(ierr);
164092f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
164192f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
164292f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1643ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(ierr);
1644e5c89e4eSSatish Balay     }
1645e5c89e4eSSatish Balay   }
164667584ceeSBarry Smith #endif
164720e2c332SMatthew G. Knepley 
16485486ca60SMatthew G. Knepley   /*
16495486ca60SMatthew G. Knepley      Close any open dynamic libraries
16505486ca60SMatthew G. Knepley   */
16515486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
16525486ca60SMatthew G. Knepley 
1653e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16544416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1655e5c89e4eSSatish Balay 
1656e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
165702c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1658e5c89e4eSSatish Balay 
1659c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1660c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
1661c2b86a48SJunchao Zhang     ierr = PetscKokkosFinalize_Private();CHKERRQ(ierr);
1662c2b86a48SJunchao Zhang     PetscBeganKokkos = PETSC_FALSE;
166345639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1664c2b86a48SJunchao Zhang   }
1665c2b86a48SJunchao Zhang #endif
1666c2b86a48SJunchao Zhang 
166771438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
166871438e86SJunchao Zhang   if (PetscBeganNvshmem) {
166971438e86SJunchao Zhang     ierr = PetscNvshmemFinalize();CHKERRQ(ierr);
167071438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
167171438e86SJunchao Zhang   }
167271438e86SJunchao Zhang #endif
1673a0e72f99SJunchao Zhang 
1674008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1675e5c89e4eSSatish Balay 
1676dbc8283eSBarry Smith   /*
1677efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1678efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1679efb80d3cSBarry Smith 
1680efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1681efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1682dbc8283eSBarry Smith  */
1683b770b1f6SSatish Balay   {
1684dbc8283eSBarry Smith     PetscCommCounter *counter;
1685dbc8283eSBarry Smith     PetscMPIInt      flg;
1686dbc8283eSBarry Smith     MPI_Comm         icomm;
1687265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
1688ffc4695bSBarry Smith     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRMPI(ierr);
1689dbc8283eSBarry Smith     if (flg) {
1690265f3f35SJed Brown       icomm = ucomm.comm;
1691ffc4695bSBarry Smith       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRMPI(ierr);
1692*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1693dbc8283eSBarry Smith 
1694ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1695ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRMPI(ierr);
1696ffc4695bSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRMPI(ierr);
1697dbc8283eSBarry Smith     }
1698ffc4695bSBarry Smith     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRMPI(ierr);
1699dbc8283eSBarry Smith     if (flg) {
1700265f3f35SJed Brown       icomm = ucomm.comm;
1701ffc4695bSBarry Smith       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRMPI(ierr);
1702*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1703dbc8283eSBarry Smith 
1704ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1705ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRMPI(ierr);
1706ffc4695bSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRMPI(ierr);
1707dbc8283eSBarry Smith     }
1708b770b1f6SSatish Balay   }
1709dbc8283eSBarry Smith 
1710ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRMPI(ierr);
1711ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1712ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRMPI(ierr);
1713ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRMPI(ierr);
1714480cf27aSJed Brown 
17155ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
17165ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
17175ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1718ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1719ef19f930SBarry Smith 
1720e5c89e4eSSatish Balay   if (PetscBeganMPI) {
172199608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
172299b1327fSBarry Smith     PetscMPIInt flag;
1723ffc4695bSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRMPI(ierr);
1724*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
172599608316SBarry Smith #endif
1726ffc4695bSBarry Smith     ierr = MPI_Finalize();CHKERRMPI(ierr);
1727e5c89e4eSSatish Balay   }
1728e5c89e4eSSatish Balay /*
1729e5c89e4eSSatish Balay 
1730e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1731e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1732e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1733e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1734e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17350e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1736e5c89e4eSSatish Balay    memory was not freed.
1737e5c89e4eSSatish Balay 
1738e5c89e4eSSatish Balay */
17391d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
174027104ee2SJacob Faibussowitsch   ierr = PetscStackReset();CHKERRQ(ierr);
1741a297a907SKarl Rupp 
17428ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_FALSE;
1743e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1744e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
174556883f60SBarry Smith #if defined(PETSC_USE_GCOV)
174656883f60SBarry Smith   /*
174756883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
174856883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
174956883f60SBarry Smith    */
175056883f60SBarry Smith   __gcov_flush();
175156883f60SBarry Smith #endif
17521724198aSStefano Zampini   /* To match PetscFunctionBegin() at the beginning of this function */
17531724198aSStefano Zampini   PetscStackClearTop;
175427104ee2SJacob Faibussowitsch   return 0;
1755e5c89e4eSSatish Balay }
1756e5c89e4eSSatish Balay 
175743db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
17588cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
175943db4dbbSBarry Smith {
176043db4dbbSBarry Smith   if (*a == *b) return 1;
176143db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
176243db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
176343db4dbbSBarry Smith   return 0;
176443db4dbbSBarry Smith }
1765a70650f6SBarry Smith #endif
176643db4dbbSBarry Smith 
176743db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
17688cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
176943db4dbbSBarry Smith {
177043db4dbbSBarry Smith   if (*a == *b) return 1;
177143db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
177243db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
177343db4dbbSBarry Smith   return 0;
177443db4dbbSBarry Smith }
177543db4dbbSBarry Smith #endif
1776