xref: /petsc/src/sys/objects/pinit.c (revision 2f613bf53f46f9356e00a2ca2bd69453be72fc31)
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 
16a0e72f99SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
17030f984aSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
189ffd0706SHong Zhang PETSC_EXTERN cudaEvent_t petsc_gputimer_begin;
199ffd0706SHong Zhang PETSC_EXTERN cudaEvent_t petsc_gputimer_end;
20a0e72f99SJunchao Zhang #endif
21a0e72f99SJunchao Zhang 
22a0e72f99SJunchao Zhang #if defined(PETSC_HAVE_HIP)
23030f984aSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
249ffd0706SHong Zhang PETSC_EXTERN hipEvent_t petsc_gputimer_begin;
259ffd0706SHong Zhang PETSC_EXTERN hipEvent_t petsc_gputimer_end;
26a0e72f99SJunchao Zhang #endif
27a0e72f99SJunchao Zhang 
2856883f60SBarry Smith #if defined(PETSC_USE_GCOV)
2956883f60SBarry Smith EXTERN_C_BEGIN
3056883f60SBarry Smith void  __gcov_flush(void);
3156883f60SBarry Smith EXTERN_C_END
3256883f60SBarry Smith #endif
338101f56cSMatthew Knepley 
3427104ee2SJacob Faibussowitsch #if PetscDefined(USE_LOG)
3595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
36a9f03627SSatish Balay #endif
37f2d66bcaSShri Abhyankar 
382d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
3995c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
402d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
412d53ad75SBarry Smith #endif
422d53ad75SBarry Smith 
4327104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
4416ad0300SBarry Smith #include <petscviewersaws.h>
45a6790183SBarry Smith #endif
46eb02082bSJunchao Zhang 
47e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
48e5c89e4eSSatish Balay 
4995c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
50e5c89e4eSSatish Balay 
5195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
5295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
5395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
5495c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
5595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
5695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
570069ddf5SShri Abhyankar 
586de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
59e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
6027104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD)
616de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
626de5d289SStefano Zampini #else
636de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
646de5d289SStefano Zampini #endif
65e5c89e4eSSatish Balay 
66480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
67480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
68480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
695f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
70480cf27aSJed Brown 
71e5c89e4eSSatish Balay /*
72e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
73e5c89e4eSSatish Balay */
7402c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
7502c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
76e5c89e4eSSatish Balay 
77ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
78ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
790f8e0872SSatish Balay 
80a2f94806SJed Brown PetscInt PetscHotRegionDepth;
81a2f94806SJed Brown 
826edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;
836edef35eSSatish Balay 
84b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
85b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
86b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
87b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
88b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
89b22622e2STadeu Manoel #endif
90b22622e2STadeu Manoel 
91e5c89e4eSSatish Balay /*
92945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
9372a42c3cSBarry Smith 
9472a42c3cSBarry Smith    Collective
9572a42c3cSBarry Smith 
9672a42c3cSBarry Smith    Level: advanced
9772a42c3cSBarry Smith 
9895452b02SPatrick Sanan     Notes:
99a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
1000f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
101a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
1020f11a792SBarry Smith 
103a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
1041ea5a559SBarry Smith 
10572a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
1060f11a792SBarry Smith */
107945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
10872a42c3cSBarry Smith {
10972a42c3cSBarry Smith   PetscErrorCode ierr;
11072a42c3cSBarry Smith   int            myargc   = argc;
11172a42c3cSBarry Smith   char           **myargs = args;
11272a42c3cSBarry Smith 
11372a42c3cSBarry Smith   PetscFunctionBegin;
11427104ee2SJacob Faibussowitsch   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) PetscFunctionReturn(ierr);
1151ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
116df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
11727104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
11872a42c3cSBarry Smith }
11972a42c3cSBarry Smith 
120f0865b08SBarry Smith /*
121a56f64adSBarry Smith       Used by Julia interface to get communicator
122f0865b08SBarry Smith */
123945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
124f0865b08SBarry Smith {
125f0865b08SBarry Smith   PetscFunctionBegin;
12627104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidPointer(comm,1);
127f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
128f0865b08SBarry Smith   PetscFunctionReturn(0);
129f0865b08SBarry Smith }
130f0865b08SBarry Smith 
131e5c89e4eSSatish Balay /*@C
132e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
133e5c89e4eSSatish Balay         the command line arguments.
134e5c89e4eSSatish Balay 
135e5c89e4eSSatish Balay    Collective
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay    Level: advanced
138e5c89e4eSSatish Balay 
139e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
140e5c89e4eSSatish Balay @*/
1417087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
142e5c89e4eSSatish Balay {
143e5c89e4eSSatish Balay   PetscErrorCode ierr;
144e5c89e4eSSatish Balay   int            argc   = 0;
14502c9f0b5SLisandro Dalcin   char           **args = NULL;
146e5c89e4eSSatish Balay 
147e5c89e4eSSatish Balay   PetscFunctionBegin;
1480298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
149e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
150e5c89e4eSSatish Balay }
151e5c89e4eSSatish Balay 
152e5c89e4eSSatish Balay /*@
153e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
154e5c89e4eSSatish Balay 
15593b6d2d1SJed Brown    Level: beginner
156e5c89e4eSSatish Balay 
157e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
158e5c89e4eSSatish Balay @*/
1597087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized)
160e5c89e4eSSatish Balay {
16127104ee2SJacob Faibussowitsch   PetscFunctionBegin;
16227104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized,1);
163e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
16427104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
165e5c89e4eSSatish Balay }
166e5c89e4eSSatish Balay 
167e5c89e4eSSatish Balay /*@
168e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
169e5c89e4eSSatish Balay 
170e5c89e4eSSatish Balay    Level: developer
171e5c89e4eSSatish Balay 
172e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
173e5c89e4eSSatish Balay @*/
1747087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
175e5c89e4eSSatish Balay {
17627104ee2SJacob Faibussowitsch   PetscFunctionBegin;
17727104ee2SJacob Faibussowitsch   if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized,1);
178e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
17927104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
180e5c89e4eSSatish Balay }
181e5c89e4eSSatish Balay 
18257171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []);
183e5c89e4eSSatish Balay 
184e5c89e4eSSatish Balay /*
185e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
186e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
187e5c89e4eSSatish Balay */
188367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
189e5c89e4eSSatish Balay 
190367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
191e5c89e4eSSatish Balay {
192e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
193e5c89e4eSSatish Balay 
194e5c89e4eSSatish Balay   PetscFunctionBegin;
195e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
196e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
19741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
198e5c89e4eSSatish Balay   }
199e5c89e4eSSatish Balay 
200e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
201e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
202e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
203e5c89e4eSSatish Balay   }
204812af9f3SBarry Smith   PetscFunctionReturnVoid();
205e5c89e4eSSatish Balay }
206e5c89e4eSSatish Balay 
207e5c89e4eSSatish Balay /*
208e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
209e5c89e4eSSatish Balay sum of the second entry.
210b693b147SBarry Smith 
21176f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
212367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
213b693b147SBarry Smith there would be no place to store the both needed results.
214e5c89e4eSSatish Balay */
21576ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
216e5c89e4eSSatish Balay {
217e5c89e4eSSatish Balay   PetscErrorCode ierr;
218e5c89e4eSSatish Balay 
219e5c89e4eSSatish Balay   PetscFunctionBegin;
220d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
221d6e4c47cSJed Brown   {
222d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
223ffc4695bSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRMPI(ierr);
224d6e4c47cSJed Brown     *max = work.max;
225d6e4c47cSJed Brown     *sum = work.sum;
226d6e4c47cSJed Brown   }
227d6e4c47cSJed Brown #else
228d6e4c47cSJed Brown   {
229d6e4c47cSJed Brown     PetscMPIInt    size,rank;
230d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
231ffc4695bSBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
232ffc4695bSBarry Smith     ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
233785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
234820f2d46SBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRMPI(ierr);
2356ac3741eSJed Brown     *max = work[rank].max;
2366ac3741eSJed Brown     *sum = work[rank].sum;
237e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
238d6e4c47cSJed Brown   }
239d6e4c47cSJed Brown #endif
240e5c89e4eSSatish Balay   PetscFunctionReturn(0);
241e5c89e4eSSatish Balay }
242e5c89e4eSSatish Balay 
243e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
244e5c89e4eSSatish Balay 
245de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
24606a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
247e5c89e4eSSatish Balay 
248092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
249e5c89e4eSSatish Balay {
250e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
251e5c89e4eSSatish Balay 
252e5c89e4eSSatish Balay   PetscFunctionBegin;
2537c2de775SJed Brown   if (*datatype == MPIU_REAL) {
254e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
255a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2567c2de775SJed Brown   }
2577c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
258c5481aeeSJose E. Roman   else if (*datatype == MPIU_COMPLEX) {
2597c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
260a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2617c2de775SJed Brown   }
2627c2de775SJed Brown #endif
2637c2de775SJed Brown   else {
2647c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
26541e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
266e2e03761SBarry Smith   }
267812af9f3SBarry Smith   PetscFunctionReturnVoid();
268e5c89e4eSSatish Balay }
269e5c89e4eSSatish Balay #endif
270e5c89e4eSSatish Balay 
271570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
272d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
273d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
274d9822059SBarry Smith 
275092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
276d9822059SBarry Smith {
277d9822059SBarry Smith   PetscInt i,count = *cnt;
278d9822059SBarry Smith 
279d9822059SBarry Smith   PetscFunctionBegin;
2807c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2818c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
282a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2837c2de775SJed Brown   }
2847c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2857c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2867c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2877c2de775SJed Brown     for (i=0; i<count; i++) {
2887c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2897c2de775SJed Brown     }
2907c2de775SJed Brown   }
2917c2de775SJed Brown #endif
2927c2de775SJed Brown   else {
2937c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
29441e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2958c764dc5SJose Roman   }
296d9822059SBarry Smith   PetscFunctionReturnVoid();
297d9822059SBarry Smith }
298d9822059SBarry Smith 
299092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
300d9822059SBarry Smith {
301d9822059SBarry Smith   PetscInt    i,count = *cnt;
302d9822059SBarry Smith 
303d9822059SBarry Smith   PetscFunctionBegin;
3047c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3058c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
306a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
3077c2de775SJed Brown   }
3087c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
3097c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3107c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
3117c2de775SJed Brown     for (i=0; i<count; i++) {
3127c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3137c2de775SJed Brown     }
3147c2de775SJed Brown   }
3157c2de775SJed Brown #endif
3167c2de775SJed Brown   else {
3178c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
31841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
3198c764dc5SJose Roman   }
320d9822059SBarry Smith   PetscFunctionReturnVoid();
321d9822059SBarry Smith }
322d9822059SBarry Smith #endif
323d9822059SBarry Smith 
324480cf27aSJed Brown /*
325480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
326480cf27aSJed Brown 
327ff0e51ddSBarry 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.
328480cf27aSJed Brown 
32912801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
330480cf27aSJed Brown 
331480cf27aSJed Brown */
33233779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
333480cf27aSJed Brown {
334480cf27aSJed Brown   PetscErrorCode   ierr;
33505342407SJunchao Zhang   PetscCommCounter *counter=(PetscCommCounter*)count_val;
336480cf27aSJed Brown 
337480cf27aSJed Brown   PetscFunctionBegin;
33802c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
33905342407SJunchao Zhang   ierr = PetscFree(counter->iflags);CHKERRMPI(ierr);
34005342407SJunchao Zhang   ierr = PetscFree(counter);CHKERRMPI(ierr);
341480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
342480cf27aSJed Brown }
343480cf27aSJed Brown 
344480cf27aSJed Brown /*
34547435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
346da3039f7SJed Brown   calls MPI_Comm_free().
347da3039f7SJed Brown 
348da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
349480cf27aSJed Brown 
350ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
351480cf27aSJed Brown 
35212801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
353480cf27aSJed Brown 
354480cf27aSJed Brown */
35533779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
356480cf27aSJed Brown {
357480cf27aSJed Brown   PetscErrorCode                    ierr;
35833779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
359480cf27aSJed Brown 
360480cf27aSJed Brown   PetscFunctionBegin;
36112801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
362ec4fadc2SJed Brown   icomm.ptr = attr_val;
36376bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
36433779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
36533779a13SJunchao Zhang     PetscMPIInt flg;
36633779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
36747435625SJed Brown     ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
36833779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
36933779a13SJunchao 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");
37033779a13SJunchao Zhang   }
37147435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
37233779a13SJunchao Zhang   ierr = PetscInfo2(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);CHKERRMPI(ierr);
373da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
374b89831e5SBarry Smith }
375da3039f7SJed Brown 
376da3039f7SJed Brown /*
37733779a13SJunchao 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.
378da3039f7SJed Brown  */
37933779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
380da3039f7SJed Brown {
381da3039f7SJed Brown   PetscErrorCode ierr;
382da3039f7SJed Brown 
383da3039f7SJed Brown   PetscFunctionBegin;
38402c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
385480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
386480cf27aSJed Brown }
387480cf27aSJed Brown 
38833779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
38942218b76SBarry Smith 
390951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3918cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3928cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3938cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
394e39fd77fSBarry Smith #endif
395e39fd77fSBarry Smith 
3960f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
39712801b39SBarry Smith 
398eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
399eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
4006ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
40102c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
402dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
403e5c89e4eSSatish Balay 
404dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
405051e4cf2SJed Brown {
406051e4cf2SJed Brown   PetscErrorCode ierr;
407051e4cf2SJed Brown 
408051e4cf2SJed Brown   PetscFunctionBegin;
409051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
410a1601671SBarry 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);
411051e4cf2SJed 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);
412051e4cf2SJed Brown   PetscFunctionReturn(0);
413051e4cf2SJed Brown }
414e5c89e4eSSatish Balay 
4152d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4162d747510SLisandro Dalcin 
4172d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4182d747510SLisandro Dalcin {
4192d747510SLisandro Dalcin   PetscErrorCode ierr;
4202d747510SLisandro Dalcin 
4212d747510SLisandro Dalcin   PetscFunctionBegin;
422589a23caSBarry Smith   ierr  = PetscStrncpy(programname,name,sizeof(programname));CHKERRQ(ierr);
4232d747510SLisandro Dalcin   PetscFunctionReturn(0);
4242d747510SLisandro Dalcin }
4252d747510SLisandro Dalcin 
4262d747510SLisandro Dalcin /*@C
4272d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4282d747510SLisandro Dalcin 
4292d747510SLisandro Dalcin     Not Collective
4302d747510SLisandro Dalcin 
4312d747510SLisandro Dalcin     Input Parameter:
4322d747510SLisandro Dalcin .   len - length of the string name
4332d747510SLisandro Dalcin 
4342d747510SLisandro Dalcin     Output Parameter:
4352d747510SLisandro Dalcin .   name - the name of the running program
4362d747510SLisandro Dalcin 
4372d747510SLisandro Dalcin    Level: advanced
4382d747510SLisandro Dalcin 
4392d747510SLisandro Dalcin     Notes:
4402d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4412d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4422d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4432d747510SLisandro Dalcin @*/
4442d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4452d747510SLisandro Dalcin {
4462d747510SLisandro Dalcin   PetscErrorCode ierr;
4472d747510SLisandro Dalcin 
4482d747510SLisandro Dalcin   PetscFunctionBegin;
4492d747510SLisandro Dalcin   ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4502d747510SLisandro Dalcin   PetscFunctionReturn(0);
4512d747510SLisandro Dalcin }
4522d747510SLisandro Dalcin 
453e5c89e4eSSatish Balay /*@C
454e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
455e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
456e5c89e4eSSatish Balay 
457e5c89e4eSSatish Balay    Not Collective
458e5c89e4eSSatish Balay 
459e5c89e4eSSatish Balay    Output Parameters:
460e5c89e4eSSatish Balay +  argc - count of number of command line arguments
461e5c89e4eSSatish Balay -  args - the command line arguments
462e5c89e4eSSatish Balay 
463e5c89e4eSSatish Balay    Level: intermediate
464e5c89e4eSSatish Balay 
465e5c89e4eSSatish Balay    Notes:
466e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
467e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
468e5c89e4eSSatish Balay 
469f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
470f177e3b1SBarry Smith 
471793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
472e5c89e4eSSatish Balay 
473e5c89e4eSSatish Balay @*/
4747087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
475e5c89e4eSSatish Balay {
476e5c89e4eSSatish Balay   PetscFunctionBegin;
47727104ee2SJacob Faibussowitsch   if (PetscUnlikely(!PetscInitializeCalled && PetscFinalizeCalled)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
478e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
479e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
480e5c89e4eSSatish Balay   PetscFunctionReturn(0);
481e5c89e4eSSatish Balay }
482e5c89e4eSSatish Balay 
483793721a6SBarry Smith /*@C
484793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
485793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
486793721a6SBarry Smith 
487793721a6SBarry Smith    Not Collective
488793721a6SBarry Smith 
489793721a6SBarry Smith    Output Parameters:
490793721a6SBarry Smith .  args - the command line arguments
491793721a6SBarry Smith 
492793721a6SBarry Smith    Level: intermediate
493793721a6SBarry Smith 
494793721a6SBarry Smith    Notes:
495793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
496793721a6SBarry Smith 
497793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
498793721a6SBarry Smith 
499793721a6SBarry Smith @*/
5007087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
501793721a6SBarry Smith {
502793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
503793721a6SBarry Smith   PetscErrorCode ierr;
504793721a6SBarry Smith 
505793721a6SBarry Smith   PetscFunctionBegin;
50627104ee2SJacob Faibussowitsch   if (PetscUnlikely(!PetscInitializeCalled && PetscFinalizeCalled)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
5072d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
508785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
509793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
510793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
511793721a6SBarry Smith   }
5122d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
513793721a6SBarry Smith   PetscFunctionReturn(0);
514793721a6SBarry Smith }
515793721a6SBarry Smith 
516793721a6SBarry Smith /*@C
517793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
518793721a6SBarry Smith 
519793721a6SBarry Smith    Not Collective
520793721a6SBarry Smith 
521793721a6SBarry Smith    Output Parameters:
522793721a6SBarry Smith .  args - the command line arguments
523793721a6SBarry Smith 
524793721a6SBarry Smith    Level: intermediate
525793721a6SBarry Smith 
526793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
527793721a6SBarry Smith 
528793721a6SBarry Smith @*/
5297087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
530793721a6SBarry Smith {
531793721a6SBarry Smith   PetscInt       i = 0;
532793721a6SBarry Smith   PetscErrorCode ierr;
533793721a6SBarry Smith 
534793721a6SBarry Smith   PetscFunctionBegin;
535a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
536793721a6SBarry Smith   while (args[i]) {
537793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
538793721a6SBarry Smith     i++;
539793721a6SBarry Smith   }
540793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
541793721a6SBarry Smith   PetscFunctionReturn(0);
542793721a6SBarry Smith }
543793721a6SBarry Smith 
54427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS)
54530befbd2SBarry Smith #include <petscconfiginfo.h>
54630befbd2SBarry Smith 
54795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
54811525c0dSBarry Smith {
54927104ee2SJacob Faibussowitsch   PetscFunctionBegin;
55011525c0dSBarry Smith   if (!PetscGlobalRank) {
55130befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
55211525c0dSBarry Smith     int            port;
553ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
55411525c0dSBarry Smith     size_t         applinelen,introlen;
55511525c0dSBarry Smith     PetscErrorCode ierr;
556ffbd1cfbSBarry Smith     char           sawsurl[256];
55711525c0dSBarry Smith 
558c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
55911525c0dSBarry Smith     if (flg) {
56011525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
56111525c0dSBarry Smith 
562589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL);CHKERRQ(ierr);
56311525c0dSBarry Smith       if (sawslog[0]) {
56411525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
56511525c0dSBarry Smith       } else {
56611525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
56711525c0dSBarry Smith       }
56811525c0dSBarry Smith     }
569589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg);CHKERRQ(ierr);
57011525c0dSBarry Smith     if (flg) {
57111525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
57211525c0dSBarry Smith     }
573c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
574ffbd1cfbSBarry Smith     if (selectport) {
575ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
576ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
577ffbd1cfbSBarry Smith     } else {
578c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
57911525c0dSBarry Smith       if (flg) {
58011525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
58111525c0dSBarry Smith       }
582ffbd1cfbSBarry Smith     }
583589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg);CHKERRQ(ierr);
58411525c0dSBarry Smith     if (flg) {
585*2f613bf5SBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
58611525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5879c1e0ce8SBarry Smith     } else {
588c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5899c1e0ce8SBarry Smith       if (flg) {
590589a23caSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root));CHKERRQ(ierr);
591*2f613bf5SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
5929c1e0ce8SBarry Smith       }
59311525c0dSBarry Smith     }
594c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
59511525c0dSBarry Smith     if (flg2) {
59611525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
59711525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
598589a23caSBarry Smith       ierr = PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root);CHKERRQ(ierr);
59911525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
60011525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
601*2f613bf5SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());
60211525c0dSBarry Smith     }
603589a23caSBarry Smith     ierr = PetscGetProgramName(programname,sizeof(programname));CHKERRQ(ierr);
60411525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
60511525c0dSBarry Smith     introlen   = 4096 + applinelen;
60630a8c9c0SSurtai Han     applinelen += 1024;
60711525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
60811525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
60911525c0dSBarry Smith 
61011525c0dSBarry Smith     if (rootlocal) {
61111525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
61211525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
61311525c0dSBarry Smith     }
61476a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
61511525c0dSBarry Smith     if (rootlocal && help) {
616928bb9adSStefano 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);
61711525c0dSBarry Smith     } else if (help) {
618928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
61911525c0dSBarry Smith     } else {
620928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
62111525c0dSBarry Smith     }
622b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
62330befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
62411525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
625a17b96a8SKyle 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"
626df62c222SBarry 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"
627928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
62843da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
62911525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
63011525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
631ffbd1cfbSBarry Smith     if (selectport) {
632aa573868SBarry Smith       PetscBool silent;
6337d812c46SBarry Smith 
6347d812c46SBarry Smith       ierr = SAWs_Initialize();
6357d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6367d812c46SBarry Smith       while (ierr) {
6377d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6387d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6397d812c46SBarry Smith         ierr = SAWs_Initialize();
6407d812c46SBarry Smith       }
6417d812c46SBarry Smith 
642aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
643aa573868SBarry Smith       if (!silent) {
644ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
645ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
646ffbd1cfbSBarry Smith       }
6477d812c46SBarry Smith     } else {
6487d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
649aa573868SBarry Smith     }
6500af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6510af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6520af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6530af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6540af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
65511525c0dSBarry Smith   }
65611525c0dSBarry Smith   PetscFunctionReturn(0);
65711525c0dSBarry Smith }
65811525c0dSBarry Smith #endif
65911525c0dSBarry Smith 
6604dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6614dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6624dfee713SSatish Balay {
6634dfee713SSatish Balay   PetscFunctionBegin;
6644dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6654dfee713SSatish Balay     /* see MPI.py for details on this bug */
6664dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6674dfee713SSatish Balay #endif
6684dfee713SSatish Balay   PetscFunctionReturn(0);
6694dfee713SSatish Balay }
6704dfee713SSatish Balay 
671a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
672a56f64adSBarry Smith #include <adios.h>
67322580e64SBarry Smith #include <adios_read.h>
6747b56e58cSSatish Balay int64_t Petsc_adios_group;
675a56f64adSBarry Smith #endif
676cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
677cd1b99a6SBarry Smith #include <omp.h>
678f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
679cd1b99a6SBarry Smith #endif
680a56f64adSBarry Smith 
68127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H)
682bc8a24ecSBarry Smith #include <dlfcn.h>
683bc8a24ecSBarry Smith #endif
684bc8a24ecSBarry Smith 
68585649d77SJunchao Zhang /*
68685649d77SJunchao Zhang   PetscInitialize_Common  - shared code between C and Fortran initialization
68785649d77SJunchao Zhang 
68885649d77SJunchao Zhang   prog:     program name
68902101c96SSatish Balay   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
69085649d77SJunchao Zhang   help:     program help message
69102101c96SSatish Balay   ftn:      is it called from Fortran initilization (petscinitializef_)?
69285649d77SJunchao Zhang   readarguments,len: used when fortran is true
69385649d77SJunchao Zhang */
69402101c96SSatish Balay PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char* prog,const char* file,const char *help,PetscBool ftn,PetscBool readarguments,PetscInt len)
69585649d77SJunchao Zhang {
69685649d77SJunchao Zhang   PetscErrorCode ierr;
69785649d77SJunchao Zhang   PetscMPIInt    size;
69885649d77SJunchao Zhang   PetscBool      flg = PETSC_TRUE;
69985649d77SJunchao Zhang   char           hostname[256];
70085649d77SJunchao Zhang 
70127104ee2SJacob Faibussowitsch   PetscFunctionBegin;
70227104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
70385649d77SJunchao Zhang   /*
70485649d77SJunchao Zhang       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
70585649d77SJunchao Zhang       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
70685649d77SJunchao Zhang         MPICH v3.1 (Released February 2014)
70785649d77SJunchao Zhang         IBM MPI v2.1 (December 2014)
70885649d77SJunchao Zhang         Intel MPI Library v5.0 (2014)
70985649d77SJunchao Zhang         Cray MPT v7.0.0 (June 2014)
71085649d77SJunchao Zhang       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
71185649d77SJunchao Zhang       listed above and since that time are compatible.
71285649d77SJunchao Zhang 
71385649d77SJunchao Zhang       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
71485649d77SJunchao Zhang       at compile time or runtime. Thus we will need to systematically track the allowed versions
71585649d77SJunchao Zhang       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
71685649d77SJunchao Zhang       to perform the checking.
71785649d77SJunchao Zhang 
71885649d77SJunchao Zhang       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
71985649d77SJunchao Zhang 
72085649d77SJunchao Zhang       Questions:
72185649d77SJunchao Zhang 
72285649d77SJunchao Zhang         Should the checks for ABI incompatibility be only on the major version number below?
72385649d77SJunchao Zhang         Presumably the output to stderr will be removed before a release.
72485649d77SJunchao Zhang   */
72585649d77SJunchao Zhang 
72685649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
72785649d77SJunchao Zhang   {
72885649d77SJunchao Zhang     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
72985649d77SJunchao Zhang     PetscMPIInt mpilibraryversionlength;
73027104ee2SJacob Faibussowitsch     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);
73127104ee2SJacob Faibussowitsch     if (ierr) PetscFunctionReturn(ierr);
73285649d77SJunchao Zhang     /* check for MPICH versions before MPI ABI initiative */
73385649d77SJunchao Zhang #if defined(MPICH_VERSION)
73485649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000
73585649d77SJunchao Zhang     {
73685649d77SJunchao Zhang       char      *ver,*lf;
73785649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
73827104ee2SJacob Faibussowitsch       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);
73927104ee2SJacob Faibussowitsch       if (ierr) PetscFunctionReturn(ierr);
74027104ee2SJacob Faibussowitsch       else if (ver) {
74127104ee2SJacob Faibussowitsch         ierr = PetscStrchr(ver,'\n',&lf);
74227104ee2SJacob Faibussowitsch         if (ierr) PetscFunctionReturn(ierr);
74327104ee2SJacob Faibussowitsch         else if (lf) {
74485649d77SJunchao Zhang           *lf = 0;
74527104ee2SJacob Faibussowitsch           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) PetscFunctionReturn(ierr);
74685649d77SJunchao Zhang         }
74785649d77SJunchao Zhang       }
74885649d77SJunchao Zhang       if (!flg) {
74985649d77SJunchao Zhang         PetscInfo1(NULL,"PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %.\n",mpilibraryversion,MPICH_VESION);
75085649d77SJunchao Zhang         flg = PETSC_TRUE;
75185649d77SJunchao Zhang       }
75285649d77SJunchao Zhang     }
75385649d77SJunchao Zhang #endif
75485649d77SJunchao 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?) */
75585649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION)
75685649d77SJunchao Zhang     {
75785649d77SJunchao Zhang       char *ver,bs[MPI_MAX_LIBRARY_VERSION_STRING],*bsf;
75885649d77SJunchao Zhang       PetscBool flg = PETSC_FALSE;
75985649d77SJunchao Zhang #define PSTRSZ 2
76085649d77SJunchao Zhang       char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI","FUJITSU MPI"};
76185649d77SJunchao Zhang       char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v","Library "};
76285649d77SJunchao Zhang       int i;
76385649d77SJunchao Zhang       for (i=0; i<PSTRSZ; i++) {
76427104ee2SJacob Faibussowitsch         ierr = PetscStrstr(mpilibraryversion,ompistr1[i],&ver);
76527104ee2SJacob Faibussowitsch         if (ierr) PetscFunctionReturn(ierr);
76627104ee2SJacob Faibussowitsch         else if (ver) {
76785649d77SJunchao Zhang           PetscSNPrintf(bs,MPI_MAX_LIBRARY_VERSION_STRING,"%s%d.%d",ompistr2[i],OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
76827104ee2SJacob Faibussowitsch           ierr = PetscStrstr(ver,bs,&bsf);
76927104ee2SJacob Faibussowitsch           if (ierr) PetscFunctionReturn(ierr);
77027104ee2SJacob Faibussowitsch           else if (bsf) flg = PETSC_TRUE;
77185649d77SJunchao Zhang           break;
77285649d77SJunchao Zhang         }
77385649d77SJunchao Zhang       }
77485649d77SJunchao Zhang       if (!flg) {
77585649d77SJunchao Zhang         PetscInfo3(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);
77685649d77SJunchao Zhang         flg = PETSC_TRUE;
77785649d77SJunchao Zhang       }
77885649d77SJunchao Zhang     }
77985649d77SJunchao Zhang #endif
78085649d77SJunchao Zhang   }
78185649d77SJunchao Zhang #endif
78285649d77SJunchao Zhang 
78385649d77SJunchao Zhang #if defined(PETSC_HAVE_DLSYM)
78485649d77SJunchao 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 */
78527104ee2SJacob Faibussowitsch   if (PetscUnlikely(dlsym(RTLD_DEFAULT,"ompi_mpi_init") && dlsym(RTLD_DEFAULT,"MPID_Abort"))) {
78685649d77SJunchao Zhang     fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
78727104ee2SJacob Faibussowitsch     ierr = PetscStackView(stderr);CHKERRQ(ierr);
78827104ee2SJacob Faibussowitsch     PetscFunctionReturn(PETSC_ERR_MPI_LIB_INCOMP);
78985649d77SJunchao Zhang   }
79085649d77SJunchao Zhang #endif
79185649d77SJunchao Zhang 
79285649d77SJunchao Zhang   /* these must be initialized in a routine, not as a constant declaration*/
79385649d77SJunchao Zhang   PETSC_STDOUT = stdout;
79485649d77SJunchao Zhang   PETSC_STDERR = stderr;
79585649d77SJunchao Zhang 
79685649d77SJunchao Zhang   /*CHKERRQ can be used from now */
79785649d77SJunchao Zhang   PetscErrorHandlingInitialized = PETSC_TRUE;
79885649d77SJunchao Zhang 
79985649d77SJunchao Zhang   /* on Windows - set printf to default to printing 2 digit exponents */
80085649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
80185649d77SJunchao Zhang   _set_output_format(_TWO_DIGIT_EXPONENT);
80285649d77SJunchao Zhang #endif
80385649d77SJunchao Zhang 
80485649d77SJunchao Zhang   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
80585649d77SJunchao Zhang 
80685649d77SJunchao Zhang   PetscFinalizeCalled = PETSC_FALSE;
80785649d77SJunchao Zhang 
80885649d77SJunchao Zhang   ierr = PetscSetProgramName(prog);CHKERRQ(ierr);
80985649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
81085649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
81185649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
81285649d77SJunchao Zhang   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
81385649d77SJunchao Zhang 
81485649d77SJunchao Zhang   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
81585649d77SJunchao Zhang   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRMPI(ierr);
81685649d77SJunchao Zhang 
81785649d77SJunchao Zhang   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
81885649d77SJunchao Zhang     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRMPI(ierr);
81985649d77SJunchao Zhang     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRMPI(ierr);
82085649d77SJunchao Zhang   }
82185649d77SJunchao Zhang 
82285649d77SJunchao Zhang   /* Done after init due to a bug in MPICH-GM? */
82385649d77SJunchao Zhang   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
82485649d77SJunchao Zhang 
82585649d77SJunchao Zhang   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRMPI(ierr);
82685649d77SJunchao Zhang   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRMPI(ierr);
82785649d77SJunchao Zhang 
82885649d77SJunchao Zhang   MPIU_BOOL = MPI_INT;
82985649d77SJunchao Zhang   MPIU_ENUM = MPI_INT;
83085649d77SJunchao Zhang   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
83185649d77SJunchao Zhang   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
83285649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
83385649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG)
83485649d77SJunchao Zhang   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
83585649d77SJunchao Zhang #endif
83627104ee2SJacob Faibussowitsch   else {
83727104ee2SJacob Faibussowitsch     (*PetscErrorPrintf)("PetscInitialize_Common: Could not find MPI type for size_t\n");
83827104ee2SJacob Faibussowitsch     PetscFunctionReturn(PETSC_ERR_SUP_SYS);
83927104ee2SJacob Faibussowitsch   }
84085649d77SJunchao Zhang 
84185649d77SJunchao Zhang   /*
84285649d77SJunchao Zhang      Initialized the global complex variable; this is because with
84385649d77SJunchao Zhang      shared libraries the constructors for global variables
84485649d77SJunchao Zhang      are not called; at least on IRIX.
84585649d77SJunchao Zhang   */
84685649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
84785649d77SJunchao Zhang   {
84885649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
84985649d77SJunchao Zhang     PetscComplex ic(0.0,1.0);
85085649d77SJunchao Zhang     PETSC_i = ic;
85185649d77SJunchao Zhang #else
85285649d77SJunchao Zhang     PETSC_i = _Complex_I;
85385649d77SJunchao Zhang #endif
85485649d77SJunchao Zhang   }
85585649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */
85685649d77SJunchao Zhang 
85785649d77SJunchao Zhang   /*
85885649d77SJunchao Zhang      Create the PETSc MPI reduction operator that sums of the first
85985649d77SJunchao Zhang      half of the entries and maxes the second half.
86085649d77SJunchao Zhang   */
86185649d77SJunchao Zhang   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRMPI(ierr);
86285649d77SJunchao Zhang 
86385649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128)
86485649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRMPI(ierr);
86585649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRMPI(ierr);
86685649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX)
86785649d77SJunchao Zhang   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRMPI(ierr);
86885649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRMPI(ierr);
86985649d77SJunchao Zhang #endif
87085649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRMPI(ierr);
87185649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRMPI(ierr);
87285649d77SJunchao Zhang #elif defined(PETSC_USE_REAL___FP16)
87385649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRMPI(ierr);
87485649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRMPI(ierr);
87585649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRMPI(ierr);
87685649d77SJunchao Zhang   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRMPI(ierr);
87785649d77SJunchao Zhang #endif
87885649d77SJunchao Zhang 
87985649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
88085649d77SJunchao Zhang   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRMPI(ierr);
88185649d77SJunchao Zhang #endif
88285649d77SJunchao Zhang 
88385649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRMPI(ierr);
88485649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRMPI(ierr);
88585649d77SJunchao Zhang 
88685649d77SJunchao Zhang   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
88785649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI)
88885649d77SJunchao Zhang   {
88985649d77SJunchao Zhang     struct PetscRealInt { PetscReal v; PetscInt i; };
89085649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
89185649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscRealInt,v),offsetof(struct PetscRealInt,i)};
89285649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_REAL,MPIU_INT}, tmpStruct;
89385649d77SJunchao Zhang 
89485649d77SJunchao Zhang     ierr = MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct);CHKERRMPI(ierr);
89585649d77SJunchao Zhang     ierr = MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscRealInt),&MPIU_REAL_INT);CHKERRMPI(ierr);
89685649d77SJunchao Zhang     ierr = MPI_Type_free(&tmpStruct);CHKERRMPI(ierr);
89785649d77SJunchao Zhang     ierr = MPI_Type_commit(&MPIU_REAL_INT);CHKERRMPI(ierr);
89885649d77SJunchao Zhang   }
89985649d77SJunchao Zhang   {
90085649d77SJunchao Zhang     struct PetscScalarInt { PetscScalar v; PetscInt i; };
90185649d77SJunchao Zhang     PetscMPIInt  blockSizes[2] = {1,1};
90285649d77SJunchao Zhang     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscScalarInt,v),offsetof(struct PetscScalarInt,i)};
90385649d77SJunchao Zhang     MPI_Datatype blockTypes[2] = {MPIU_SCALAR,MPIU_INT}, tmpStruct;
90485649d77SJunchao Zhang 
90585649d77SJunchao Zhang     ierr = MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct);CHKERRMPI(ierr);
90685649d77SJunchao Zhang     ierr = MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscScalarInt),&MPIU_SCALAR_INT);CHKERRMPI(ierr);
90785649d77SJunchao Zhang     ierr = MPI_Type_free(&tmpStruct);CHKERRMPI(ierr);
90885649d77SJunchao Zhang     ierr = MPI_Type_commit(&MPIU_SCALAR_INT);CHKERRMPI(ierr);
90985649d77SJunchao Zhang   }
91085649d77SJunchao Zhang #endif
91185649d77SJunchao Zhang 
91285649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES)
91385649d77SJunchao Zhang   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRMPI(ierr);
91485649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRMPI(ierr);
91585649d77SJunchao Zhang #endif
91685649d77SJunchao Zhang   ierr = MPI_Type_contiguous(4,MPI_INT,&MPI_4INT);CHKERRMPI(ierr);
91785649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPI_4INT);CHKERRMPI(ierr);
91885649d77SJunchao Zhang   ierr = MPI_Type_contiguous(4,MPIU_INT,&MPIU_4INT);CHKERRMPI(ierr);
91985649d77SJunchao Zhang   ierr = MPI_Type_commit(&MPIU_4INT);CHKERRMPI(ierr);
92085649d77SJunchao Zhang 
92185649d77SJunchao Zhang   /*
92285649d77SJunchao Zhang      Attributes to be set on PETSc communicators
92385649d77SJunchao Zhang   */
92485649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);CHKERRMPI(ierr);
92585649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);CHKERRMPI(ierr);
92685649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);CHKERRMPI(ierr);
92785649d77SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);CHKERRMPI(ierr);
92885649d77SJunchao Zhang 
92985649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN)
93002101c96SSatish Balay   if (ftn) {ierr = PetscInitFortran_Private(readarguments,file,len);CHKERRQ(ierr);}
93185649d77SJunchao Zhang   else
93285649d77SJunchao Zhang #endif
93385649d77SJunchao Zhang   {ierr = PetscOptionsInsert(NULL,&PetscGlobalArgc,&PetscGlobalArgs,file);CHKERRQ(ierr);}
93485649d77SJunchao Zhang 
93585649d77SJunchao Zhang   /* call a second time so it can look in the options database */
93685649d77SJunchao Zhang   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
93785649d77SJunchao Zhang 
93885649d77SJunchao Zhang   /*
93985649d77SJunchao Zhang      Check system options and print help
94085649d77SJunchao Zhang   */
94185649d77SJunchao Zhang   ierr = PetscOptionsCheckInitial_Private(help);CHKERRQ(ierr);
94285649d77SJunchao Zhang 
94385649d77SJunchao Zhang   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
94485649d77SJunchao Zhang 
94585649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS)
94602101c96SSatish Balay   ierr = PetscInitializeSAWs(ftn ? NULL : help);CHKERRQ(ierr);
94727104ee2SJacob Faibussowitsch   flg = PETSC_FALSE;
94827104ee2SJacob Faibussowitsch   ierr = PetscOptionsHasName(NULL,NULL,"-stack_view",&flg);CHKERRQ(ierr);
94927104ee2SJacob Faibussowitsch   if (flg) PetscStackViewSAWs();
95085649d77SJunchao Zhang #endif
95185649d77SJunchao Zhang 
95285649d77SJunchao Zhang   /*
95385649d77SJunchao Zhang      Load the dynamic libraries (on machines that support them), this registers all
95485649d77SJunchao Zhang      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
95585649d77SJunchao Zhang   */
95685649d77SJunchao Zhang   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
95785649d77SJunchao Zhang 
95885649d77SJunchao Zhang   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
95985649d77SJunchao Zhang   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
96085649d77SJunchao Zhang   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
96185649d77SJunchao Zhang   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
96285649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP)
96385649d77SJunchao Zhang   {
96485649d77SJunchao Zhang     PetscBool omp_view_flag;
96585649d77SJunchao Zhang     char      *threads = getenv("OMP_NUM_THREADS");
96685649d77SJunchao Zhang 
96785649d77SJunchao Zhang     if (threads) {
9682f840973SStefano Zampini       ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
96985649d77SJunchao Zhang       (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
97085649d77SJunchao Zhang     } else {
9712f840973SStefano Zampini       PetscNumOMPThreads = (PetscInt) omp_get_max_threads();
9722f840973SStefano Zampini       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (as given by omp_get_max_threads())\n",PetscNumOMPThreads);CHKERRQ(ierr);
97385649d77SJunchao Zhang     }
97485649d77SJunchao Zhang     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
97585649d77SJunchao 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);
97685649d77SJunchao Zhang     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
97785649d77SJunchao Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
97885649d77SJunchao Zhang     if (flg) {
97985649d77SJunchao Zhang       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
98085649d77SJunchao Zhang       omp_set_num_threads((int)PetscNumOMPThreads);
98185649d77SJunchao Zhang     }
98285649d77SJunchao Zhang     if (omp_view_flag) {
98385649d77SJunchao Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
98485649d77SJunchao Zhang     }
98585649d77SJunchao Zhang   }
98685649d77SJunchao Zhang #endif
98785649d77SJunchao Zhang 
98885649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
98985649d77SJunchao Zhang   /*
99085649d77SJunchao Zhang       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
99185649d77SJunchao Zhang 
99285649d77SJunchao Zhang       Currently not used because it is not supported by MPICH.
99385649d77SJunchao Zhang   */
99485649d77SJunchao Zhang   if (!PetscBinaryBigEndian()) {
99585649d77SJunchao Zhang     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRMPI(ierr);
99685649d77SJunchao Zhang   }
99785649d77SJunchao Zhang #endif
99885649d77SJunchao Zhang 
99985649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS)
100085649d77SJunchao Zhang   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
100185649d77SJunchao Zhang #endif
100285649d77SJunchao Zhang 
100385649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC)
100485649d77SJunchao Zhang   {
100585649d77SJunchao Zhang     PetscViewer viewer;
100685649d77SJunchao Zhang     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
100785649d77SJunchao Zhang     if (flg) {
100885649d77SJunchao Zhang       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
100985649d77SJunchao Zhang       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
101085649d77SJunchao Zhang     }
101185649d77SJunchao Zhang   }
101285649d77SJunchao Zhang #endif
101385649d77SJunchao Zhang 
101485649d77SJunchao Zhang   flg  = PETSC_TRUE;
101585649d77SJunchao Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
101685649d77SJunchao Zhang   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);}
101785649d77SJunchao Zhang 
101885649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS)
101985649d77SJunchao Zhang   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
102085649d77SJunchao Zhang   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
102185649d77SJunchao Zhang   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
102285649d77SJunchao Zhang   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
102385649d77SJunchao Zhang #endif
102485649d77SJunchao Zhang 
102585649d77SJunchao Zhang #if defined(__VALGRIND_H)
102685649d77SJunchao Zhang   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND? PETSC_TRUE: PETSC_FALSE;
102785649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
102885649d77SJunchao Zhang   if (PETSC_RUNNING_ON_VALGRIND) {
102985649d77SJunchao 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);
103085649d77SJunchao Zhang     }
103185649d77SJunchao Zhang #endif
103285649d77SJunchao Zhang #endif
103385649d77SJunchao Zhang 
103485649d77SJunchao Zhang #if (defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)) && defined(PETSC_EXPERIMENTAL)
103585649d77SJunchao Zhang   ierr = PetscDeviceInitializeDefaultDevices_Internal();CHKERRQ(ierr);
103685649d77SJunchao Zhang   ierr = PetscDeviceContextInitializeRootContext_Internal(PETSC_COMM_WORLD,NULL);CHKERRQ(ierr);
103785649d77SJunchao Zhang #endif
103885649d77SJunchao Zhang 
103985649d77SJunchao Zhang   /*
104085649d77SJunchao Zhang       Set flag that we are completely initialized
104185649d77SJunchao Zhang   */
104285649d77SJunchao Zhang   PetscInitializeCalled = PETSC_TRUE;
104385649d77SJunchao Zhang 
104485649d77SJunchao Zhang   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
104585649d77SJunchao Zhang   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
104627104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
104785649d77SJunchao Zhang }
104885649d77SJunchao Zhang 
1049e5c89e4eSSatish Balay /*@C
1050e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
1051e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
1052e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
1053e5c89e4eSSatish Balay    your program -- usually the very first line!
1054e5c89e4eSSatish Balay 
1055e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
1056e5c89e4eSSatish Balay 
1057e5c89e4eSSatish Balay    Input Parameters:
1058e5c89e4eSSatish Balay +  argc - count of number of command line arguments
1059e5c89e4eSSatish Balay .  args - the command line arguments
1060be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1061be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
1062be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
1063c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
10640298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
1065e5c89e4eSSatish Balay 
106605827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
106705827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
106805827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
106905827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
107005827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
1071e5c89e4eSSatish Balay 
1072e5c89e4eSSatish Balay    Options Database Keys:
10737ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
10747ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1075e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
10768a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
10778a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
10788a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
10798a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
10808a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1081bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1082e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1083e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
1084e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
108579dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
108679dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
108779dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
1088aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
108992f119d6SBarry 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
109092f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
109192f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1092608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
109392f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
1094e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
1095e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
1096e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
1097e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
1098e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
10990841954dSBarry Smith -  -memory_view - Print memory usage at end of run
1100e5c89e4eSSatish Balay 
1101c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
1102c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1103c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
1104c5b5d8d5SVaclav Hapla -  -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet()
1105c5b5d8d5SVaclav Hapla 
1106c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
1107c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
1108c5b5d8d5SVaclav Hapla    They can be used also in option files.
1109c5b5d8d5SVaclav Hapla 
1110c5b5d8d5SVaclav Hapla    See PetscOptionsMonitorSet() to do monitoring programmatically.
1111c5b5d8d5SVaclav Hapla 
1112e5c89e4eSSatish Balay    Options Database Keys for Profiling:
1113a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
1114fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
1115217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1116217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
1117495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1118e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
11199a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
112079dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
11219a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1122495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
1123f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1124495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
1125495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
1126c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
112787c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
1128c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
1129495fc317SBarry Smith 
1130609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
1131e5c89e4eSSatish Balay 
1132ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
1133ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1134ffbd1cfbSBarry 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
1135ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
1136ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
1137ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1138ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
1139ffbd1cfbSBarry Smith 
1140e5c89e4eSSatish Balay    Environmental Variables:
1141e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
1142e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
1143e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
11444a971ea4SToby Isaac .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
11454a971ea4SToby Isaac .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1146e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
1147e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
1148e5c89e4eSSatish Balay 
1149e5c89e4eSSatish Balay    Level: beginner
1150e5c89e4eSSatish Balay 
1151e5c89e4eSSatish Balay    Notes:
1152e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
1153e5c89e4eSSatish Balay    it before PetscInitialize().
1154e5c89e4eSSatish Balay 
1155e5c89e4eSSatish Balay    Fortran Version:
1156e5c89e4eSSatish Balay    In Fortran this routine has the format
1157e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
1158e5c89e4eSSatish Balay 
1159e5c89e4eSSatish Balay +  ierr - error return code
1160c5b5d8d5SVaclav Hapla -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
1161c5b5d8d5SVaclav Hapla           Use PETSC_NULL_CHARACTER to not check for code specific file.
1162c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
1163e5c89e4eSSatish Balay 
1164e5c89e4eSSatish Balay    Important Fortran Note:
11650eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
11660298fd71SBarry Smith    null character string; you CANNOT just use NULL as
1167a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
1168e5c89e4eSSatish Balay 
116901cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
117001cb0274SBarry Smith    calling PetscInitialize().
1171e5c89e4eSSatish Balay 
117201cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
1173e5c89e4eSSatish Balay 
1174e5c89e4eSSatish Balay @*/
11757087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
1176e5c89e4eSSatish Balay {
1177e5c89e4eSSatish Balay   PetscErrorCode ierr;
117885649d77SJunchao Zhang   PetscMPIInt    flag;
117985649d77SJunchao Zhang   const char     *prog = "Unknown Name";
1180e5c89e4eSSatish Balay 
118127104ee2SJacob Faibussowitsch   PetscFunctionBegin;
118227104ee2SJacob Faibussowitsch   if (PetscInitializeCalled) PetscFunctionReturn(0);
1183ffc4695bSBarry Smith   ierr = MPI_Initialized(&flag);CHKERRMPI(ierr);
1184e5c89e4eSSatish Balay   if (!flag) {
1185e32f2f54SBarry Smith     if (PETSC_COMM_WORLD != MPI_COMM_NULL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You cannot set PETSC_COMM_WORLD if you have not initialized MPI first");
11864dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
11875e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
11885e765c61SJed Brown     {
11895e765c61SJed Brown       PetscMPIInt provided;
1190ffc4695bSBarry Smith       ierr = MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided);CHKERRMPI(ierr);
11915e765c61SJed Brown     }
11925e765c61SJed Brown #else
1193ffc4695bSBarry Smith     ierr = MPI_Init(argc,args);CHKERRMPI(ierr);
11945e765c61SJed Brown #endif
1195e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
1196e5c89e4eSSatish Balay   }
11974dfee713SSatish Balay 
119885649d77SJunchao Zhang   if (argc && *argc) prog = **args;
1199e5c89e4eSSatish Balay   if (argc && args) {
1200e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
1201e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
1202e5c89e4eSSatish Balay   }
120385649d77SJunchao Zhang   ierr = PetscInitialize_Common(prog,file,help,PETSC_FALSE/*C*/,PETSC_FALSE,0);CHKERRQ(ierr);
120427104ee2SJacob Faibussowitsch   PetscFunctionReturn(0);
1205e5c89e4eSSatish Balay }
1206e5c89e4eSSatish Balay 
12074097062eSBarry Smith #if defined(PETSC_USE_LOG)
120895c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1209ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1210ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
121195c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
12124097062eSBarry Smith #endif
1213e5c89e4eSSatish Balay 
1214008a6e76SBarry Smith /*
1215008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1216008a6e76SBarry Smith */
1217008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1218008a6e76SBarry Smith {
1219008a6e76SBarry Smith   PetscErrorCode ierr;
1220008a6e76SBarry Smith 
1221008a6e76SBarry Smith   PetscFunctionBegin;
1222008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1223ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRMPI(ierr);
1224008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1225ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRMPI(ierr);
1226008a6e76SBarry Smith #endif
1227ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRMPI(ierr);
1228ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRMPI(ierr);
1229008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1230ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRMPI(ierr);
1231ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRMPI(ierr);
1232ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRMPI(ierr);
1233008a6e76SBarry Smith #endif
1234008a6e76SBarry Smith 
1235de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1236ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRMPI(ierr);
1237008a6e76SBarry Smith #endif
1238008a6e76SBarry Smith 
1239ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRMPI(ierr);
1240092991acSStefano Zampini   ierr = MPI_Type_free(&MPIU_REAL_INT);CHKERRMPI(ierr);
1241092991acSStefano Zampini   ierr = MPI_Type_free(&MPIU_SCALAR_INT);CHKERRMPI(ierr);
124240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1243ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRMPI(ierr);
1244008a6e76SBarry Smith #endif
1245b5a892a1SMatthew G. Knepley   ierr = MPI_Type_free(&MPI_4INT);CHKERRMPI(ierr);
1246b5a892a1SMatthew G. Knepley   ierr = MPI_Type_free(&MPIU_4INT);CHKERRMPI(ierr);
1247ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRMPI(ierr);
1248008a6e76SBarry Smith   PetscFunctionReturn(0);
1249008a6e76SBarry Smith }
1250008a6e76SBarry Smith 
1251e5c89e4eSSatish Balay /*@C
1252e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1253e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1254e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1255e5c89e4eSSatish Balay 
1256e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1257e5c89e4eSSatish Balay 
1258e5c89e4eSSatish Balay    Options Database Keys:
125926a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1260e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
12617eb1d149SBarry 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
1262e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
126392f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1264e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
126592f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1266e5c89e4eSSatish Balay 
1267e5c89e4eSSatish Balay    Level: beginner
1268e5c89e4eSSatish Balay 
1269e5c89e4eSSatish Balay    Note:
1270e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1271e5c89e4eSSatish Balay 
127288c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1273e5c89e4eSSatish Balay @*/
12747087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1275e5c89e4eSSatish Balay {
1276e5c89e4eSSatish Balay   PetscErrorCode ierr;
12774bb5149bSJed Brown   PetscMPIInt    rank;
1278a8d2bbe5SBarry Smith   PetscInt       nopt;
12792bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1280dff31646SBarry Smith   PetscBool      flg;
128110463e74SBarry Smith #if defined(PETSC_USE_LOG)
128210463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
128310463e74SBarry Smith #endif
1284e5c89e4eSSatish Balay 
12853cbbc5ffSBarry Smith   PetscFunctionBegin;
128627104ee2SJacob Faibussowitsch   if (PetscUnlikely(!PetscInitializeCalled)) {
128727104ee2SJacob Faibussowitsch     fprintf(stderr,"PetscInitialize() must be called before PetscFinalize()\n");
128827104ee2SJacob Faibussowitsch     ierr = PetscStackView(stderr);CHKERRQ(ierr);
128927104ee2SJacob Faibussowitsch     return PETSC_ERR_ARG_WRONGSTATE;
129027104ee2SJacob Faibussowitsch   }
12910298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1292b022a5c1SBarry Smith 
1293ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
1294a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
129522580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1296a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1297a56f64adSBarry Smith #endif
1298c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1299dff31646SBarry Smith   if (flg) {
13001f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
13011f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
13021f817a21SBarry Smith 
1303589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL);CHKERRQ(ierr);
13041f817a21SBarry Smith     if (filename[0]) {
13051f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
13061f817a21SBarry Smith     }
1307dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1308dff31646SBarry Smith     cits[0] = 0;
1309dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
13101f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
13111f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
13121f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
13131f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
13141f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1315dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1316dff31646SBarry Smith   }
1317dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1318dff31646SBarry Smith 
1319c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
132004102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
132104102261SBarry Smith   {
132204102261SBarry Smith     PetscInt nmax = 2;
132304102261SBarry Smith     char     **buffs;
132404102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1325c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
132604102261SBarry Smith     if (flg1) {
132727104ee2SJacob Faibussowitsch       if (PetscUnlikely(!nmax)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
132804102261SBarry Smith       if (nmax == 1) {
132904102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
133004102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
133104102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
133204102261SBarry Smith       }
133304102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
133404102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
133504102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
133604102261SBarry Smith     }
133704102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
133804102261SBarry Smith   }
1339763ec7b1SBarry Smith   {
1340763ec7b1SBarry Smith     PetscInt nmax = 2;
1341763ec7b1SBarry Smith     char     **buffs;
1342763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1343763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1344763ec7b1SBarry Smith     if (flg1) {
1345763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1346763ec7b1SBarry Smith       if (nmax == 1) {
1347763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1348763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1349763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1350763ec7b1SBarry Smith       }
1351763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1352763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1353763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1354763ec7b1SBarry Smith     }
1355763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1356763ec7b1SBarry Smith   }
135704102261SBarry Smith #endif
135804102261SBarry Smith 
13592d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
13602d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
13612d53ad75SBarry Smith #endif
13622d53ad75SBarry Smith 
1363e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1364dff31646SBarry Smith   flg = PETSC_FALSE;
1365c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1366d5649816SBarry Smith   if (flg) {
1367e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1368d5649816SBarry Smith   }
1369d5649816SBarry Smith #endif
1370d5649816SBarry Smith 
1371681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1372681455b2SBarry Smith   flg1 = PETSC_FALSE;
1373c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1374681455b2SBarry Smith   if (flg1) {
1375681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1376681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1377681455b2SBarry Smith   }
1378681455b2SBarry Smith #endif
1379681455b2SBarry Smith 
138067584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1381c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1382e5c89e4eSSatish Balay   if (!flg2) {
138390d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1384c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1385e5c89e4eSSatish Balay   }
1386e5c89e4eSSatish Balay   if (flg2) {
13870841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1388e5c89e4eSSatish Balay   }
138967584ceeSBarry Smith #endif
1390e5c89e4eSSatish Balay 
1391e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
139290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1393c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1394e5c89e4eSSatish Balay   if (flg1) {
1395e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1396ffc4695bSBarry Smith     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRMPI(ierr);
1397e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1398e5c89e4eSSatish Balay   }
1399e5c89e4eSSatish Balay #endif
1400e5c89e4eSSatish Balay 
1401e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1402e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1403e5c89e4eSSatish Balay   mname[0] = 0;
1404589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1405e5c89e4eSSatish Balay   if (flg1) {
1406e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1407e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1408e5c89e4eSSatish Balay   }
1409e5c89e4eSSatish Balay #endif
1410356e5837SBarry Smith #endif
1411a297a907SKarl Rupp 
1412dd710f27SBarry Smith   /*
1413dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1414dd710f27SBarry Smith   */
1415dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1416dd710f27SBarry Smith 
1417356e5837SBarry Smith #if defined(PETSC_USE_LOG)
141887c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1419f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
142087c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
142187c3beb6SLisandro Dalcin 
1422356e5837SBarry Smith   mname[0] = 0;
1423589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1424e5c89e4eSSatish Balay   if (flg1) {
142591eabc43SBarry Smith     PetscViewer viewer;
142620a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
142791eabc43SBarry Smith     if (mname[0]) {
142891eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
142991eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
14306bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
143133f85c2fSBarry Smith     } else {
143233f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
14339a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
143433f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
14359a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
143633f85c2fSBarry Smith     }
1437e5c89e4eSSatish Balay   }
1438a297a907SKarl Rupp 
1439dd710f27SBarry Smith   /*
1440dd710f27SBarry Smith      Free any objects created by the last block of code.
1441dd710f27SBarry Smith   */
1442dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1443dd710f27SBarry Smith 
1444dd710f27SBarry Smith   mname[0] = 0;
1445589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1446589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2);CHKERRQ(ierr);
14477ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1448e5c89e4eSSatish Balay #endif
144910463e74SBarry Smith 
145090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1451c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1452e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
145390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1454c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1455e5c89e4eSSatish Balay   if (flg1) {
1456e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1457e5c89e4eSSatish Balay   }
145890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
145990d69ab7SBarry Smith   flg2 = PETSC_FALSE;
14608bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1461c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1462c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1463c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1464e4c476e2SSatish Balay 
1465e5c89e4eSSatish Balay   if (flg2) {
1466be56827dSJed Brown     PetscViewer viewer;
146702ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
146802ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1469c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1470be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1471e5c89e4eSSatish Balay   }
1472e5c89e4eSSatish Balay 
1473e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1474c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1475c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1476e5c89e4eSSatish Balay 
147733fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1478c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
14799245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1480e5c89e4eSSatish Balay   if (flg3) {
14813de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1482be56827dSJed Brown       PetscViewer viewer;
148302ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
148402ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1485c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1486be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1487e5c89e4eSSatish Balay     }
14883de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
14893de2bfdfSBarry Smith     if (nopt) {
14903de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
14913de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
14923de2bfdfSBarry Smith       if (nopt == 1) {
1493e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1494e5c89e4eSSatish Balay       } else {
14957582186dSLisandro Dalcin         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1496e5c89e4eSSatish Balay       }
14973de2bfdfSBarry Smith     } else if (flg3 && flg1) {
14983de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1499df12ba86SBarry Smith     }
1500c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1501e5c89e4eSSatish Balay   }
1502e5c89e4eSSatish Balay 
1503e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1504d45a07a7SBarry Smith   if (!PetscGlobalRank) {
150587f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
150616ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1507d45a07a7SBarry Smith   }
1508ec957eceSBarry Smith #endif
1509ec957eceSBarry Smith 
15104097062eSBarry Smith #if defined(PETSC_USE_LOG)
151110463e74SBarry Smith   /*
1512dbc8283eSBarry Smith        List all objects the user may have forgot to free
15132eff7a51SBarry Smith   */
151405df10baSBarry Smith   if (PetscObjectsLog) {
1515c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1516a64a8e02SBarry Smith     if (flg1) {
1517a64a8e02SBarry Smith       MPI_Comm local_comm;
15187eb1d149SBarry Smith       char     string[64];
1519a64a8e02SBarry Smith 
1520589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL);CHKERRQ(ierr);
1521ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(ierr);
1522a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
15237eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1524a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1525ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(ierr);
15260a1571b3SBarry Smith     }
152705df10baSBarry Smith   }
15284097062eSBarry Smith #endif
15294097062eSBarry Smith 
15304097062eSBarry Smith #if defined(PETSC_USE_LOG)
1531dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1532dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1533a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
15344097062eSBarry Smith #endif
15352eff7a51SBarry Smith 
153633f85c2fSBarry Smith   /*
153733f85c2fSBarry Smith      Destroy any packages that registered a finalize
153833f85c2fSBarry Smith   */
153933f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
154033f85c2fSBarry Smith 
1541101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1542fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1543101409b8SToby Isaac #endif
1544101409b8SToby Isaac 
154533f85c2fSBarry Smith   /*
154648dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
154748dd1dffSBarry Smith 
154837e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
154948dd1dffSBarry Smith   */
155037e93019SBarry Smith 
15514028d114SSatish Balay   if (petsc_history) {
1552f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
155302c9f0b5SLisandro Dalcin     petsc_history = NULL;
1554e5c89e4eSSatish Balay   }
15559de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1556e94e781bSJacob Faibussowitsch   ierr = PetscInfoDestroy();CHKERRQ(ierr);
1557e5c89e4eSSatish Balay 
155867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
155992f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1560e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
156192f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1562e5c89e4eSSatish Balay     FILE *fd;
1563ed9cf6e9SBarry Smith     int  err;
1564e5c89e4eSSatish Balay 
1565dc92acbaSJed Brown     flg2 = PETSC_FALSE;
156692f119d6SBarry Smith     flg3 = PETSC_FALSE;
1567cf9c20a2SJed Brown     if (PetscDefined(USE_DEBUG)) {ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);}
156892f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
156992f119d6SBarry Smith     fname[0] = 0;
1570589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1571e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1572e5c89e4eSSatish Balay 
1573589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1574e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1575e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1576ed9cf6e9SBarry Smith       err  = fclose(fd);
1577e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
157892f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1579e5c89e4eSSatish Balay       MPI_Comm local_comm;
1580e5c89e4eSSatish Balay 
1581ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(ierr);
1582e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1583e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1584e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1585ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(ierr);
1586e5c89e4eSSatish Balay     }
1587e5c89e4eSSatish Balay     fname[0] = 0;
1588589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1589e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1590e5c89e4eSSatish Balay 
1591589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
159292f119d6SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
159392f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1594ed9cf6e9SBarry Smith       err  = fclose(fd);
1595e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
159692f119d6SBarry Smith     } else if (flg1) {
159792f119d6SBarry Smith       MPI_Comm local_comm;
159892f119d6SBarry Smith 
1599ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(ierr);
160092f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
160192f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
160292f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1603ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(ierr);
1604e5c89e4eSSatish Balay     }
1605e5c89e4eSSatish Balay   }
160667584ceeSBarry Smith #endif
160720e2c332SMatthew G. Knepley 
16085486ca60SMatthew G. Knepley   /*
16095486ca60SMatthew G. Knepley      Close any open dynamic libraries
16105486ca60SMatthew G. Knepley   */
16115486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
16125486ca60SMatthew G. Knepley 
1613e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
16144416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1615e5c89e4eSSatish Balay 
1616e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
161702c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1618e5c89e4eSSatish Balay 
1619c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1620c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
1621c2b86a48SJunchao Zhang     ierr = PetscKokkosFinalize_Private();CHKERRQ(ierr);
1622c2b86a48SJunchao Zhang     PetscBeganKokkos = PETSC_FALSE;
162345639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1624c2b86a48SJunchao Zhang   }
1625c2b86a48SJunchao Zhang #endif
1626c2b86a48SJunchao Zhang 
162771438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM)
162871438e86SJunchao Zhang   if (PetscBeganNvshmem) {
162971438e86SJunchao Zhang     ierr = PetscNvshmemFinalize();CHKERRQ(ierr);
163071438e86SJunchao Zhang     PetscBeganNvshmem = PETSC_FALSE;
163171438e86SJunchao Zhang   }
163271438e86SJunchao Zhang #endif
1633a0e72f99SJunchao Zhang 
1634a0e72f99SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1635a0e72f99SJunchao Zhang   if (PetscDefaultCudaStream) {cudaError_t cerr = cudaStreamDestroy(PetscDefaultCudaStream);CHKERRCUDA(cerr);}
16369ffd0706SHong Zhang   if (petsc_gputimer_begin) {
16379ffd0706SHong Zhang     cudaError_t cerr = cudaEventDestroy(petsc_gputimer_begin);CHKERRCUDA(cerr);
16389ffd0706SHong Zhang   }
16399ffd0706SHong Zhang   if (petsc_gputimer_end) {
16409ffd0706SHong Zhang     cudaError_t cerr = cudaEventDestroy(petsc_gputimer_end);CHKERRCUDA(cerr);
16419ffd0706SHong Zhang   }
1642a0e72f99SJunchao Zhang #endif
1643a0e72f99SJunchao Zhang 
1644a0e72f99SJunchao Zhang #if defined(PETSC_HAVE_HIP)
1645a0e72f99SJunchao Zhang   if (PetscDefaultHipStream)  {hipError_t cerr  = hipStreamDestroy(PetscDefaultHipStream);CHKERRHIP(cerr);}
16469ffd0706SHong Zhang   if (petsc_gputimer_begin) {
16479ffd0706SHong Zhang     hipError_t cerr = hipEventDestroy(petsc_gputimer_begin);CHKERRHIP(cerr);
16489ffd0706SHong Zhang   }
16499ffd0706SHong Zhang   if (petsc_gputimer_end) {
16509ffd0706SHong Zhang     hipError_t cerr = hipEventDestroy(petsc_gputimer_end);CHKERRHIP(cerr);
16519ffd0706SHong Zhang   }
1652a0e72f99SJunchao Zhang #endif
1653a0e72f99SJunchao Zhang 
1654008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1655e5c89e4eSSatish Balay 
1656dbc8283eSBarry Smith   /*
1657efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1658efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1659efb80d3cSBarry Smith 
1660efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1661efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1662dbc8283eSBarry Smith  */
1663b770b1f6SSatish Balay   {
1664dbc8283eSBarry Smith     PetscCommCounter *counter;
1665dbc8283eSBarry Smith     PetscMPIInt      flg;
1666dbc8283eSBarry Smith     MPI_Comm         icomm;
1667265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
1668ffc4695bSBarry Smith     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRMPI(ierr);
1669dbc8283eSBarry Smith     if (flg) {
1670265f3f35SJed Brown       icomm = ucomm.comm;
1671ffc4695bSBarry Smith       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRMPI(ierr);
167227104ee2SJacob Faibussowitsch       if (PetscUnlikely(!flg)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1673dbc8283eSBarry Smith 
1674ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1675ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRMPI(ierr);
1676ffc4695bSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRMPI(ierr);
1677dbc8283eSBarry Smith     }
1678ffc4695bSBarry Smith     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRMPI(ierr);
1679dbc8283eSBarry Smith     if (flg) {
1680265f3f35SJed Brown       icomm = ucomm.comm;
1681ffc4695bSBarry Smith       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRMPI(ierr);
168227104ee2SJacob Faibussowitsch       if (PetscUnlikely(!flg)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1683dbc8283eSBarry Smith 
1684ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1685ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRMPI(ierr);
1686ffc4695bSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRMPI(ierr);
1687dbc8283eSBarry Smith     }
1688b770b1f6SSatish Balay   }
1689dbc8283eSBarry Smith 
1690ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRMPI(ierr);
1691ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1692ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRMPI(ierr);
1693ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRMPI(ierr);
1694480cf27aSJed Brown 
16955ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
16965ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
16975ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1698ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1699ef19f930SBarry Smith 
1700e5c89e4eSSatish Balay   if (PetscBeganMPI) {
170199608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
170299b1327fSBarry Smith     PetscMPIInt flag;
1703ffc4695bSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRMPI(ierr);
170427104ee2SJacob Faibussowitsch     if (PetscUnlikely(flag)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
170599608316SBarry Smith #endif
1706ffc4695bSBarry Smith     ierr = MPI_Finalize();CHKERRMPI(ierr);
1707e5c89e4eSSatish Balay   }
1708e5c89e4eSSatish Balay /*
1709e5c89e4eSSatish Balay 
1710e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1711e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1712e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1713e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1714e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
17150e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1716e5c89e4eSSatish Balay    memory was not freed.
1717e5c89e4eSSatish Balay 
1718e5c89e4eSSatish Balay */
17191d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
172027104ee2SJacob Faibussowitsch   ierr = PetscStackReset();CHKERRQ(ierr);
1721a297a907SKarl Rupp 
17228ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_FALSE;
1723e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1724e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
172556883f60SBarry Smith #if defined(PETSC_USE_GCOV)
172656883f60SBarry Smith   /*
172756883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
172856883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
172956883f60SBarry Smith    */
173056883f60SBarry Smith   __gcov_flush();
173156883f60SBarry Smith #endif
173227104ee2SJacob Faibussowitsch   return 0;
1733e5c89e4eSSatish Balay }
1734e5c89e4eSSatish Balay 
173543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
17368cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
173743db4dbbSBarry Smith {
173843db4dbbSBarry Smith   if (*a == *b) return 1;
173943db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
174043db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
174143db4dbbSBarry Smith   return 0;
174243db4dbbSBarry Smith }
1743a70650f6SBarry Smith #endif
174443db4dbbSBarry Smith 
174543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
17468cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
174743db4dbbSBarry Smith {
174843db4dbbSBarry Smith   if (*a == *b) return 1;
174943db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
175043db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
175143db4dbbSBarry Smith   return 0;
175243db4dbbSBarry Smith }
175343db4dbbSBarry Smith #endif
1754