xref: /petsc/src/sys/objects/pinit.c (revision 9245e7496f19085098a102c852c40f8ba5bf1cee)
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 
336de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
34e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
356de5d289SStefano Zampini #if defined(PETSC_HAVE_MPI_INIT_THREAD)
366de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
376de5d289SStefano Zampini #else
386de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
396de5d289SStefano Zampini #endif
40e5c89e4eSSatish Balay 
41480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
42480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
43480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
445f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
45480cf27aSJed Brown 
46e5c89e4eSSatish Balay /*
47e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
48e5c89e4eSSatish Balay */
4902c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
5002c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
51e5c89e4eSSatish Balay 
52ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
53ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
540f8e0872SSatish Balay 
55a2f94806SJed Brown PetscInt PetscHotRegionDepth;
56a2f94806SJed Brown 
57b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
58b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
59b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
60b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
61b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
62b22622e2STadeu Manoel #endif
63b22622e2STadeu Manoel 
64e5c89e4eSSatish Balay /*
65945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
6672a42c3cSBarry Smith 
6772a42c3cSBarry Smith    Collective
6872a42c3cSBarry Smith 
6972a42c3cSBarry Smith    Level: advanced
7072a42c3cSBarry Smith 
7195452b02SPatrick Sanan     Notes:
72a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
730f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
74a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
750f11a792SBarry Smith 
76a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
771ea5a559SBarry Smith 
7872a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
790f11a792SBarry Smith */
80945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
8172a42c3cSBarry Smith {
8272a42c3cSBarry Smith   PetscErrorCode ierr;
8372a42c3cSBarry Smith   int            myargc   = argc;
8472a42c3cSBarry Smith   char           **myargs = args;
8572a42c3cSBarry Smith 
8672a42c3cSBarry Smith   PetscFunctionBegin;
87c52ac268SRichard Tran Mills   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr;
881ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
89df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
9072a42c3cSBarry Smith   PetscFunctionReturn(ierr);
9172a42c3cSBarry Smith }
9272a42c3cSBarry Smith 
93f0865b08SBarry Smith /*
94a56f64adSBarry Smith       Used by Julia interface to get communicator
95f0865b08SBarry Smith */
96945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
97f0865b08SBarry Smith {
98f0865b08SBarry Smith   PetscFunctionBegin;
99f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
100f0865b08SBarry Smith   PetscFunctionReturn(0);
101f0865b08SBarry Smith }
102f0865b08SBarry Smith 
103e5c89e4eSSatish Balay /*@C
104e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
105e5c89e4eSSatish Balay         the command line arguments.
106e5c89e4eSSatish Balay 
107e5c89e4eSSatish Balay    Collective
108e5c89e4eSSatish Balay 
109e5c89e4eSSatish Balay    Level: advanced
110e5c89e4eSSatish Balay 
111e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
112e5c89e4eSSatish Balay @*/
1137087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
114e5c89e4eSSatish Balay {
115e5c89e4eSSatish Balay   PetscErrorCode ierr;
116e5c89e4eSSatish Balay   int            argc   = 0;
11702c9f0b5SLisandro Dalcin   char           **args = NULL;
118e5c89e4eSSatish Balay 
119e5c89e4eSSatish Balay   PetscFunctionBegin;
1200298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
121e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
122e5c89e4eSSatish Balay }
123e5c89e4eSSatish Balay 
124e5c89e4eSSatish Balay /*@
125e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
126e5c89e4eSSatish Balay 
12793b6d2d1SJed Brown    Level: beginner
128e5c89e4eSSatish Balay 
129e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
130e5c89e4eSSatish Balay @*/
1317087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized)
132e5c89e4eSSatish Balay {
133e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
13493b6d2d1SJed Brown   return 0;
135e5c89e4eSSatish Balay }
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay /*@
138e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
139e5c89e4eSSatish Balay 
140e5c89e4eSSatish Balay    Level: developer
141e5c89e4eSSatish Balay 
142e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
143e5c89e4eSSatish Balay @*/
1447087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
145e5c89e4eSSatish Balay {
146e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
14793b6d2d1SJed Brown   return 0;
148e5c89e4eSSatish Balay }
149e5c89e4eSSatish Balay 
15057171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []);
151e5c89e4eSSatish Balay 
152e5c89e4eSSatish Balay /*
153e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
154e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
155e5c89e4eSSatish Balay */
156367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
157e5c89e4eSSatish Balay 
158367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
159e5c89e4eSSatish Balay {
160e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
161e5c89e4eSSatish Balay 
162e5c89e4eSSatish Balay   PetscFunctionBegin;
163e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
164e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
16541e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
166e5c89e4eSSatish Balay   }
167e5c89e4eSSatish Balay 
168e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
169e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
170e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
171e5c89e4eSSatish Balay   }
172812af9f3SBarry Smith   PetscFunctionReturnVoid();
173e5c89e4eSSatish Balay }
174e5c89e4eSSatish Balay 
175e5c89e4eSSatish Balay /*
176e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
177e5c89e4eSSatish Balay sum of the second entry.
178b693b147SBarry Smith 
17976f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
180367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
181b693b147SBarry Smith there would be no place to store the both needed results.
182e5c89e4eSSatish Balay */
18376ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
184e5c89e4eSSatish Balay {
185e5c89e4eSSatish Balay   PetscErrorCode ierr;
186e5c89e4eSSatish Balay 
187e5c89e4eSSatish Balay   PetscFunctionBegin;
188d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
189d6e4c47cSJed Brown   {
190d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
191367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
192d6e4c47cSJed Brown     *max = work.max;
193d6e4c47cSJed Brown     *sum = work.sum;
194d6e4c47cSJed Brown   }
195d6e4c47cSJed Brown #else
196d6e4c47cSJed Brown   {
197d6e4c47cSJed Brown     PetscMPIInt    size,rank;
198d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
199e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
200e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
201785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
202367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2036ac3741eSJed Brown     *max = work[rank].max;
2046ac3741eSJed Brown     *sum = work[rank].sum;
205e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
206d6e4c47cSJed Brown   }
207d6e4c47cSJed Brown #endif
208e5c89e4eSSatish Balay   PetscFunctionReturn(0);
209e5c89e4eSSatish Balay }
210e5c89e4eSSatish Balay 
211e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
212e5c89e4eSSatish Balay 
213570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
21406a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
215e5c89e4eSSatish Balay 
2168cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
217e5c89e4eSSatish Balay {
218e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
219e5c89e4eSSatish Balay 
220e5c89e4eSSatish Balay   PetscFunctionBegin;
2217c2de775SJed Brown   if (*datatype == MPIU_REAL) {
222e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
223a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2247c2de775SJed Brown   }
2257c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2267c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2277c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
228a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2297c2de775SJed Brown   }
2307c2de775SJed Brown #endif
2317c2de775SJed Brown   else {
2327c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
23341e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
234e2e03761SBarry Smith   }
235812af9f3SBarry Smith   PetscFunctionReturnVoid();
236e5c89e4eSSatish Balay }
237e5c89e4eSSatish Balay #endif
238e5c89e4eSSatish Balay 
239570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
240d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
241d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
242d9822059SBarry Smith 
2438cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
244d9822059SBarry Smith {
245d9822059SBarry Smith   PetscInt i,count = *cnt;
246d9822059SBarry Smith 
247d9822059SBarry Smith   PetscFunctionBegin;
2487c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2498c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
250a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2517c2de775SJed Brown   }
2527c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2537c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2547c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2557c2de775SJed Brown     for (i=0; i<count; i++) {
2567c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2577c2de775SJed Brown     }
2587c2de775SJed Brown   }
2597c2de775SJed Brown #endif
2607c2de775SJed Brown   else {
2617c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
26241e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2638c764dc5SJose Roman   }
264d9822059SBarry Smith   PetscFunctionReturnVoid();
265d9822059SBarry Smith }
266d9822059SBarry Smith 
2678cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
268d9822059SBarry Smith {
269d9822059SBarry Smith   PetscInt    i,count = *cnt;
270d9822059SBarry Smith 
271d9822059SBarry Smith   PetscFunctionBegin;
2727c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2738c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
274a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2757c2de775SJed Brown   }
2767c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2777c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2787c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2797c2de775SJed Brown     for (i=0; i<count; i++) {
2807c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2817c2de775SJed Brown     }
2827c2de775SJed Brown   }
2837c2de775SJed Brown #endif
2847c2de775SJed Brown   else {
2858c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
28641e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2878c764dc5SJose Roman   }
288d9822059SBarry Smith   PetscFunctionReturnVoid();
289d9822059SBarry Smith }
290d9822059SBarry Smith #endif
291d9822059SBarry Smith 
292480cf27aSJed Brown /*
293480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
294480cf27aSJed Brown 
295ff0e51ddSBarry 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.
296480cf27aSJed Brown 
29712801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
298480cf27aSJed Brown 
299480cf27aSJed Brown */
30033779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
301480cf27aSJed Brown {
302480cf27aSJed Brown   PetscErrorCode   ierr;
30305342407SJunchao Zhang   PetscCommCounter *counter=(PetscCommCounter*)count_val;
304480cf27aSJed Brown 
305480cf27aSJed Brown   PetscFunctionBegin;
30602c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
30705342407SJunchao Zhang   ierr = PetscFree(counter->iflags);CHKERRMPI(ierr);
30805342407SJunchao Zhang   ierr = PetscFree(counter);CHKERRMPI(ierr);
309480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
310480cf27aSJed Brown }
311480cf27aSJed Brown 
312480cf27aSJed Brown /*
31347435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
314da3039f7SJed Brown   calls MPI_Comm_free().
315da3039f7SJed Brown 
316da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
317480cf27aSJed Brown 
318ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
319480cf27aSJed Brown 
32012801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
321480cf27aSJed Brown 
322480cf27aSJed Brown */
32333779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
324480cf27aSJed Brown {
325480cf27aSJed Brown   PetscErrorCode                    ierr;
32633779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
327480cf27aSJed Brown 
328480cf27aSJed Brown   PetscFunctionBegin;
32912801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
330ec4fadc2SJed Brown   icomm.ptr = attr_val;
33176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
33233779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
33333779a13SJunchao Zhang     PetscMPIInt flg;
33433779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
33547435625SJed Brown     ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
33633779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
33733779a13SJunchao 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");
33833779a13SJunchao Zhang   }
33947435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
34033779a13SJunchao Zhang   ierr = PetscInfo2(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);CHKERRMPI(ierr);
341da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
342b89831e5SBarry Smith }
343da3039f7SJed Brown 
344da3039f7SJed Brown /*
34533779a13SJunchao 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.
346da3039f7SJed Brown  */
34733779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
348da3039f7SJed Brown {
349da3039f7SJed Brown   PetscErrorCode ierr;
350da3039f7SJed Brown 
351da3039f7SJed Brown   PetscFunctionBegin;
35202c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
353480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
354480cf27aSJed Brown }
355480cf27aSJed Brown 
35633779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
35742218b76SBarry Smith 
358951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3598cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3608cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3618cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
362e39fd77fSBarry Smith #endif
363e39fd77fSBarry Smith 
3640f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
36512801b39SBarry Smith 
366eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
367eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3686ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
36902c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
370dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
371e5c89e4eSSatish Balay 
372dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
373051e4cf2SJed Brown {
374051e4cf2SJed Brown   PetscErrorCode ierr;
375051e4cf2SJed Brown 
376051e4cf2SJed Brown   PetscFunctionBegin;
377051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
378a1601671SBarry 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);
379051e4cf2SJed 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);
380051e4cf2SJed Brown   PetscFunctionReturn(0);
381051e4cf2SJed Brown }
382e5c89e4eSSatish Balay 
3832d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3842d747510SLisandro Dalcin 
3852d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
3862d747510SLisandro Dalcin {
3872d747510SLisandro Dalcin   PetscErrorCode ierr;
3882d747510SLisandro Dalcin 
3892d747510SLisandro Dalcin   PetscFunctionBegin;
390589a23caSBarry Smith   ierr  = PetscStrncpy(programname,name,sizeof(programname));CHKERRQ(ierr);
3912d747510SLisandro Dalcin   PetscFunctionReturn(0);
3922d747510SLisandro Dalcin }
3932d747510SLisandro Dalcin 
3942d747510SLisandro Dalcin /*@C
3952d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
3962d747510SLisandro Dalcin 
3972d747510SLisandro Dalcin     Not Collective
3982d747510SLisandro Dalcin 
3992d747510SLisandro Dalcin     Input Parameter:
4002d747510SLisandro Dalcin .   len - length of the string name
4012d747510SLisandro Dalcin 
4022d747510SLisandro Dalcin     Output Parameter:
4032d747510SLisandro Dalcin .   name - the name of the running program
4042d747510SLisandro Dalcin 
4052d747510SLisandro Dalcin    Level: advanced
4062d747510SLisandro Dalcin 
4072d747510SLisandro Dalcin     Notes:
4082d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4092d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4102d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4112d747510SLisandro Dalcin @*/
4122d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4132d747510SLisandro Dalcin {
4142d747510SLisandro Dalcin   PetscErrorCode ierr;
4152d747510SLisandro Dalcin 
4162d747510SLisandro Dalcin   PetscFunctionBegin;
4172d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4182d747510SLisandro Dalcin   PetscFunctionReturn(0);
4192d747510SLisandro Dalcin }
4202d747510SLisandro Dalcin 
421e5c89e4eSSatish Balay /*@C
422e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
423e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
424e5c89e4eSSatish Balay 
425e5c89e4eSSatish Balay    Not Collective
426e5c89e4eSSatish Balay 
427e5c89e4eSSatish Balay    Output Parameters:
428e5c89e4eSSatish Balay +  argc - count of number of command line arguments
429e5c89e4eSSatish Balay -  args - the command line arguments
430e5c89e4eSSatish Balay 
431e5c89e4eSSatish Balay    Level: intermediate
432e5c89e4eSSatish Balay 
433e5c89e4eSSatish Balay    Notes:
434e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
435e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
436e5c89e4eSSatish Balay 
437f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
438f177e3b1SBarry Smith 
439793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
440e5c89e4eSSatish Balay 
441e5c89e4eSSatish Balay @*/
4427087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
443e5c89e4eSSatish Balay {
444e5c89e4eSSatish Balay   PetscFunctionBegin;
44517186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
446e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
447e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
448e5c89e4eSSatish Balay   PetscFunctionReturn(0);
449e5c89e4eSSatish Balay }
450e5c89e4eSSatish Balay 
451793721a6SBarry Smith /*@C
452793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
453793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
454793721a6SBarry Smith 
455793721a6SBarry Smith    Not Collective
456793721a6SBarry Smith 
457793721a6SBarry Smith    Output Parameters:
458793721a6SBarry Smith .  args - the command line arguments
459793721a6SBarry Smith 
460793721a6SBarry Smith    Level: intermediate
461793721a6SBarry Smith 
462793721a6SBarry Smith    Notes:
463793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
464793721a6SBarry Smith 
465793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
466793721a6SBarry Smith 
467793721a6SBarry Smith @*/
4687087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
469793721a6SBarry Smith {
470793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
471793721a6SBarry Smith   PetscErrorCode ierr;
472793721a6SBarry Smith 
473793721a6SBarry Smith   PetscFunctionBegin;
47417186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4752d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
476785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
477793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
478793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
479793721a6SBarry Smith   }
4802d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
481793721a6SBarry Smith   PetscFunctionReturn(0);
482793721a6SBarry Smith }
483793721a6SBarry Smith 
484793721a6SBarry Smith /*@C
485793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
486793721a6SBarry Smith 
487793721a6SBarry Smith    Not Collective
488793721a6SBarry Smith 
489793721a6SBarry Smith    Output Parameters:
490793721a6SBarry Smith .  args - the command line arguments
491793721a6SBarry Smith 
492793721a6SBarry Smith    Level: intermediate
493793721a6SBarry Smith 
494793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
495793721a6SBarry Smith 
496793721a6SBarry Smith @*/
4977087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
498793721a6SBarry Smith {
499793721a6SBarry Smith   PetscInt       i = 0;
500793721a6SBarry Smith   PetscErrorCode ierr;
501793721a6SBarry Smith 
502793721a6SBarry Smith   PetscFunctionBegin;
503a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
504793721a6SBarry Smith   while (args[i]) {
505793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
506793721a6SBarry Smith     i++;
507793721a6SBarry Smith   }
508793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
509793721a6SBarry Smith   PetscFunctionReturn(0);
510793721a6SBarry Smith }
511793721a6SBarry Smith 
51211525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
51330befbd2SBarry Smith #include <petscconfiginfo.h>
51430befbd2SBarry Smith 
51595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
51611525c0dSBarry Smith {
51711525c0dSBarry Smith   if (!PetscGlobalRank) {
51830befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
51911525c0dSBarry Smith     int            port;
520ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
52111525c0dSBarry Smith     size_t         applinelen,introlen;
52211525c0dSBarry Smith     PetscErrorCode ierr;
523ffbd1cfbSBarry Smith     char           sawsurl[256];
52411525c0dSBarry Smith 
525c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
52611525c0dSBarry Smith     if (flg) {
52711525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
52811525c0dSBarry Smith 
529589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL);CHKERRQ(ierr);
53011525c0dSBarry Smith       if (sawslog[0]) {
53111525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
53211525c0dSBarry Smith       } else {
53311525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
53411525c0dSBarry Smith       }
53511525c0dSBarry Smith     }
536589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg);CHKERRQ(ierr);
53711525c0dSBarry Smith     if (flg) {
53811525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
53911525c0dSBarry Smith     }
540c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
541ffbd1cfbSBarry Smith     if (selectport) {
542ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
543ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
544ffbd1cfbSBarry Smith     } else {
545c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
54611525c0dSBarry Smith       if (flg) {
54711525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
54811525c0dSBarry Smith       }
549ffbd1cfbSBarry Smith     }
550589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg);CHKERRQ(ierr);
55111525c0dSBarry Smith     if (flg) {
55211525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
55311525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5549c1e0ce8SBarry Smith     } else {
555c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5569c1e0ce8SBarry Smith       if (flg) {
557589a23caSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root));CHKERRQ(ierr);
5589c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5599c1e0ce8SBarry Smith       }
56011525c0dSBarry Smith     }
561c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
56211525c0dSBarry Smith     if (flg2) {
56311525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
56411525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
565589a23caSBarry Smith       ierr = PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root);CHKERRQ(ierr);
56611525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
56711525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
56843da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
56911525c0dSBarry Smith     }
570589a23caSBarry Smith     ierr = PetscGetProgramName(programname,sizeof(programname));CHKERRQ(ierr);
57111525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
57211525c0dSBarry Smith     introlen   = 4096 + applinelen;
57330a8c9c0SSurtai Han     applinelen += 1024;
57411525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
57511525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
57611525c0dSBarry Smith 
57711525c0dSBarry Smith     if (rootlocal) {
57811525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
57911525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
58011525c0dSBarry Smith     }
58176a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
58211525c0dSBarry Smith     if (rootlocal && help) {
583928bb9adSStefano 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);
58411525c0dSBarry Smith     } else if (help) {
585928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
58611525c0dSBarry Smith     } else {
587928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
58811525c0dSBarry Smith     }
589b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
59030befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
59111525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
592a8d69d7bSBarry 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"
593df62c222SBarry 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"
594928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
59543da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
59611525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
59711525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
598ffbd1cfbSBarry Smith     if (selectport) {
599aa573868SBarry Smith       PetscBool silent;
6007d812c46SBarry Smith 
6017d812c46SBarry Smith       ierr = SAWs_Initialize();
6027d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6037d812c46SBarry Smith       while (ierr) {
6047d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6057d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6067d812c46SBarry Smith         ierr = SAWs_Initialize();
6077d812c46SBarry Smith       }
6087d812c46SBarry Smith 
609aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
610aa573868SBarry Smith       if (!silent) {
611ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
612ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
613ffbd1cfbSBarry Smith       }
6147d812c46SBarry Smith     } else {
6157d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
616aa573868SBarry Smith     }
6170af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6180af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6190af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6200af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6210af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
62211525c0dSBarry Smith   }
62311525c0dSBarry Smith   PetscFunctionReturn(0);
62411525c0dSBarry Smith }
62511525c0dSBarry Smith #endif
62611525c0dSBarry Smith 
6274dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6284dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6294dfee713SSatish Balay {
6304dfee713SSatish Balay   PetscFunctionBegin;
6314dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6324dfee713SSatish Balay     /* see MPI.py for details on this bug */
6334dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6344dfee713SSatish Balay #endif
6354dfee713SSatish Balay   PetscFunctionReturn(0);
6364dfee713SSatish Balay }
6374dfee713SSatish Balay 
638a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
639a56f64adSBarry Smith #include <adios.h>
64022580e64SBarry Smith #include <adios_read.h>
6417b56e58cSSatish Balay int64_t Petsc_adios_group;
642a56f64adSBarry Smith #endif
64355d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
64455d657eeSBarry Smith #include <adios2_c.h>
64555d657eeSBarry Smith #endif
646cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
647cd1b99a6SBarry Smith #include <omp.h>
648f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
649cd1b99a6SBarry Smith #endif
650a56f64adSBarry Smith 
651bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H)
652bc8a24ecSBarry Smith #include <dlfcn.h>
653bc8a24ecSBarry Smith #endif
654bc8a24ecSBarry Smith 
655e5c89e4eSSatish Balay /*@C
656e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
657e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
658e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
659e5c89e4eSSatish Balay    your program -- usually the very first line!
660e5c89e4eSSatish Balay 
661e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
662e5c89e4eSSatish Balay 
663e5c89e4eSSatish Balay    Input Parameters:
664e5c89e4eSSatish Balay +  argc - count of number of command line arguments
665e5c89e4eSSatish Balay .  args - the command line arguments
666c5b5d8d5SVaclav Hapla .  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
667c5b5d8d5SVaclav Hapla           Use NULL to not check for code specific file.
668c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
6690298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
670e5c89e4eSSatish Balay 
67105827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
67205827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
67305827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
67405827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
67505827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
676e5c89e4eSSatish Balay 
677e5c89e4eSSatish Balay    Options Database Keys:
6787ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6797ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
680e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
6818a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
6828a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
6838a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
6848a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
6858a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
686e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
687e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
688e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
689e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
69079dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
69179dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
69279dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
693aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
69492f119d6SBarry 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
69592f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
69692f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
69792f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
698e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
699e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
700e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
701e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
702e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7030841954dSBarry Smith -  -memory_view - Print memory usage at end of run
704e5c89e4eSSatish Balay 
705c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
706c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
707c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
708c5b5d8d5SVaclav Hapla -  -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet()
709c5b5d8d5SVaclav Hapla 
710c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
711c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
712c5b5d8d5SVaclav Hapla    They can be used also in option files.
713c5b5d8d5SVaclav Hapla 
714c5b5d8d5SVaclav Hapla    See PetscOptionsMonitorSet() to do monitoring programmatically.
715c5b5d8d5SVaclav Hapla 
716e5c89e4eSSatish Balay    Options Database Keys for Profiling:
717a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
718fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
719217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
720217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
721495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
722e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7239a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
72479dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
7259a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
726495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
727f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
728495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
729495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
730c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
73187c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
732c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
733495fc317SBarry Smith 
734609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
735e5c89e4eSSatish Balay 
736ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
737ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
738ffbd1cfbSBarry 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
739ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
740ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
741ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
742ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
743ffbd1cfbSBarry Smith 
744e5c89e4eSSatish Balay    Environmental Variables:
745e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
746e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
747e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
7484a971ea4SToby Isaac .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
7494a971ea4SToby Isaac .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
750e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
751e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
752e5c89e4eSSatish Balay 
753e5c89e4eSSatish Balay 
754e5c89e4eSSatish Balay    Level: beginner
755e5c89e4eSSatish Balay 
756e5c89e4eSSatish Balay    Notes:
757e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
758e5c89e4eSSatish Balay    it before PetscInitialize().
759e5c89e4eSSatish Balay 
760e5c89e4eSSatish Balay    Fortran Version:
761e5c89e4eSSatish Balay    In Fortran this routine has the format
762e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
763e5c89e4eSSatish Balay 
764e5c89e4eSSatish Balay +   ierr - error return code
765c5b5d8d5SVaclav Hapla -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
766c5b5d8d5SVaclav Hapla           Use PETSC_NULL_CHARACTER to not check for code specific file.
767c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
768e5c89e4eSSatish Balay 
769e5c89e4eSSatish Balay    Important Fortran Note:
7700eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7710298fd71SBarry Smith    null character string; you CANNOT just use NULL as
772a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
773e5c89e4eSSatish Balay 
77401cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
77501cb0274SBarry Smith    calling PetscInitialize().
776e5c89e4eSSatish Balay 
77701cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
778e5c89e4eSSatish Balay 
779e5c89e4eSSatish Balay @*/
7807087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
781e5c89e4eSSatish Balay {
782e5c89e4eSSatish Balay   PetscErrorCode ierr;
7834bb5149bSJed Brown   PetscMPIInt    flag, size;
78419c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
785e5c89e4eSSatish Balay   char           hostname[256];
786e5c89e4eSSatish Balay 
787e5c89e4eSSatish Balay   PetscFunctionBegin;
788e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7893d96e996SBarry Smith   /*
7903d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7913d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7923d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7933d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7943d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7953d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7963d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7973d96e996SBarry Smith       listed above and since that time are compatible.
7983d96e996SBarry Smith 
7993d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
8003d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
8013d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
8023d96e996SBarry Smith       to perform the checking.
8033d96e996SBarry Smith 
8043d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
8053d96e996SBarry Smith 
8063d96e996SBarry Smith       Questions:
8073d96e996SBarry Smith 
8083d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
8093d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
8103d96e996SBarry Smith   */
8113d96e996SBarry Smith 
81219c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
81319c5658aSBarry Smith   {
81419c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
81519c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
81619c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8173d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
81819c5658aSBarry Smith #if defined(MPICH_VERSION)
8193d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
82019c5658aSBarry Smith     {
82119c5658aSBarry Smith       char *ver,*lf;
82219c5658aSBarry Smith       flg = PETSC_FALSE;
82319c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
82419c5658aSBarry Smith       if (ver) {
82519c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
82619c5658aSBarry Smith         if (lf) {
82719c5658aSBarry Smith           *lf = 0;
82819c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
82919c5658aSBarry Smith         }
83019c5658aSBarry Smith       }
83119c5658aSBarry Smith       if (!flg) {
83219c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8333d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
83419c5658aSBarry Smith       }
83519c5658aSBarry Smith     }
8363d96e996SBarry Smith #endif
8373d96e996SBarry 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?) */
83819c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
83919c5658aSBarry Smith     {
84019c5658aSBarry Smith       char *ver,bs[32],*bsf;
84119c5658aSBarry Smith       flg = PETSC_FALSE;
84219c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
84319c5658aSBarry Smith       if (ver) {
8442e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
84519c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
84619c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
84719c5658aSBarry Smith       }
84819c5658aSBarry Smith       if (!flg) {
84919c5658aSBarry 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);
8503d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
85119c5658aSBarry Smith       }
85219c5658aSBarry Smith     }
85319c5658aSBarry Smith #endif
85419c5658aSBarry Smith   }
85519c5658aSBarry Smith #endif
85619c5658aSBarry Smith 
857bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM)
858bc8a24ecSBarry Smith   {
859bc8a24ecSBarry Smith     PetscInt cnt = 0;
860bc8a24ecSBarry 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 */
861bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++;
862bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++;
863bc8a24ecSBarry Smith     if (cnt > 1) {
864bc8a24ecSBarry Smith       fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
865bc8a24ecSBarry Smith       return PETSC_ERR_MPI_LIB_INCOMP;
866bc8a24ecSBarry Smith     }
867bc8a24ecSBarry Smith   }
868bc8a24ecSBarry Smith #endif
869bc8a24ecSBarry Smith 
870ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
871d89683f4Sbcordonn   PETSC_STDOUT = stdout;
872ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
873e5c89e4eSSatish Balay 
8748ad20175SVaclav Hapla   /* CHKERRQ can be used from now */
8758ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_TRUE;
8768ad20175SVaclav Hapla 
8770c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8780c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8790c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8800c30907bSSatish Balay #endif
8810c30907bSSatish Balay 
8824416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
883e5c89e4eSSatish Balay 
884e5c89e4eSSatish Balay   /*
885e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
886e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
887e5c89e4eSSatish Balay   */
888e5c89e4eSSatish Balay   if (argc && *argc) {
889e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
890e5c89e4eSSatish Balay   } else {
891e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
892e5c89e4eSSatish Balay   }
893e5c89e4eSSatish Balay 
894e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
895e5c89e4eSSatish Balay   if (!flag) {
896e32f2f54SBarry 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");
8974dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
8985e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8995e765c61SJed Brown     {
9005e765c61SJed Brown       PetscMPIInt provided;
9016de5d289SStefano Zampini       ierr = MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided);CHKERRQ(ierr);
9025e765c61SJed Brown     }
9035e765c61SJed Brown #else
904e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
9055e765c61SJed Brown #endif
906e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
907e5c89e4eSSatish Balay   }
9084dfee713SSatish Balay 
909e5c89e4eSSatish Balay   if (argc && args) {
910e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
911e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
912e5c89e4eSSatish Balay   }
913e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
9145ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
9155ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
9165ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
917ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
918e5c89e4eSSatish Balay 
919a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
920d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
921e5c89e4eSSatish Balay 
9220f9be574SPeter Hill   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
92312801b39SBarry Smith     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
92412801b39SBarry Smith     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
9250f9be574SPeter Hill   }
92612801b39SBarry Smith 
927e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
928e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
929e5c89e4eSSatish Balay 
930e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
931e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
932e5c89e4eSSatish Balay 
9338ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9348ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9357cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
936e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
937e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
938e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
939e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
940e316c87fSJed Brown #endif
941e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9428ad47952SJed Brown 
943e5c89e4eSSatish Balay   /*
944e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
945e5c89e4eSSatish Balay      shared libraries the constructors for global variables
946e5c89e4eSSatish Balay      are not called; at least on IRIX.
947e5c89e4eSSatish Balay   */
948886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
949e5c89e4eSSatish Balay   {
950252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
95150f81f78SJed Brown     PetscComplex ic(0.0,1.0);
952e5c89e4eSSatish Balay     PETSC_i = ic;
953252ecd31SSatish Balay #else
95450f81f78SJed Brown     PETSC_i = _Complex_I;
955b7940d39SSatish Balay #endif
956762437b8SSatish Balay   }
957762437b8SSatish Balay 
9582c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
959e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
960500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
961500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
962500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9632c876bd9SBarry Smith #endif
964886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
965e5c89e4eSSatish Balay 
966e5c89e4eSSatish Balay   /*
967e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
968e5c89e4eSSatish Balay      half of the entries and maxes the second half.
969e5c89e4eSSatish Balay   */
970367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
971e5c89e4eSSatish Balay 
972ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
973c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
974c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9757c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9768c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9778c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9788c764dc5SJose Roman #endif
979d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
980d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
981570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
982570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
983570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
984570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
985570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
986c90a1750SBarry Smith #endif
987c90a1750SBarry Smith 
988570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
989cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
990cca4cb22SSatish Balay #endif
991cca4cb22SSatish Balay 
992e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
993e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
994e5c89e4eSSatish Balay 
99540df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
996e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
997e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
99844041f26SJed Brown #endif
999e5c89e4eSSatish Balay 
1000e5c89e4eSSatish Balay   /*
1001480cf27aSJed Brown      Attributes to be set on PETSc communicators
1002480cf27aSJed Brown   */
100333779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
100433779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
100533779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
100633779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
1007480cf27aSJed Brown 
1008480cf27aSJed Brown   /*
1009e8fb0fc0SBarry Smith      Build the options database
1010e5c89e4eSSatish Balay   */
1011c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
1012e5c89e4eSSatish Balay 
1013703f3eceSBarry Smith   /* call a second time so it can look in the options database */
1014703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
10156dc8fec2Sbcordonn 
1016e5c89e4eSSatish Balay   /*
101757171095SVaclav Hapla      Check system options and print help
1018e5c89e4eSSatish Balay   */
101957171095SVaclav Hapla   ierr = PetscOptionsCheckInitial_Private(help);CHKERRQ(ierr);
1020e5c89e4eSSatish Balay 
1021d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
1022d45a07a7SBarry Smith 
1023e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
102411525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1025f4202a44SBarry Smith #endif
1026f4202a44SBarry Smith 
1027e5c89e4eSSatish Balay   /*
1028e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1029e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1030e5c89e4eSSatish Balay   */
1031e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1032e5c89e4eSSatish Balay 
1033e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
103402c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1035e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
103602c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1037cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1038101491d0SBarry Smith   {
1039101491d0SBarry Smith     PetscBool omp_view_flag;
1040cd1b99a6SBarry Smith     char      *threads = getenv("OMP_NUM_THREADS");
1041e5c89e4eSSatish Balay 
1042cd1b99a6SBarry Smith    if (threads) {
104302c9f0b5SLisandro Dalcin      ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
1044f90b075cSBarry Smith      (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
1045cd1b99a6SBarry Smith    } else {
1046cd1b99a6SBarry Smith #define NMAX  10000
1047101491d0SBarry Smith      int          i;
1048cd1b99a6SBarry Smith       PetscScalar *x;
1049cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1050cd1b99a6SBarry Smith #pragma omp parallel for
1051cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1052cd1b99a6SBarry Smith         x[i] = 0.0;
1053f90b075cSBarry Smith         PetscNumOMPThreads  = (PetscInt) omp_get_num_threads();
1054cd1b99a6SBarry Smith       }
1055cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
105602c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1057101491d0SBarry Smith     }
1058101491d0SBarry Smith     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
1059f90b075cSBarry 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);
1060101491d0SBarry Smith     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
1061101491d0SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1062101491d0SBarry Smith     if (flg) {
106302c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1064f90b075cSBarry Smith       omp_set_num_threads((int)PetscNumOMPThreads);
1065101491d0SBarry Smith     }
1066101491d0SBarry Smith     if (omp_view_flag) {
1067f90b075cSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
1068cd1b99a6SBarry Smith     }
1069cd1b99a6SBarry Smith   }
1070cd1b99a6SBarry Smith #endif
1071ef6c6fedSBoyana Norris 
1072951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1073e39fd77fSBarry Smith   /*
1074e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1075e39fd77fSBarry Smith 
1076e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1077e39fd77fSBarry Smith   */
107830815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
10790298fd71SBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
108030815ce0SLisandro Dalcin   }
1081951e3c8eSBarry Smith #endif
1082e39fd77fSBarry Smith 
108341c0b4b3SShri Abhyankar   /*
108441c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
108541c0b4b3SShri Abhyankar   */
1086ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1087e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1088e1167bb9SShri Abhyankar #endif
1089e1167bb9SShri Abhyankar 
10902d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10912d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10922d53ad75SBarry Smith #endif
10932d53ad75SBarry Smith 
10945e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
109587c3beb6SLisandro Dalcin   {
109687c3beb6SLisandro Dalcin     PetscViewer viewer;
109722e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
10985e71baefSBarry Smith     if (flg) {
10995e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
110087c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
110187c3beb6SLisandro Dalcin     }
11025e71baefSBarry Smith   }
11035e71baefSBarry Smith #endif
1104dff31646SBarry Smith 
110587c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
110687c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
110787c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
110887c3beb6SLisandro Dalcin 
1109a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1110a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
11117b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1112a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
11139fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1114a56f64adSBarry Smith #endif
111555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
111655d657eeSBarry Smith #endif
1117a56f64adSBarry Smith 
1118301d30feSBarry Smith   /*
111957171095SVaclav Hapla       Set flag that we are completely initialized
1120301d30feSBarry Smith   */
1121301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
11222db0e300SLisandro Dalcin 
11232db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
11242db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1125301d30feSBarry Smith   PetscFunctionReturn(0);
1126e5c89e4eSSatish Balay }
1127e5c89e4eSSatish Balay 
11284097062eSBarry Smith #if defined(PETSC_USE_LOG)
112995c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1130ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1131ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
113295c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11334097062eSBarry Smith #endif
1134e5c89e4eSSatish Balay 
1135008a6e76SBarry Smith /*
1136008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1137008a6e76SBarry Smith */
1138008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1139008a6e76SBarry Smith {
1140008a6e76SBarry Smith   PetscErrorCode ierr;
1141008a6e76SBarry Smith 
1142008a6e76SBarry Smith   PetscFunctionBegin;
1143008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1144008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1145008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1146008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1147008a6e76SBarry Smith #endif
1148008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1149008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1150008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1151008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1152008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1153008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1154008a6e76SBarry Smith #endif
1155008a6e76SBarry Smith 
1156008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1157008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1158008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1159008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1160008a6e76SBarry Smith #endif
1161008a6e76SBarry Smith #endif
1162008a6e76SBarry Smith 
1163008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1164008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1165008a6e76SBarry Smith #endif
1166008a6e76SBarry Smith 
1167008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
116840df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1169008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1170008a6e76SBarry Smith #endif
1171008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1172008a6e76SBarry Smith   PetscFunctionReturn(0);
1173008a6e76SBarry Smith }
1174008a6e76SBarry Smith 
1175e5c89e4eSSatish Balay /*@C
1176e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1177e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1178e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1179e5c89e4eSSatish Balay 
1180e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1181e5c89e4eSSatish Balay 
1182e5c89e4eSSatish Balay    Options Database Keys:
118326a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1184e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11857eb1d149SBarry 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
1186e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
118792f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1188e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
118992f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1190e5c89e4eSSatish Balay 
1191e5c89e4eSSatish Balay    Level: beginner
1192e5c89e4eSSatish Balay 
1193e5c89e4eSSatish Balay    Note:
1194e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1195e5c89e4eSSatish Balay 
119688c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1197e5c89e4eSSatish Balay @*/
11987087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1199e5c89e4eSSatish Balay {
1200e5c89e4eSSatish Balay   PetscErrorCode ierr;
12014bb5149bSJed Brown   PetscMPIInt    rank;
1202a8d2bbe5SBarry Smith   PetscInt       nopt;
12032bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1204dff31646SBarry Smith   PetscBool      flg;
120510463e74SBarry Smith #if defined(PETSC_USE_LOG)
120610463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
120710463e74SBarry Smith #endif
1208e5c89e4eSSatish Balay 
1209e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
12104b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
12113cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1212e5c89e4eSSatish Balay   }
12133cbbc5ffSBarry Smith   PetscFunctionBegin;
12140298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1215b022a5c1SBarry Smith 
12161f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1217a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
121822580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1219a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1220a56f64adSBarry Smith #endif
122155d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
122255d657eeSBarry Smith #endif
1223c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1224dff31646SBarry Smith   if (flg) {
12251f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
12261f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
12271f817a21SBarry Smith 
1228589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL);CHKERRQ(ierr);
12291f817a21SBarry Smith     if (filename[0]) {
12301f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
12311f817a21SBarry Smith     }
1232dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1233dff31646SBarry Smith     cits[0] = 0;
1234dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
12351f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
12361f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12371f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
12381f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12391f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1240dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1241dff31646SBarry Smith   }
1242dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1243dff31646SBarry Smith 
1244c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
124504102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
124604102261SBarry Smith   {
124704102261SBarry Smith     PetscInt nmax = 2;
124804102261SBarry Smith     char     **buffs;
124904102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1250c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
125104102261SBarry Smith     if (flg1) {
125204102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
125304102261SBarry Smith       if (nmax == 1) {
125404102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
125504102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
125604102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
125704102261SBarry Smith       }
125804102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
125904102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
126004102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
126104102261SBarry Smith     }
126204102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
126304102261SBarry Smith   }
1264763ec7b1SBarry Smith   {
1265763ec7b1SBarry Smith     PetscInt nmax = 2;
1266763ec7b1SBarry Smith     char     **buffs;
1267763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1268763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1269763ec7b1SBarry Smith     if (flg1) {
1270763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1271763ec7b1SBarry Smith       if (nmax == 1) {
1272763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1273763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1274763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1275763ec7b1SBarry Smith       }
1276763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1277763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1278763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1279763ec7b1SBarry Smith     }
1280763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1281763ec7b1SBarry Smith   }
128204102261SBarry Smith #endif
128304102261SBarry Smith 
12842d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12852d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12862d53ad75SBarry Smith #endif
12872d53ad75SBarry Smith 
1288e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1289dff31646SBarry Smith   flg = PETSC_FALSE;
1290c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1291d5649816SBarry Smith   if (flg) {
1292e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1293d5649816SBarry Smith   }
1294d5649816SBarry Smith #endif
1295d5649816SBarry Smith 
1296681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1297681455b2SBarry Smith   flg1 = PETSC_FALSE;
1298c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1299681455b2SBarry Smith   if (flg1) {
1300681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1301681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1302681455b2SBarry Smith   }
1303681455b2SBarry Smith #endif
1304681455b2SBarry Smith 
130567584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1306c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1307e5c89e4eSSatish Balay   if (!flg2) {
130890d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1309c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1310e5c89e4eSSatish Balay   }
1311e5c89e4eSSatish Balay   if (flg2) {
13120841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1313e5c89e4eSSatish Balay   }
131467584ceeSBarry Smith #endif
1315e5c89e4eSSatish Balay 
1316e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
131790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1318c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1319e5c89e4eSSatish Balay   if (flg1) {
1320e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1321205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1322e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1323e5c89e4eSSatish Balay   }
1324e5c89e4eSSatish Balay #endif
1325e5c89e4eSSatish Balay 
1326e5c89e4eSSatish Balay 
1327e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1328e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1329e5c89e4eSSatish Balay   mname[0] = 0;
1330589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1331e5c89e4eSSatish Balay   if (flg1) {
1332e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1333e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1334e5c89e4eSSatish Balay   }
1335e5c89e4eSSatish Balay #endif
1336356e5837SBarry Smith #endif
1337a297a907SKarl Rupp 
1338dd710f27SBarry Smith   /*
1339dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1340dd710f27SBarry Smith   */
1341dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1342dd710f27SBarry Smith 
1343356e5837SBarry Smith #if defined(PETSC_USE_LOG)
134487c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1345f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
134687c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
134787c3beb6SLisandro Dalcin 
1348356e5837SBarry Smith   mname[0] = 0;
1349589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1350e5c89e4eSSatish Balay   if (flg1) {
135191eabc43SBarry Smith     PetscViewer viewer;
135220a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
135391eabc43SBarry Smith     if (mname[0]) {
135491eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
135591eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13566bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
135733f85c2fSBarry Smith     } else {
135833f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13599a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
136033f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13619a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
136233f85c2fSBarry Smith     }
1363e5c89e4eSSatish Balay   }
1364a297a907SKarl Rupp 
1365dd710f27SBarry Smith   /*
1366dd710f27SBarry Smith      Free any objects created by the last block of code.
1367dd710f27SBarry Smith   */
1368dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1369dd710f27SBarry Smith 
1370dd710f27SBarry Smith   mname[0] = 0;
1371589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1372589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2);CHKERRQ(ierr);
13737ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1374e5c89e4eSSatish Balay #endif
137510463e74SBarry Smith 
137633f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
137710463e74SBarry Smith 
137890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1379c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1380e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
138190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1382c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1383e5c89e4eSSatish Balay   if (flg1) {
1384e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1385e5c89e4eSSatish Balay   }
138690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
138790d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13888bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1389c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1390c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1391c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1392e4c476e2SSatish Balay 
1393e5c89e4eSSatish Balay   if (flg2) {
1394be56827dSJed Brown     PetscViewer viewer;
139502ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
139602ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1397c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1398be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1399e5c89e4eSSatish Balay   }
1400e5c89e4eSSatish Balay 
1401e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1402c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1403c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1404e5c89e4eSSatish Balay 
140533fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1406c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
1407*9245e749SBarry Smith   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1408e5c89e4eSSatish Balay   if (flg3) {
14093de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1410be56827dSJed Brown       PetscViewer viewer;
141102ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
141202ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1413c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1414be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1415e5c89e4eSSatish Balay     }
14163de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
14173de2bfdfSBarry Smith     if (nopt) {
14183de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
14193de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
14203de2bfdfSBarry Smith       if (nopt == 1) {
1421e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1422e5c89e4eSSatish Balay       } else {
14237582186dSLisandro Dalcin         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1424e5c89e4eSSatish Balay       }
14253de2bfdfSBarry Smith     } else if (flg3 && flg1) {
14263de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1427df12ba86SBarry Smith     }
1428c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1429e5c89e4eSSatish Balay   }
1430e5c89e4eSSatish Balay 
1431e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1432d45a07a7SBarry Smith   if (!PetscGlobalRank) {
143387f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
143416ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1435d45a07a7SBarry Smith   }
1436ec957eceSBarry Smith #endif
1437ec957eceSBarry Smith 
14384097062eSBarry Smith #if defined(PETSC_USE_LOG)
143910463e74SBarry Smith   /*
1440dbc8283eSBarry Smith        List all objects the user may have forgot to free
14412eff7a51SBarry Smith   */
144205df10baSBarry Smith   if (PetscObjectsLog) {
1443c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1444a64a8e02SBarry Smith     if (flg1) {
1445a64a8e02SBarry Smith       MPI_Comm local_comm;
14467eb1d149SBarry Smith       char     string[64];
1447a64a8e02SBarry Smith 
1448589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL);CHKERRQ(ierr);
1449a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1450a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14517eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1452a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1453a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14540a1571b3SBarry Smith     }
145505df10baSBarry Smith   }
14564097062eSBarry Smith #endif
14574097062eSBarry Smith 
14584097062eSBarry Smith #if defined(PETSC_USE_LOG)
1459dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1460dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1461a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14624097062eSBarry Smith #endif
14632eff7a51SBarry Smith 
146433f85c2fSBarry Smith   /*
146533f85c2fSBarry Smith      Destroy any packages that registered a finalize
146633f85c2fSBarry Smith   */
146733f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
146833f85c2fSBarry Smith 
1469101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1470fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1471101409b8SToby Isaac #endif
1472101409b8SToby Isaac 
147333f85c2fSBarry Smith   /*
147448dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
147548dd1dffSBarry Smith 
147637e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
147748dd1dffSBarry Smith   */
147837e93019SBarry Smith 
14794028d114SSatish Balay   if (petsc_history) {
1480f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
148102c9f0b5SLisandro Dalcin     petsc_history = NULL;
1482e5c89e4eSSatish Balay   }
14839de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1484e94e781bSJacob Faibussowitsch   ierr = PetscInfoDestroy();CHKERRQ(ierr);
1485e5c89e4eSSatish Balay 
148667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
148792f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1488e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
148992f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1490e5c89e4eSSatish Balay     FILE *fd;
1491ed9cf6e9SBarry Smith     int  err;
1492e5c89e4eSSatish Balay 
1493dc92acbaSJed Brown     flg2 = PETSC_FALSE;
149492f119d6SBarry Smith     flg3 = PETSC_FALSE;
1495cf9c20a2SJed Brown     if (PetscDefined(USE_DEBUG)) {ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);}
149692f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
149792f119d6SBarry Smith     fname[0] = 0;
1498589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1499e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1500e5c89e4eSSatish Balay 
1501589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1502e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1503e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1504ed9cf6e9SBarry Smith       err  = fclose(fd);
1505e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
150692f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1507e5c89e4eSSatish Balay       MPI_Comm local_comm;
1508e5c89e4eSSatish Balay 
1509e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1510e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1511e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1512e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1513e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1514e5c89e4eSSatish Balay     }
1515e5c89e4eSSatish Balay     fname[0] = 0;
1516589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1517e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1518e5c89e4eSSatish Balay 
1519589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
152092f119d6SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
152192f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1522ed9cf6e9SBarry Smith       err  = fclose(fd);
1523e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
152492f119d6SBarry Smith     } else if (flg1) {
152592f119d6SBarry Smith       MPI_Comm local_comm;
152692f119d6SBarry Smith 
152792f119d6SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
152892f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
152992f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
153092f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
153192f119d6SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1532e5c89e4eSSatish Balay     }
1533e5c89e4eSSatish Balay   }
153467584ceeSBarry Smith #endif
153520e2c332SMatthew G. Knepley 
15365486ca60SMatthew G. Knepley   /*
15375486ca60SMatthew G. Knepley      Close any open dynamic libraries
15385486ca60SMatthew G. Knepley   */
15395486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15405486ca60SMatthew G. Knepley 
1541e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15424416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1543e5c89e4eSSatish Balay 
1544e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
154502c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1546e5c89e4eSSatish Balay 
1547008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1548e5c89e4eSSatish Balay 
1549dbc8283eSBarry Smith   /*
1550efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1551efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1552efb80d3cSBarry Smith 
1553efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1554efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1555dbc8283eSBarry Smith  */
1556b770b1f6SSatish Balay   {
1557dbc8283eSBarry Smith     PetscCommCounter *counter;
1558dbc8283eSBarry Smith     PetscMPIInt      flg;
1559dbc8283eSBarry Smith     MPI_Comm         icomm;
1560265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
156147435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1562dbc8283eSBarry Smith     if (flg) {
1563265f3f35SJed Brown       icomm = ucomm.comm;
156447435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1565dbc8283eSBarry 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");
1566dbc8283eSBarry Smith 
156747435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
156847435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1569efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1570dbc8283eSBarry Smith     }
157147435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1572dbc8283eSBarry Smith     if (flg) {
1573265f3f35SJed Brown       icomm = ucomm.comm;
157447435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1575dbc8283eSBarry 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");
1576dbc8283eSBarry Smith 
157747435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
157847435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1579efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1580dbc8283eSBarry Smith     }
1581b770b1f6SSatish Balay   }
1582dbc8283eSBarry Smith 
158347435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
158447435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
158547435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
15865f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1587480cf27aSJed Brown 
15885ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15895ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15905ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1591ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1592ef19f930SBarry Smith 
1593e5c89e4eSSatish Balay   if (PetscBeganMPI) {
159499608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
159599b1327fSBarry Smith     PetscMPIInt flag;
159699b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1597e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
159899608316SBarry Smith #endif
1599e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1600e5c89e4eSSatish Balay   }
1601e5c89e4eSSatish Balay /*
1602e5c89e4eSSatish Balay 
1603e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1604e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1605e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1606e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1607e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16080e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1609e5c89e4eSSatish Balay    memory was not freed.
1610e5c89e4eSSatish Balay 
1611e5c89e4eSSatish Balay */
16121d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1613a297a907SKarl Rupp 
16148ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_FALSE;
1615e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1616e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
16173db9a53dSBarry Smith   PetscFunctionReturn(0);
1618e5c89e4eSSatish Balay }
1619e5c89e4eSSatish Balay 
162043db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16218cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
162243db4dbbSBarry Smith {
162343db4dbbSBarry Smith   if (*a == *b) return 1;
162443db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
162543db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
162643db4dbbSBarry Smith   return 0;
162743db4dbbSBarry Smith }
1628a70650f6SBarry Smith #endif
162943db4dbbSBarry Smith 
163043db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16318cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
163243db4dbbSBarry Smith {
163343db4dbbSBarry Smith   if (*a == *b) return 1;
163443db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
163543db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
163643db4dbbSBarry Smith   return 0;
163743db4dbbSBarry Smith }
163843db4dbbSBarry Smith #endif
1639