xref: /petsc/src/sys/objects/pinit.c (revision cf9c20a2d58da010f7c4712defbcdf61cc8f72b5)
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*/
6022afb99SBarry Smith #include <petscvalgrind.h>
7665c2dedSJed Brown #include <petscviewer.h>
88101f56cSMatthew Knepley 
9a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
11a9f03627SSatish Balay #endif
12f2d66bcaSShri Abhyankar 
132d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1495c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
152d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
162d53ad75SBarry Smith #endif
172d53ad75SBarry Smith 
18a6790183SBarry Smith #if defined(PETSC_HAVE_SAWS)
1916ad0300SBarry Smith #include <petscviewersaws.h>
20a6790183SBarry Smith #endif
21eb02082bSJunchao Zhang 
22e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
23e5c89e4eSSatish Balay 
2495c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
25e5c89e4eSSatish Balay 
2695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
2795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
2895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
2995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
3195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
320069ddf5SShri Abhyankar 
33e5c89e4eSSatish Balay /* user may set this BEFORE calling PetscInitialize() */
34e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
35e5c89e4eSSatish Balay 
36480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
37480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
38480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
395f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
40480cf27aSJed Brown 
41e5c89e4eSSatish Balay /*
42e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
43e5c89e4eSSatish Balay */
4402c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
4502c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
46e5c89e4eSSatish Balay 
47ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
48ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
490f8e0872SSatish Balay 
50a2f94806SJed Brown PetscInt PetscHotRegionDepth;
51a2f94806SJed Brown 
52b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
53b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
54b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
55b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
56b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
57b22622e2STadeu Manoel #endif
58b22622e2STadeu Manoel 
59e5c89e4eSSatish Balay /*
60945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
6172a42c3cSBarry Smith 
6272a42c3cSBarry Smith    Collective
6372a42c3cSBarry Smith 
6472a42c3cSBarry Smith    Level: advanced
6572a42c3cSBarry Smith 
6695452b02SPatrick Sanan     Notes:
67a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
680f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
69a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
700f11a792SBarry Smith 
71a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
721ea5a559SBarry Smith 
7372a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
740f11a792SBarry Smith */
75945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
7672a42c3cSBarry Smith {
7772a42c3cSBarry Smith   PetscErrorCode ierr;
7872a42c3cSBarry Smith   int            myargc   = argc;
7972a42c3cSBarry Smith   char           **myargs = args;
8072a42c3cSBarry Smith 
8172a42c3cSBarry Smith   PetscFunctionBegin;
82c52ac268SRichard Tran Mills   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr;
831ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
84df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
8572a42c3cSBarry Smith   PetscFunctionReturn(ierr);
8672a42c3cSBarry Smith }
8772a42c3cSBarry Smith 
88f0865b08SBarry Smith /*
89a56f64adSBarry Smith       Used by Julia interface to get communicator
90f0865b08SBarry Smith */
91945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
92f0865b08SBarry Smith {
93f0865b08SBarry Smith   PetscFunctionBegin;
94f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
95f0865b08SBarry Smith   PetscFunctionReturn(0);
96f0865b08SBarry Smith }
97f0865b08SBarry Smith 
98e5c89e4eSSatish Balay /*@C
99e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
100e5c89e4eSSatish Balay         the command line arguments.
101e5c89e4eSSatish Balay 
102e5c89e4eSSatish Balay    Collective
103e5c89e4eSSatish Balay 
104e5c89e4eSSatish Balay    Level: advanced
105e5c89e4eSSatish Balay 
106e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
107e5c89e4eSSatish Balay @*/
1087087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
109e5c89e4eSSatish Balay {
110e5c89e4eSSatish Balay   PetscErrorCode ierr;
111e5c89e4eSSatish Balay   int            argc   = 0;
11202c9f0b5SLisandro Dalcin   char           **args = NULL;
113e5c89e4eSSatish Balay 
114e5c89e4eSSatish Balay   PetscFunctionBegin;
1150298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
116e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
117e5c89e4eSSatish Balay }
118e5c89e4eSSatish Balay 
119e5c89e4eSSatish Balay /*@
120e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
121e5c89e4eSSatish Balay 
12293b6d2d1SJed Brown    Level: beginner
123e5c89e4eSSatish Balay 
124e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
125e5c89e4eSSatish Balay @*/
1267087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool  *isInitialized)
127e5c89e4eSSatish Balay {
128e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
12993b6d2d1SJed Brown   return 0;
130e5c89e4eSSatish Balay }
131e5c89e4eSSatish Balay 
132e5c89e4eSSatish Balay /*@
133e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
134e5c89e4eSSatish Balay 
135e5c89e4eSSatish Balay    Level: developer
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
138e5c89e4eSSatish Balay @*/
1397087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
140e5c89e4eSSatish Balay {
141e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
14293b6d2d1SJed Brown   return 0;
143e5c89e4eSSatish Balay }
144e5c89e4eSSatish Balay 
14595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void);
146e5c89e4eSSatish Balay 
147e5c89e4eSSatish Balay /*
148e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
149e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
150e5c89e4eSSatish Balay */
151367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
152e5c89e4eSSatish Balay 
153367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
154e5c89e4eSSatish Balay {
155e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
156e5c89e4eSSatish Balay 
157e5c89e4eSSatish Balay   PetscFunctionBegin;
158e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
159e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
16041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
161e5c89e4eSSatish Balay   }
162e5c89e4eSSatish Balay 
163e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
164e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
165e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
166e5c89e4eSSatish Balay   }
167812af9f3SBarry Smith   PetscFunctionReturnVoid();
168e5c89e4eSSatish Balay }
169e5c89e4eSSatish Balay 
170e5c89e4eSSatish Balay /*
171e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
172e5c89e4eSSatish Balay sum of the second entry.
173b693b147SBarry Smith 
17476f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
175367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
176b693b147SBarry Smith there would be no place to store the both needed results.
177e5c89e4eSSatish Balay */
17876ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
179e5c89e4eSSatish Balay {
180e5c89e4eSSatish Balay   PetscErrorCode ierr;
181e5c89e4eSSatish Balay 
182e5c89e4eSSatish Balay   PetscFunctionBegin;
183d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
184d6e4c47cSJed Brown   {
185d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
186367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
187d6e4c47cSJed Brown     *max = work.max;
188d6e4c47cSJed Brown     *sum = work.sum;
189d6e4c47cSJed Brown   }
190d6e4c47cSJed Brown #else
191d6e4c47cSJed Brown   {
192d6e4c47cSJed Brown     PetscMPIInt    size,rank;
193d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
194e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
195e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
196785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
197367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
1986ac3741eSJed Brown     *max = work[rank].max;
1996ac3741eSJed Brown     *sum = work[rank].sum;
200e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
201d6e4c47cSJed Brown   }
202d6e4c47cSJed Brown #endif
203e5c89e4eSSatish Balay   PetscFunctionReturn(0);
204e5c89e4eSSatish Balay }
205e5c89e4eSSatish Balay 
206e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
207e5c89e4eSSatish Balay 
208570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
20906a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
210e5c89e4eSSatish Balay 
2118cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
212e5c89e4eSSatish Balay {
213e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
214e5c89e4eSSatish Balay 
215e5c89e4eSSatish Balay   PetscFunctionBegin;
2167c2de775SJed Brown   if (*datatype == MPIU_REAL) {
217e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
218a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2197c2de775SJed Brown   }
2207c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2217c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2227c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
223a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2247c2de775SJed Brown   }
2257c2de775SJed Brown #endif
2267c2de775SJed Brown   else {
2277c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
22841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
229e2e03761SBarry Smith   }
230812af9f3SBarry Smith   PetscFunctionReturnVoid();
231e5c89e4eSSatish Balay }
232e5c89e4eSSatish Balay #endif
233e5c89e4eSSatish Balay 
234570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
235d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
236d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
237d9822059SBarry Smith 
2388cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
239d9822059SBarry Smith {
240d9822059SBarry Smith   PetscInt i,count = *cnt;
241d9822059SBarry Smith 
242d9822059SBarry Smith   PetscFunctionBegin;
2437c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2448c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
245a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2467c2de775SJed Brown   }
2477c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2487c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2497c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2507c2de775SJed Brown     for (i=0; i<count; i++) {
2517c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2527c2de775SJed Brown     }
2537c2de775SJed Brown   }
2547c2de775SJed Brown #endif
2557c2de775SJed Brown   else {
2567c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
25741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2588c764dc5SJose Roman   }
259d9822059SBarry Smith   PetscFunctionReturnVoid();
260d9822059SBarry Smith }
261d9822059SBarry Smith 
2628cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
263d9822059SBarry Smith {
264d9822059SBarry Smith   PetscInt    i,count = *cnt;
265d9822059SBarry Smith 
266d9822059SBarry Smith   PetscFunctionBegin;
2677c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2688c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
269a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2707c2de775SJed Brown   }
2717c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2727c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2737c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2747c2de775SJed Brown     for (i=0; i<count; i++) {
2757c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2767c2de775SJed Brown     }
2777c2de775SJed Brown   }
2787c2de775SJed Brown #endif
2797c2de775SJed Brown   else {
2808c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
28141e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2828c764dc5SJose Roman   }
283d9822059SBarry Smith   PetscFunctionReturnVoid();
284d9822059SBarry Smith }
285d9822059SBarry Smith #endif
286d9822059SBarry Smith 
287480cf27aSJed Brown /*
288480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
289480cf27aSJed Brown 
290ff0e51ddSBarry 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.
291480cf27aSJed Brown 
29212801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
293480cf27aSJed Brown 
294480cf27aSJed Brown */
29533779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
296480cf27aSJed Brown {
297480cf27aSJed Brown   PetscErrorCode   ierr;
29805342407SJunchao Zhang   PetscCommCounter *counter=(PetscCommCounter*)count_val;
299480cf27aSJed Brown 
300480cf27aSJed Brown   PetscFunctionBegin;
30102c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
30205342407SJunchao Zhang   ierr = PetscFree(counter->iflags);CHKERRMPI(ierr);
30305342407SJunchao Zhang   ierr = PetscFree(counter);CHKERRMPI(ierr);
304480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
305480cf27aSJed Brown }
306480cf27aSJed Brown 
307480cf27aSJed Brown /*
30847435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
309da3039f7SJed Brown   calls MPI_Comm_free().
310da3039f7SJed Brown 
311da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
312480cf27aSJed Brown 
313ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
314480cf27aSJed Brown 
31512801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
316480cf27aSJed Brown 
317480cf27aSJed Brown */
31833779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
319480cf27aSJed Brown {
320480cf27aSJed Brown   PetscErrorCode                    ierr;
32133779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
322480cf27aSJed Brown 
323480cf27aSJed Brown   PetscFunctionBegin;
32412801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
325ec4fadc2SJed Brown   icomm.ptr = attr_val;
32676bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
32733779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
32833779a13SJunchao Zhang     PetscMPIInt flg;
32933779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
33047435625SJed Brown     ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
33133779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
33233779a13SJunchao 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");
33333779a13SJunchao Zhang   }
33447435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
33533779a13SJunchao Zhang   ierr = PetscInfo2(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);CHKERRMPI(ierr);
336da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
337b89831e5SBarry Smith }
338da3039f7SJed Brown 
339da3039f7SJed Brown /*
34033779a13SJunchao 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.
341da3039f7SJed Brown  */
34233779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
343da3039f7SJed Brown {
344da3039f7SJed Brown   PetscErrorCode ierr;
345da3039f7SJed Brown 
346da3039f7SJed Brown   PetscFunctionBegin;
34702c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
348480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
349480cf27aSJed Brown }
350480cf27aSJed Brown 
35133779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
35242218b76SBarry Smith 
353951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3548cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3558cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3568cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
357e39fd77fSBarry Smith #endif
358e39fd77fSBarry Smith 
3590f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
36012801b39SBarry Smith 
361eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
362eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3636ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
36402c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
365dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
366e5c89e4eSSatish Balay 
367dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
368051e4cf2SJed Brown {
369051e4cf2SJed Brown   PetscErrorCode ierr;
370051e4cf2SJed Brown 
371051e4cf2SJed Brown   PetscFunctionBegin;
372051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
373a1601671SBarry 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);
374051e4cf2SJed 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);
375051e4cf2SJed Brown   PetscFunctionReturn(0);
376051e4cf2SJed Brown }
377e5c89e4eSSatish Balay 
3782d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3792d747510SLisandro Dalcin 
3802d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
3812d747510SLisandro Dalcin {
3822d747510SLisandro Dalcin   PetscErrorCode ierr;
3832d747510SLisandro Dalcin 
3842d747510SLisandro Dalcin   PetscFunctionBegin;
3852d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
3862d747510SLisandro Dalcin   PetscFunctionReturn(0);
3872d747510SLisandro Dalcin }
3882d747510SLisandro Dalcin 
3892d747510SLisandro Dalcin /*@C
3902d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
3912d747510SLisandro Dalcin 
3922d747510SLisandro Dalcin     Not Collective
3932d747510SLisandro Dalcin 
3942d747510SLisandro Dalcin     Input Parameter:
3952d747510SLisandro Dalcin .   len - length of the string name
3962d747510SLisandro Dalcin 
3972d747510SLisandro Dalcin     Output Parameter:
3982d747510SLisandro Dalcin .   name - the name of the running program
3992d747510SLisandro Dalcin 
4002d747510SLisandro Dalcin    Level: advanced
4012d747510SLisandro Dalcin 
4022d747510SLisandro Dalcin     Notes:
4032d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4042d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4052d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4062d747510SLisandro Dalcin @*/
4072d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4082d747510SLisandro Dalcin {
4092d747510SLisandro Dalcin   PetscErrorCode ierr;
4102d747510SLisandro Dalcin 
4112d747510SLisandro Dalcin   PetscFunctionBegin;
4122d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4132d747510SLisandro Dalcin   PetscFunctionReturn(0);
4142d747510SLisandro Dalcin }
4152d747510SLisandro Dalcin 
416e5c89e4eSSatish Balay /*@C
417e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
418e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
419e5c89e4eSSatish Balay 
420e5c89e4eSSatish Balay    Not Collective
421e5c89e4eSSatish Balay 
422e5c89e4eSSatish Balay    Output Parameters:
423e5c89e4eSSatish Balay +  argc - count of number of command line arguments
424e5c89e4eSSatish Balay -  args - the command line arguments
425e5c89e4eSSatish Balay 
426e5c89e4eSSatish Balay    Level: intermediate
427e5c89e4eSSatish Balay 
428e5c89e4eSSatish Balay    Notes:
429e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
430e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
431e5c89e4eSSatish Balay 
432f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
433f177e3b1SBarry Smith 
434793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
435e5c89e4eSSatish Balay 
436e5c89e4eSSatish Balay @*/
4377087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
438e5c89e4eSSatish Balay {
439e5c89e4eSSatish Balay   PetscFunctionBegin;
44017186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
441e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
442e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
443e5c89e4eSSatish Balay   PetscFunctionReturn(0);
444e5c89e4eSSatish Balay }
445e5c89e4eSSatish Balay 
446793721a6SBarry Smith /*@C
447793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
448793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
449793721a6SBarry Smith 
450793721a6SBarry Smith    Not Collective
451793721a6SBarry Smith 
452793721a6SBarry Smith    Output Parameters:
453793721a6SBarry Smith .  args - the command line arguments
454793721a6SBarry Smith 
455793721a6SBarry Smith    Level: intermediate
456793721a6SBarry Smith 
457793721a6SBarry Smith    Notes:
458793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
459793721a6SBarry Smith 
460793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
461793721a6SBarry Smith 
462793721a6SBarry Smith @*/
4637087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
464793721a6SBarry Smith {
465793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
466793721a6SBarry Smith   PetscErrorCode ierr;
467793721a6SBarry Smith 
468793721a6SBarry Smith   PetscFunctionBegin;
46917186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4702d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
471785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
472793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
473793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
474793721a6SBarry Smith   }
4752d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
476793721a6SBarry Smith   PetscFunctionReturn(0);
477793721a6SBarry Smith }
478793721a6SBarry Smith 
479793721a6SBarry Smith /*@C
480793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
481793721a6SBarry Smith 
482793721a6SBarry Smith    Not Collective
483793721a6SBarry Smith 
484793721a6SBarry Smith    Output Parameters:
485793721a6SBarry Smith .  args - the command line arguments
486793721a6SBarry Smith 
487793721a6SBarry Smith    Level: intermediate
488793721a6SBarry Smith 
489793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
490793721a6SBarry Smith 
491793721a6SBarry Smith @*/
4927087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
493793721a6SBarry Smith {
494793721a6SBarry Smith   PetscInt       i = 0;
495793721a6SBarry Smith   PetscErrorCode ierr;
496793721a6SBarry Smith 
497793721a6SBarry Smith   PetscFunctionBegin;
498a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
499793721a6SBarry Smith   while (args[i]) {
500793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
501793721a6SBarry Smith     i++;
502793721a6SBarry Smith   }
503793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
504793721a6SBarry Smith   PetscFunctionReturn(0);
505793721a6SBarry Smith }
506793721a6SBarry Smith 
50711525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
50830befbd2SBarry Smith #include <petscconfiginfo.h>
50930befbd2SBarry Smith 
51095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
51111525c0dSBarry Smith {
51211525c0dSBarry Smith   if (!PetscGlobalRank) {
51330befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
51411525c0dSBarry Smith     int            port;
515ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
51611525c0dSBarry Smith     size_t         applinelen,introlen;
51711525c0dSBarry Smith     PetscErrorCode ierr;
518ffbd1cfbSBarry Smith     char           sawsurl[256];
51911525c0dSBarry Smith 
520c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
52111525c0dSBarry Smith     if (flg) {
52211525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
52311525c0dSBarry Smith 
524c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
52511525c0dSBarry Smith       if (sawslog[0]) {
52611525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
52711525c0dSBarry Smith       } else {
52811525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
52911525c0dSBarry Smith       }
53011525c0dSBarry Smith     }
531c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
53211525c0dSBarry Smith     if (flg) {
53311525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
53411525c0dSBarry Smith     }
535c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
536ffbd1cfbSBarry Smith     if (selectport) {
537ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
538ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
539ffbd1cfbSBarry Smith     } else {
540c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
54111525c0dSBarry Smith       if (flg) {
54211525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
54311525c0dSBarry Smith       }
544ffbd1cfbSBarry Smith     }
545c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
54611525c0dSBarry Smith     if (flg) {
54711525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
54811525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5499c1e0ce8SBarry Smith     } else {
550c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5519c1e0ce8SBarry Smith       if (flg) {
5523c01dfcfSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
5539c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5549c1e0ce8SBarry Smith       }
55511525c0dSBarry Smith     }
556c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
55711525c0dSBarry Smith     if (flg2) {
55811525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
55911525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
56011525c0dSBarry Smith       ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr);
56111525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
56211525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
56343da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
56411525c0dSBarry Smith     }
56511525c0dSBarry Smith     ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr);
56611525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
56711525c0dSBarry Smith     introlen   = 4096 + applinelen;
56830a8c9c0SSurtai Han     applinelen += 1024;
56911525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
57011525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
57111525c0dSBarry Smith 
57211525c0dSBarry Smith     if (rootlocal) {
57311525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
57411525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
57511525c0dSBarry Smith     }
57676a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
57711525c0dSBarry Smith     if (rootlocal && help) {
578928bb9adSStefano 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);
57911525c0dSBarry Smith     } else if (help) {
580928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
58111525c0dSBarry Smith     } else {
582928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
58311525c0dSBarry Smith     }
584b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
58530befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
58611525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
587a8d69d7bSBarry Smith                                     "<center><h2> <a href=\"https://www.mcs.anl.gov/petsc\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n"
588df62c222SBarry 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"
589928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
59043da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
59111525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
59211525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
593ffbd1cfbSBarry Smith     if (selectport) {
594aa573868SBarry Smith       PetscBool silent;
5957d812c46SBarry Smith 
5967d812c46SBarry Smith       ierr = SAWs_Initialize();
5977d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
5987d812c46SBarry Smith       while (ierr) {
5997d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6007d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6017d812c46SBarry Smith         ierr = SAWs_Initialize();
6027d812c46SBarry Smith       }
6037d812c46SBarry Smith 
604aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
605aa573868SBarry Smith       if (!silent) {
606ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
607ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
608ffbd1cfbSBarry Smith       }
6097d812c46SBarry Smith     } else {
6107d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
611aa573868SBarry Smith     }
6120af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6130af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6140af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6150af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6160af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
61711525c0dSBarry Smith   }
61811525c0dSBarry Smith   PetscFunctionReturn(0);
61911525c0dSBarry Smith }
62011525c0dSBarry Smith #endif
62111525c0dSBarry Smith 
6224dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6234dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6244dfee713SSatish Balay {
6254dfee713SSatish Balay   PetscFunctionBegin;
6264dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6274dfee713SSatish Balay     /* see MPI.py for details on this bug */
6284dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6294dfee713SSatish Balay #endif
6304dfee713SSatish Balay   PetscFunctionReturn(0);
6314dfee713SSatish Balay }
6324dfee713SSatish Balay 
633a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
634a56f64adSBarry Smith #include <adios.h>
63522580e64SBarry Smith #include <adios_read.h>
6367b56e58cSSatish Balay int64_t Petsc_adios_group;
637a56f64adSBarry Smith #endif
63855d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
63955d657eeSBarry Smith #include <adios2_c.h>
64055d657eeSBarry Smith #endif
641cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
642cd1b99a6SBarry Smith #include <omp.h>
643f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
644cd1b99a6SBarry Smith #endif
645a56f64adSBarry Smith 
646bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H)
647bc8a24ecSBarry Smith #include <dlfcn.h>
648bc8a24ecSBarry Smith #endif
649bc8a24ecSBarry Smith 
650e5c89e4eSSatish Balay /*@C
651e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
652e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
653e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
654e5c89e4eSSatish Balay    your program -- usually the very first line!
655e5c89e4eSSatish Balay 
656e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
657e5c89e4eSSatish Balay 
658e5c89e4eSSatish Balay    Input Parameters:
659e5c89e4eSSatish Balay +  argc - count of number of command line arguments
660e5c89e4eSSatish Balay .  args - the command line arguments
6610298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
662fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6630298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
664e5c89e4eSSatish Balay 
66505827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
66605827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
66705827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
66805827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
66905827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
670e5c89e4eSSatish Balay 
671e5c89e4eSSatish Balay    Options Database Keys:
6727ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6737ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
674e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
6758a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
6768a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
6778a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
6788a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
6798a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
680e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
681e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
682e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
683e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
68479dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
68579dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
68679dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
687aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
68892f119d6SBarry 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
68992f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
69092f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
69192f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
692e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
693e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
694e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
695e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
696e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
6970841954dSBarry Smith -  -memory_view - Print memory usage at end of run
698e5c89e4eSSatish Balay 
699e5c89e4eSSatish Balay    Options Database Keys for Profiling:
700a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
701fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
702217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
703217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
704495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
705e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7069a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
70779dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
7089a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
709495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
710f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
711495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
712495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
713c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
71487c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
715c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
716495fc317SBarry Smith 
717609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
718e5c89e4eSSatish Balay 
719ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
720ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
721ffbd1cfbSBarry 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
722ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
723ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
724ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
725ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
726ffbd1cfbSBarry Smith 
727e5c89e4eSSatish Balay    Environmental Variables:
728e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
729e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
730e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
731e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
732e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
733e5c89e4eSSatish Balay 
734e5c89e4eSSatish Balay 
735e5c89e4eSSatish Balay    Level: beginner
736e5c89e4eSSatish Balay 
737e5c89e4eSSatish Balay    Notes:
738e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
739e5c89e4eSSatish Balay    it before PetscInitialize().
740e5c89e4eSSatish Balay 
741e5c89e4eSSatish Balay    Fortran Version:
742e5c89e4eSSatish Balay    In Fortran this routine has the format
743e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
744e5c89e4eSSatish Balay 
745e5c89e4eSSatish Balay +   ierr - error return code
7460eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
747fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
748e5c89e4eSSatish Balay 
749e5c89e4eSSatish Balay    Important Fortran Note:
7500eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7510298fd71SBarry Smith    null character string; you CANNOT just use NULL as
752a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
753e5c89e4eSSatish Balay 
75401cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
75501cb0274SBarry Smith    calling PetscInitialize().
756e5c89e4eSSatish Balay 
75701cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
758e5c89e4eSSatish Balay 
759e5c89e4eSSatish Balay @*/
7607087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
761e5c89e4eSSatish Balay {
762e5c89e4eSSatish Balay   PetscErrorCode ierr;
7634bb5149bSJed Brown   PetscMPIInt    flag, size;
76419c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
765e5c89e4eSSatish Balay   char           hostname[256];
766e5c89e4eSSatish Balay 
767e5c89e4eSSatish Balay   PetscFunctionBegin;
768e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7693d96e996SBarry Smith   /*
7703d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7713d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7723d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7733d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7743d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7753d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7763d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7773d96e996SBarry Smith       listed above and since that time are compatible.
7783d96e996SBarry Smith 
7793d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
7803d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
7813d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
7823d96e996SBarry Smith       to perform the checking.
7833d96e996SBarry Smith 
7843d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
7853d96e996SBarry Smith 
7863d96e996SBarry Smith       Questions:
7873d96e996SBarry Smith 
7883d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
7893d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
7903d96e996SBarry Smith   */
7913d96e996SBarry Smith 
79219c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
79319c5658aSBarry Smith   {
79419c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
79519c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
79619c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
7973d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
79819c5658aSBarry Smith #if defined(MPICH_VERSION)
7993d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
80019c5658aSBarry Smith     {
80119c5658aSBarry Smith       char *ver,*lf;
80219c5658aSBarry Smith       flg = PETSC_FALSE;
80319c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
80419c5658aSBarry Smith       if (ver) {
80519c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
80619c5658aSBarry Smith         if (lf) {
80719c5658aSBarry Smith           *lf = 0;
80819c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
80919c5658aSBarry Smith         }
81019c5658aSBarry Smith       }
81119c5658aSBarry Smith       if (!flg) {
81219c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8133d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
81419c5658aSBarry Smith       }
81519c5658aSBarry Smith     }
8163d96e996SBarry Smith #endif
8173d96e996SBarry Smith     /* check for OpenMPI version, it is not part of the MPI ABI initiative (is it part of another initiative that needs to be handled?) */
81819c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
81919c5658aSBarry Smith     {
82019c5658aSBarry Smith       char *ver,bs[32],*bsf;
82119c5658aSBarry Smith       flg = PETSC_FALSE;
82219c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
82319c5658aSBarry Smith       if (ver) {
8242e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
82519c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
82619c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
82719c5658aSBarry Smith       }
82819c5658aSBarry Smith       if (!flg) {
82919c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d, aborting\n",mpilibraryversion,OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
8303d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
83119c5658aSBarry Smith       }
83219c5658aSBarry Smith     }
83319c5658aSBarry Smith #endif
83419c5658aSBarry Smith   }
83519c5658aSBarry Smith #endif
83619c5658aSBarry Smith 
837bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM)
838bc8a24ecSBarry Smith   {
839bc8a24ecSBarry Smith     PetscInt cnt = 0;
840bc8a24ecSBarry Smith     /* 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 */
841bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++;
842bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++;
843bc8a24ecSBarry Smith     if (cnt > 1) {
844bc8a24ecSBarry Smith       fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
845bc8a24ecSBarry Smith       return PETSC_ERR_MPI_LIB_INCOMP;
846bc8a24ecSBarry Smith     }
847bc8a24ecSBarry Smith   }
848bc8a24ecSBarry Smith #endif
849bc8a24ecSBarry Smith 
850ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
851d89683f4Sbcordonn   PETSC_STDOUT = stdout;
852ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
853e5c89e4eSSatish Balay 
8540c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8550c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8560c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8570c30907bSSatish Balay #endif
8580c30907bSSatish Balay 
8594416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
860e5c89e4eSSatish Balay 
861e5c89e4eSSatish Balay   /*
862e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
863e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
864e5c89e4eSSatish Balay   */
865e5c89e4eSSatish Balay   if (argc && *argc) {
866e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
867e5c89e4eSSatish Balay   } else {
868e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
869e5c89e4eSSatish Balay   }
870e5c89e4eSSatish Balay 
871e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
872e5c89e4eSSatish Balay   if (!flag) {
873e32f2f54SBarry 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");
8744dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
8755e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8765e765c61SJed Brown     {
8775e765c61SJed Brown       PetscMPIInt provided;
8785e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
8795e765c61SJed Brown     }
8805e765c61SJed Brown #else
881e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
8825e765c61SJed Brown #endif
883e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
884e5c89e4eSSatish Balay   }
8854dfee713SSatish Balay 
886e5c89e4eSSatish Balay   if (argc && args) {
887e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
888e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
889e5c89e4eSSatish Balay   }
890e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
8915ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
8925ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
8935ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
894ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
895e5c89e4eSSatish Balay 
896a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
897d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
898e5c89e4eSSatish Balay 
8990f9be574SPeter Hill   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
90012801b39SBarry Smith     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
90112801b39SBarry Smith     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
9020f9be574SPeter Hill   }
90312801b39SBarry Smith 
904e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
905e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
906e5c89e4eSSatish Balay 
907e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
908e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
909e5c89e4eSSatish Balay 
9108ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9118ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9127cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
913e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
914e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
915e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
916e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
917e316c87fSJed Brown #endif
918e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9198ad47952SJed Brown 
920e5c89e4eSSatish Balay   /*
921e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
922e5c89e4eSSatish Balay      shared libraries the constructors for global variables
923e5c89e4eSSatish Balay      are not called; at least on IRIX.
924e5c89e4eSSatish Balay   */
925886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
926e5c89e4eSSatish Balay   {
927252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
92850f81f78SJed Brown     PetscComplex ic(0.0,1.0);
929e5c89e4eSSatish Balay     PETSC_i = ic;
930252ecd31SSatish Balay #else
93150f81f78SJed Brown     PETSC_i = _Complex_I;
932b7940d39SSatish Balay #endif
933762437b8SSatish Balay   }
934762437b8SSatish Balay 
9352c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
936e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
937500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
938500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
939500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9402c876bd9SBarry Smith #endif
941886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
942e5c89e4eSSatish Balay 
943e5c89e4eSSatish Balay   /*
944e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
945e5c89e4eSSatish Balay      half of the entries and maxes the second half.
946e5c89e4eSSatish Balay   */
947367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
948e5c89e4eSSatish Balay 
949ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
950c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
951c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9527c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9538c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9548c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9558c764dc5SJose Roman #endif
956d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
957d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
958570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
959570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
960570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
961570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
962570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
963c90a1750SBarry Smith #endif
964c90a1750SBarry Smith 
965570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
966cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
967cca4cb22SSatish Balay #endif
968cca4cb22SSatish Balay 
969e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
970e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
971e5c89e4eSSatish Balay 
97240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
973e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
974e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
97544041f26SJed Brown #endif
976e5c89e4eSSatish Balay 
977ec957eceSBarry Smith 
978e5c89e4eSSatish Balay   /*
979480cf27aSJed Brown      Attributes to be set on PETSc communicators
980480cf27aSJed Brown   */
98133779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
98233779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
98333779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
98433779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
985480cf27aSJed Brown 
986480cf27aSJed Brown   /*
987e8fb0fc0SBarry Smith      Build the options database
988e5c89e4eSSatish Balay   */
989c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
990e5c89e4eSSatish Balay 
991703f3eceSBarry Smith   /* call a second time so it can look in the options database */
992703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
9936dc8fec2Sbcordonn 
994e5c89e4eSSatish Balay   /*
995e5c89e4eSSatish Balay      Print main application help message
996e5c89e4eSSatish Balay   */
9972d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
998e5c89e4eSSatish Balay   if (help && flg) {
999e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
1000c0bb3764SVaclav Hapla     ierr = PetscPrintf(PETSC_COMM_WORLD,"----------------------------------------\n");CHKERRQ(ierr);
1001e5c89e4eSSatish Balay   }
1002e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
1003e5c89e4eSSatish Balay 
1004d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
1005d45a07a7SBarry Smith 
1006e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
100711525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1008f4202a44SBarry Smith #endif
1009f4202a44SBarry Smith 
1010e5c89e4eSSatish Balay   /*
1011e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1012e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1013e5c89e4eSSatish Balay   */
1014e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1015e5c89e4eSSatish Balay 
1016e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
101702c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1018e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
101902c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1020cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1021101491d0SBarry Smith   {
1022101491d0SBarry Smith     PetscBool omp_view_flag;
1023cd1b99a6SBarry Smith     char      *threads = getenv("OMP_NUM_THREADS");
1024e5c89e4eSSatish Balay 
1025cd1b99a6SBarry Smith    if (threads) {
102602c9f0b5SLisandro Dalcin      ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
1027f90b075cSBarry Smith      (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
1028cd1b99a6SBarry Smith    } else {
1029cd1b99a6SBarry Smith #define NMAX  10000
1030101491d0SBarry Smith      int          i;
1031cd1b99a6SBarry Smith       PetscScalar *x;
1032cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1033cd1b99a6SBarry Smith #pragma omp parallel for
1034cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1035cd1b99a6SBarry Smith         x[i] = 0.0;
1036f90b075cSBarry Smith         PetscNumOMPThreads  = (PetscInt) omp_get_num_threads();
1037cd1b99a6SBarry Smith       }
1038cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
103902c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1040101491d0SBarry Smith     }
1041101491d0SBarry Smith     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
1042f90b075cSBarry Smith     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);
1043101491d0SBarry Smith     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
1044101491d0SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1045101491d0SBarry Smith     if (flg) {
104602c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1047f90b075cSBarry Smith       omp_set_num_threads((int)PetscNumOMPThreads);
1048101491d0SBarry Smith     }
1049101491d0SBarry Smith     if (omp_view_flag) {
1050f90b075cSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
1051cd1b99a6SBarry Smith     }
1052cd1b99a6SBarry Smith   }
1053cd1b99a6SBarry Smith #endif
1054ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1055c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1056ef6c6fedSBoyana Norris 
1057951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1058e39fd77fSBarry Smith   /*
1059e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1060e39fd77fSBarry Smith 
1061e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1062e39fd77fSBarry Smith   */
106330815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
10640298fd71SBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
106530815ce0SLisandro Dalcin   }
1066951e3c8eSBarry Smith #endif
1067e39fd77fSBarry Smith 
106841c0b4b3SShri Abhyankar   /*
106941c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
107041c0b4b3SShri Abhyankar   */
1071ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1072e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1073e1167bb9SShri Abhyankar #endif
1074e1167bb9SShri Abhyankar 
10752d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10762d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10772d53ad75SBarry Smith #endif
10782d53ad75SBarry Smith 
10795e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
108087c3beb6SLisandro Dalcin   {
108187c3beb6SLisandro Dalcin     PetscViewer viewer;
108222e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
10835e71baefSBarry Smith     if (flg) {
10845e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
108587c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
108687c3beb6SLisandro Dalcin     }
10875e71baefSBarry Smith   }
10885e71baefSBarry Smith #endif
1089dff31646SBarry Smith 
109087c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
109187c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
109287c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
109387c3beb6SLisandro Dalcin 
1094a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1095a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
10967b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1097a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
10989fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1099a56f64adSBarry Smith #endif
110055d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
110155d657eeSBarry Smith #endif
1102a56f64adSBarry Smith 
1103301d30feSBarry Smith   /*
1104301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1105301d30feSBarry Smith   */
1106301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
11072db0e300SLisandro Dalcin 
11082db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
11092db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1110301d30feSBarry Smith   PetscFunctionReturn(0);
1111e5c89e4eSSatish Balay }
1112e5c89e4eSSatish Balay 
11134097062eSBarry Smith #if defined(PETSC_USE_LOG)
111495c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1115ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1116ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
111795c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11184097062eSBarry Smith #endif
1119e5c89e4eSSatish Balay 
1120008a6e76SBarry Smith /*
1121008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1122008a6e76SBarry Smith */
1123008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1124008a6e76SBarry Smith {
1125008a6e76SBarry Smith   PetscErrorCode ierr;
1126008a6e76SBarry Smith 
1127008a6e76SBarry Smith   PetscFunctionBegin;
1128008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1129008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1130008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1131008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1132008a6e76SBarry Smith #endif
1133008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1134008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1135008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1136008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1137008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1138008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1139008a6e76SBarry Smith #endif
1140008a6e76SBarry Smith 
1141008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1142008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1143008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1144008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1145008a6e76SBarry Smith #endif
1146008a6e76SBarry Smith #endif
1147008a6e76SBarry Smith 
1148008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1149008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1150008a6e76SBarry Smith #endif
1151008a6e76SBarry Smith 
1152008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
115340df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1154008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1155008a6e76SBarry Smith #endif
1156008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1157008a6e76SBarry Smith   PetscFunctionReturn(0);
1158008a6e76SBarry Smith }
1159008a6e76SBarry Smith 
1160e5c89e4eSSatish Balay /*@C
1161e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1162e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1163e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1164e5c89e4eSSatish Balay 
1165e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1166e5c89e4eSSatish Balay 
1167e5c89e4eSSatish Balay    Options Database Keys:
116826a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1169e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11707eb1d149SBarry 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
1171e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
117292f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1173e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
117492f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1175e5c89e4eSSatish Balay 
1176e5c89e4eSSatish Balay    Level: beginner
1177e5c89e4eSSatish Balay 
1178e5c89e4eSSatish Balay    Note:
1179e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1180e5c89e4eSSatish Balay 
118188c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1182e5c89e4eSSatish Balay @*/
11837087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1184e5c89e4eSSatish Balay {
1185e5c89e4eSSatish Balay   PetscErrorCode ierr;
11864bb5149bSJed Brown   PetscMPIInt    rank;
1187a8d2bbe5SBarry Smith   PetscInt       nopt;
11882bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1189dff31646SBarry Smith   PetscBool      flg;
119010463e74SBarry Smith #if defined(PETSC_USE_LOG)
119110463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
119210463e74SBarry Smith #endif
1193e5c89e4eSSatish Balay 
1194e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
11954b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
11963cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1197e5c89e4eSSatish Balay   }
11983cbbc5ffSBarry Smith   PetscFunctionBegin;
11990298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1200b022a5c1SBarry Smith 
12011f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1202a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
120322580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1204a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1205a56f64adSBarry Smith #endif
120655d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
120755d657eeSBarry Smith #endif
1208c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1209dff31646SBarry Smith   if (flg) {
12101f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
12111f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
12121f817a21SBarry Smith 
1213c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
12141f817a21SBarry Smith     if (filename[0]) {
12151f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
12161f817a21SBarry Smith     }
1217dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1218dff31646SBarry Smith     cits[0] = 0;
1219dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
12201f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
12211f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12221f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
12231f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12241f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1225dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1226dff31646SBarry Smith   }
1227dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1228dff31646SBarry Smith 
1229c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
123004102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
123104102261SBarry Smith   {
123204102261SBarry Smith     PetscInt nmax = 2;
123304102261SBarry Smith     char     **buffs;
123404102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1235c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
123604102261SBarry Smith     if (flg1) {
123704102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
123804102261SBarry Smith       if (nmax == 1) {
123904102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
124004102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
124104102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
124204102261SBarry Smith       }
124304102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
124404102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
124504102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
124604102261SBarry Smith     }
124704102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
124804102261SBarry Smith   }
1249763ec7b1SBarry Smith   {
1250763ec7b1SBarry Smith     PetscInt nmax = 2;
1251763ec7b1SBarry Smith     char     **buffs;
1252763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1253763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1254763ec7b1SBarry Smith     if (flg1) {
1255763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1256763ec7b1SBarry Smith       if (nmax == 1) {
1257763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1258763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1259763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1260763ec7b1SBarry Smith       }
1261763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1262763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1263763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1264763ec7b1SBarry Smith     }
1265763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1266763ec7b1SBarry Smith   }
126704102261SBarry Smith #endif
126867234432SDmitry Karpeev   /*
126967234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
127067234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
127167234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
127267234432SDmitry Karpeev   */
127367234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
127404102261SBarry Smith 
12752d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12762d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12772d53ad75SBarry Smith #endif
12782d53ad75SBarry Smith 
12792d53ad75SBarry Smith 
1280e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1281dff31646SBarry Smith   flg = PETSC_FALSE;
1282c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1283d5649816SBarry Smith   if (flg) {
1284e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1285d5649816SBarry Smith   }
1286d5649816SBarry Smith #endif
1287d5649816SBarry Smith 
1288681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1289681455b2SBarry Smith   flg1 = PETSC_FALSE;
1290c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1291681455b2SBarry Smith   if (flg1) {
1292681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1293681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1294681455b2SBarry Smith   }
1295681455b2SBarry Smith #endif
1296681455b2SBarry Smith 
129767584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1298c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1299e5c89e4eSSatish Balay   if (!flg2) {
130090d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1301c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1302e5c89e4eSSatish Balay   }
1303e5c89e4eSSatish Balay   if (flg2) {
13040841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1305e5c89e4eSSatish Balay   }
130667584ceeSBarry Smith #endif
1307e5c89e4eSSatish Balay 
1308e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
130990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1310c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1311e5c89e4eSSatish Balay   if (flg1) {
1312e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1313205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1314e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1315e5c89e4eSSatish Balay   }
1316e5c89e4eSSatish Balay #endif
1317e5c89e4eSSatish Balay 
1318e5c89e4eSSatish Balay 
1319e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1320e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1321e5c89e4eSSatish Balay   mname[0] = 0;
1322c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1323e5c89e4eSSatish Balay   if (flg1) {
1324e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1325e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1326e5c89e4eSSatish Balay   }
1327e5c89e4eSSatish Balay #endif
1328356e5837SBarry Smith #endif
1329a297a907SKarl Rupp 
1330dd710f27SBarry Smith   /*
1331dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1332dd710f27SBarry Smith   */
1333dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1334dd710f27SBarry Smith 
1335356e5837SBarry Smith #if defined(PETSC_USE_LOG)
133687c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1337f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
133887c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
133987c3beb6SLisandro Dalcin 
1340356e5837SBarry Smith   mname[0] = 0;
1341c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1342e5c89e4eSSatish Balay   if (flg1) {
134391eabc43SBarry Smith     PetscViewer viewer;
134420a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
134591eabc43SBarry Smith     if (mname[0]) {
134691eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
134791eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13486bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
134933f85c2fSBarry Smith     } else {
135033f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13519a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
135233f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13539a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
135433f85c2fSBarry Smith     }
1355e5c89e4eSSatish Balay   }
1356a297a907SKarl Rupp 
1357dd710f27SBarry Smith   /*
1358dd710f27SBarry Smith      Free any objects created by the last block of code.
1359dd710f27SBarry Smith   */
1360dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1361dd710f27SBarry Smith 
1362dd710f27SBarry Smith   mname[0] = 0;
1363c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1364c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
13657ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1366e5c89e4eSSatish Balay #endif
136710463e74SBarry Smith 
136833f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
136910463e74SBarry Smith 
137090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1371c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1372e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
137390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1374c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1375e5c89e4eSSatish Balay   if (flg1) {
1376e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1377e5c89e4eSSatish Balay   }
137890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
137990d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13808bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1381c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1382c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1383c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1384e4c476e2SSatish Balay 
1385e5c89e4eSSatish Balay   if (flg2) {
1386be56827dSJed Brown     PetscViewer viewer;
138702ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
138802ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1389c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1390be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1391e5c89e4eSSatish Balay   }
1392e5c89e4eSSatish Balay 
1393e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1394c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1395c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1396e5c89e4eSSatish Balay 
139733fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1398c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
1399*cf9c20a2SJed Brown   if (PetscDefined(USE_DEBUG) && !flg1) flg3 = PETSC_TRUE;
1400e5c89e4eSSatish Balay   if (flg3) {
14013de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1402be56827dSJed Brown       PetscViewer viewer;
140302ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
140402ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1405c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1406be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1407e5c89e4eSSatish Balay     }
14083de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
14093de2bfdfSBarry Smith     if (nopt) {
14103de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
14113de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
14123de2bfdfSBarry Smith       if (nopt == 1) {
1413e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1414e5c89e4eSSatish Balay       } else {
14157582186dSLisandro Dalcin         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1416e5c89e4eSSatish Balay       }
14173de2bfdfSBarry Smith     } else if (flg3 && flg1) {
14183de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1419df12ba86SBarry Smith     }
1420c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1421e5c89e4eSSatish Balay   }
1422e5c89e4eSSatish Balay 
1423e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1424d45a07a7SBarry Smith   if (!PetscGlobalRank) {
142587f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
142616ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1427d45a07a7SBarry Smith   }
1428ec957eceSBarry Smith #endif
1429ec957eceSBarry Smith 
14304097062eSBarry Smith #if defined(PETSC_USE_LOG)
143110463e74SBarry Smith   /*
1432dbc8283eSBarry Smith        List all objects the user may have forgot to free
14332eff7a51SBarry Smith   */
143405df10baSBarry Smith   if (PetscObjectsLog) {
1435c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1436a64a8e02SBarry Smith     if (flg1) {
1437a64a8e02SBarry Smith       MPI_Comm local_comm;
14387eb1d149SBarry Smith       char     string[64];
1439a64a8e02SBarry Smith 
1440c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1441a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1442a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14437eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1444a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1445a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14460a1571b3SBarry Smith     }
144705df10baSBarry Smith   }
14484097062eSBarry Smith #endif
14494097062eSBarry Smith 
14504097062eSBarry Smith #if defined(PETSC_USE_LOG)
1451dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1452dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1453a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14544097062eSBarry Smith #endif
14552eff7a51SBarry Smith 
145633f85c2fSBarry Smith   /*
145733f85c2fSBarry Smith      Destroy any packages that registered a finalize
145833f85c2fSBarry Smith   */
145933f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
146033f85c2fSBarry Smith 
1461101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1462fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1463101409b8SToby Isaac #endif
1464101409b8SToby Isaac 
146533f85c2fSBarry Smith   /*
146648dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
146748dd1dffSBarry Smith 
146837e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
146948dd1dffSBarry Smith   */
147037e93019SBarry Smith 
14714028d114SSatish Balay   if (petsc_history) {
1472f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
147302c9f0b5SLisandro Dalcin     petsc_history = NULL;
1474e5c89e4eSSatish Balay   }
14759de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1476e94e781bSJacob Faibussowitsch   ierr = PetscInfoDestroy();CHKERRQ(ierr);
1477e5c89e4eSSatish Balay 
147867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
147992f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1480e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
148192f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1482e5c89e4eSSatish Balay     FILE *fd;
1483ed9cf6e9SBarry Smith     int  err;
1484e5c89e4eSSatish Balay 
1485dc92acbaSJed Brown     flg2 = PETSC_FALSE;
148692f119d6SBarry Smith     flg3 = PETSC_FALSE;
1487*cf9c20a2SJed Brown     if (PetscDefined(USE_DEBUG)) {ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);}
148892f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
148992f119d6SBarry Smith     fname[0] = 0;
149092f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1491e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1492e5c89e4eSSatish Balay 
14932e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1494e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1495e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1496ed9cf6e9SBarry Smith       err  = fclose(fd);
1497e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
149892f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1499e5c89e4eSSatish Balay       MPI_Comm local_comm;
1500e5c89e4eSSatish Balay 
1501e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1502e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1503e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1504e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1505e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1506e5c89e4eSSatish Balay     }
1507e5c89e4eSSatish Balay     fname[0] = 0;
150892f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,250,&flg1);CHKERRQ(ierr);
1509e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1510e5c89e4eSSatish Balay 
151192f119d6SBarry Smith       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
151292f119d6SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
151392f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1514ed9cf6e9SBarry Smith       err  = fclose(fd);
1515e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
151692f119d6SBarry Smith     } else if (flg1) {
151792f119d6SBarry Smith       MPI_Comm local_comm;
151892f119d6SBarry Smith 
151992f119d6SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
152092f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
152192f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
152292f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
152392f119d6SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1524e5c89e4eSSatish Balay     }
1525e5c89e4eSSatish Balay   }
152667584ceeSBarry Smith #endif
152720e2c332SMatthew G. Knepley 
15285486ca60SMatthew G. Knepley   /*
15295486ca60SMatthew G. Knepley      Close any open dynamic libraries
15305486ca60SMatthew G. Knepley   */
15315486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15325486ca60SMatthew G. Knepley 
1533e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15344416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1535e5c89e4eSSatish Balay 
1536e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
153702c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1538e5c89e4eSSatish Balay 
1539008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1540e5c89e4eSSatish Balay 
1541dbc8283eSBarry Smith   /*
1542efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1543efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1544efb80d3cSBarry Smith 
1545efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1546efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1547dbc8283eSBarry Smith  */
1548b770b1f6SSatish Balay   {
1549dbc8283eSBarry Smith     PetscCommCounter *counter;
1550dbc8283eSBarry Smith     PetscMPIInt      flg;
1551dbc8283eSBarry Smith     MPI_Comm         icomm;
1552265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
155347435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1554dbc8283eSBarry Smith     if (flg) {
1555265f3f35SJed Brown       icomm = ucomm.comm;
155647435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1557dbc8283eSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1558dbc8283eSBarry Smith 
155947435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
156047435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1561efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1562dbc8283eSBarry Smith     }
156347435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1564dbc8283eSBarry Smith     if (flg) {
1565265f3f35SJed Brown       icomm = ucomm.comm;
156647435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1567dbc8283eSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");
1568dbc8283eSBarry Smith 
156947435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
157047435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1571efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1572dbc8283eSBarry Smith     }
1573b770b1f6SSatish Balay   }
1574dbc8283eSBarry Smith 
157547435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
157647435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
157747435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
15785f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1579480cf27aSJed Brown 
15805ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15815ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15825ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1583ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1584ef19f930SBarry Smith 
1585e5c89e4eSSatish Balay   if (PetscBeganMPI) {
158699608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
158799b1327fSBarry Smith     PetscMPIInt flag;
158899b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1589e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
159099608316SBarry Smith #endif
1591e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1592e5c89e4eSSatish Balay   }
1593e5c89e4eSSatish Balay /*
1594e5c89e4eSSatish Balay 
1595e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1596e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1597e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1598e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1599e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16000e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1601e5c89e4eSSatish Balay    memory was not freed.
1602e5c89e4eSSatish Balay 
1603e5c89e4eSSatish Balay */
16041d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1605a297a907SKarl Rupp 
1606e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1607e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
16083db9a53dSBarry Smith   PetscFunctionReturn(0);
1609e5c89e4eSSatish Balay }
1610e5c89e4eSSatish Balay 
161143db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16128cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
161343db4dbbSBarry Smith {
161443db4dbbSBarry Smith   if (*a == *b) return 1;
161543db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
161643db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
161743db4dbbSBarry Smith   return 0;
161843db4dbbSBarry Smith }
1619a70650f6SBarry Smith #endif
162043db4dbbSBarry Smith 
162143db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16228cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
162343db4dbbSBarry Smith {
162443db4dbbSBarry Smith   if (*a == *b) return 1;
162543db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
162643db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
162743db4dbbSBarry Smith   return 0;
162843db4dbbSBarry Smith }
162943db4dbbSBarry Smith #endif
1630