xref: /petsc/src/sys/objects/pinit.c (revision de272c7ab8a1acfc6d767cecdf5f3b7cfc1e1d55)
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>
856883f60SBarry Smith #if defined(PETSC_USE_GCOV)
956883f60SBarry Smith EXTERN_C_BEGIN
1056883f60SBarry Smith void  __gcov_flush(void);
1156883f60SBarry Smith EXTERN_C_END
1256883f60SBarry Smith #endif
138101f56cSMatthew Knepley 
14a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
16a9f03627SSatish Balay #endif
17f2d66bcaSShri Abhyankar 
182d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1995c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
202d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
212d53ad75SBarry Smith #endif
222d53ad75SBarry Smith 
23a6790183SBarry Smith #if defined(PETSC_HAVE_SAWS)
2416ad0300SBarry Smith #include <petscviewersaws.h>
25a6790183SBarry Smith #endif
26eb02082bSJunchao Zhang 
27e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
28e5c89e4eSSatish Balay 
2995c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
30e5c89e4eSSatish Balay 
3195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
3395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
3495c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
3695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
370069ddf5SShri Abhyankar 
386de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
39e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
406de5d289SStefano Zampini #if defined(PETSC_HAVE_MPI_INIT_THREAD)
416de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
426de5d289SStefano Zampini #else
436de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
446de5d289SStefano Zampini #endif
45e5c89e4eSSatish Balay 
46480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
47480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
48480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
495f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
50480cf27aSJed Brown 
51e5c89e4eSSatish Balay /*
52e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
53e5c89e4eSSatish Balay */
5402c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
5502c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
56e5c89e4eSSatish Balay 
57ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
58ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
590f8e0872SSatish Balay 
60a2f94806SJed Brown PetscInt PetscHotRegionDepth;
61a2f94806SJed Brown 
62b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
63b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
64b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
65b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
66b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
67b22622e2STadeu Manoel #endif
68b22622e2STadeu Manoel 
69e5c89e4eSSatish Balay /*
70945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
7172a42c3cSBarry Smith 
7272a42c3cSBarry Smith    Collective
7372a42c3cSBarry Smith 
7472a42c3cSBarry Smith    Level: advanced
7572a42c3cSBarry Smith 
7695452b02SPatrick Sanan     Notes:
77a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
780f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
79a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
800f11a792SBarry Smith 
81a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
821ea5a559SBarry Smith 
8372a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
840f11a792SBarry Smith */
85945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
8672a42c3cSBarry Smith {
8772a42c3cSBarry Smith   PetscErrorCode ierr;
8872a42c3cSBarry Smith   int            myargc   = argc;
8972a42c3cSBarry Smith   char           **myargs = args;
9072a42c3cSBarry Smith 
9172a42c3cSBarry Smith   PetscFunctionBegin;
92c52ac268SRichard Tran Mills   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr;
931ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
94df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
9572a42c3cSBarry Smith   PetscFunctionReturn(ierr);
9672a42c3cSBarry Smith }
9772a42c3cSBarry Smith 
98f0865b08SBarry Smith /*
99a56f64adSBarry Smith       Used by Julia interface to get communicator
100f0865b08SBarry Smith */
101945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
102f0865b08SBarry Smith {
103f0865b08SBarry Smith   PetscFunctionBegin;
104f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
105f0865b08SBarry Smith   PetscFunctionReturn(0);
106f0865b08SBarry Smith }
107f0865b08SBarry Smith 
108e5c89e4eSSatish Balay /*@C
109e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
110e5c89e4eSSatish Balay         the command line arguments.
111e5c89e4eSSatish Balay 
112e5c89e4eSSatish Balay    Collective
113e5c89e4eSSatish Balay 
114e5c89e4eSSatish Balay    Level: advanced
115e5c89e4eSSatish Balay 
116e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
117e5c89e4eSSatish Balay @*/
1187087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
119e5c89e4eSSatish Balay {
120e5c89e4eSSatish Balay   PetscErrorCode ierr;
121e5c89e4eSSatish Balay   int            argc   = 0;
12202c9f0b5SLisandro Dalcin   char           **args = NULL;
123e5c89e4eSSatish Balay 
124e5c89e4eSSatish Balay   PetscFunctionBegin;
1250298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
126e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
127e5c89e4eSSatish Balay }
128e5c89e4eSSatish Balay 
129e5c89e4eSSatish Balay /*@
130e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
131e5c89e4eSSatish Balay 
13293b6d2d1SJed Brown    Level: beginner
133e5c89e4eSSatish Balay 
134e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
135e5c89e4eSSatish Balay @*/
1367087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized)
137e5c89e4eSSatish Balay {
138e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
13993b6d2d1SJed Brown   return 0;
140e5c89e4eSSatish Balay }
141e5c89e4eSSatish Balay 
142e5c89e4eSSatish Balay /*@
143e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
144e5c89e4eSSatish Balay 
145e5c89e4eSSatish Balay    Level: developer
146e5c89e4eSSatish Balay 
147e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
148e5c89e4eSSatish Balay @*/
1497087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
150e5c89e4eSSatish Balay {
151e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
15293b6d2d1SJed Brown   return 0;
153e5c89e4eSSatish Balay }
154e5c89e4eSSatish Balay 
15557171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []);
156e5c89e4eSSatish Balay 
157e5c89e4eSSatish Balay /*
158e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
159e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
160e5c89e4eSSatish Balay */
161367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
162e5c89e4eSSatish Balay 
163367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
164e5c89e4eSSatish Balay {
165e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
166e5c89e4eSSatish Balay 
167e5c89e4eSSatish Balay   PetscFunctionBegin;
168e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
169e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
17041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
171e5c89e4eSSatish Balay   }
172e5c89e4eSSatish Balay 
173e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
174e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
175e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
176e5c89e4eSSatish Balay   }
177812af9f3SBarry Smith   PetscFunctionReturnVoid();
178e5c89e4eSSatish Balay }
179e5c89e4eSSatish Balay 
180e5c89e4eSSatish Balay /*
181e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
182e5c89e4eSSatish Balay sum of the second entry.
183b693b147SBarry Smith 
18476f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
185367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
186b693b147SBarry Smith there would be no place to store the both needed results.
187e5c89e4eSSatish Balay */
18876ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
189e5c89e4eSSatish Balay {
190e5c89e4eSSatish Balay   PetscErrorCode ierr;
191e5c89e4eSSatish Balay 
192e5c89e4eSSatish Balay   PetscFunctionBegin;
193d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
194d6e4c47cSJed Brown   {
195d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
196ffc4695bSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRMPI(ierr);
197d6e4c47cSJed Brown     *max = work.max;
198d6e4c47cSJed Brown     *sum = work.sum;
199d6e4c47cSJed Brown   }
200d6e4c47cSJed Brown #else
201d6e4c47cSJed Brown   {
202d6e4c47cSJed Brown     PetscMPIInt    size,rank;
203d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
204ffc4695bSBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
205ffc4695bSBarry Smith     ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
206785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
207367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2086ac3741eSJed Brown     *max = work[rank].max;
2096ac3741eSJed Brown     *sum = work[rank].sum;
210e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
211d6e4c47cSJed Brown   }
212d6e4c47cSJed Brown #endif
213e5c89e4eSSatish Balay   PetscFunctionReturn(0);
214e5c89e4eSSatish Balay }
215e5c89e4eSSatish Balay 
216e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
217e5c89e4eSSatish Balay 
218*de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
21906a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
220e5c89e4eSSatish Balay 
2218cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
222e5c89e4eSSatish Balay {
223e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
224e5c89e4eSSatish Balay 
225e5c89e4eSSatish Balay   PetscFunctionBegin;
2267c2de775SJed Brown   if (*datatype == MPIU_REAL) {
227e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
228a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2297c2de775SJed Brown   }
2307c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
231*de272c7aSSatish Balay   else if (*datatype == MPI_COMPLEX) {
2327c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
233a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2347c2de775SJed Brown   }
2357c2de775SJed Brown #endif
2367c2de775SJed Brown   else {
2377c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
23841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
239e2e03761SBarry Smith   }
240812af9f3SBarry Smith   PetscFunctionReturnVoid();
241e5c89e4eSSatish Balay }
242e5c89e4eSSatish Balay #endif
243e5c89e4eSSatish Balay 
244570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
245d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
246d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
247d9822059SBarry Smith 
2488cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
249d9822059SBarry Smith {
250d9822059SBarry Smith   PetscInt i,count = *cnt;
251d9822059SBarry Smith 
252d9822059SBarry Smith   PetscFunctionBegin;
2537c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2548c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
255a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2567c2de775SJed Brown   }
2577c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2587c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2597c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2607c2de775SJed Brown     for (i=0; i<count; i++) {
2617c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2627c2de775SJed Brown     }
2637c2de775SJed Brown   }
2647c2de775SJed Brown #endif
2657c2de775SJed Brown   else {
2667c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
26741e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2688c764dc5SJose Roman   }
269d9822059SBarry Smith   PetscFunctionReturnVoid();
270d9822059SBarry Smith }
271d9822059SBarry Smith 
2728cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
273d9822059SBarry Smith {
274d9822059SBarry Smith   PetscInt    i,count = *cnt;
275d9822059SBarry Smith 
276d9822059SBarry Smith   PetscFunctionBegin;
2777c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2788c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
279a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2807c2de775SJed Brown   }
2817c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2827c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2837c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2847c2de775SJed Brown     for (i=0; i<count; i++) {
2857c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2867c2de775SJed Brown     }
2877c2de775SJed Brown   }
2887c2de775SJed Brown #endif
2897c2de775SJed Brown   else {
2908c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
29141e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2928c764dc5SJose Roman   }
293d9822059SBarry Smith   PetscFunctionReturnVoid();
294d9822059SBarry Smith }
295d9822059SBarry Smith #endif
296d9822059SBarry Smith 
297480cf27aSJed Brown /*
298480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
299480cf27aSJed Brown 
300ff0e51ddSBarry 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.
301480cf27aSJed Brown 
30212801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
303480cf27aSJed Brown 
304480cf27aSJed Brown */
30533779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
306480cf27aSJed Brown {
307480cf27aSJed Brown   PetscErrorCode   ierr;
30805342407SJunchao Zhang   PetscCommCounter *counter=(PetscCommCounter*)count_val;
309480cf27aSJed Brown 
310480cf27aSJed Brown   PetscFunctionBegin;
31102c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
31205342407SJunchao Zhang   ierr = PetscFree(counter->iflags);CHKERRMPI(ierr);
31305342407SJunchao Zhang   ierr = PetscFree(counter);CHKERRMPI(ierr);
314480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
315480cf27aSJed Brown }
316480cf27aSJed Brown 
317480cf27aSJed Brown /*
31847435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
319da3039f7SJed Brown   calls MPI_Comm_free().
320da3039f7SJed Brown 
321da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
322480cf27aSJed Brown 
323ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
324480cf27aSJed Brown 
32512801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
326480cf27aSJed Brown 
327480cf27aSJed Brown */
32833779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
329480cf27aSJed Brown {
330480cf27aSJed Brown   PetscErrorCode                    ierr;
33133779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
332480cf27aSJed Brown 
333480cf27aSJed Brown   PetscFunctionBegin;
33412801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
335ec4fadc2SJed Brown   icomm.ptr = attr_val;
33676bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
33733779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
33833779a13SJunchao Zhang     PetscMPIInt flg;
33933779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
34047435625SJed Brown     ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
34133779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
34233779a13SJunchao 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");
34333779a13SJunchao Zhang   }
34447435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
34533779a13SJunchao Zhang   ierr = PetscInfo2(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);CHKERRMPI(ierr);
346da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
347b89831e5SBarry Smith }
348da3039f7SJed Brown 
349da3039f7SJed Brown /*
35033779a13SJunchao 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.
351da3039f7SJed Brown  */
35233779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
353da3039f7SJed Brown {
354da3039f7SJed Brown   PetscErrorCode ierr;
355da3039f7SJed Brown 
356da3039f7SJed Brown   PetscFunctionBegin;
35702c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
358480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
359480cf27aSJed Brown }
360480cf27aSJed Brown 
36133779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
36242218b76SBarry Smith 
363951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3648cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3658cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3668cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
367e39fd77fSBarry Smith #endif
368e39fd77fSBarry Smith 
3690f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
37012801b39SBarry Smith 
371eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
372eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3736ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
37402c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
375dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
376e5c89e4eSSatish Balay 
377dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
378051e4cf2SJed Brown {
379051e4cf2SJed Brown   PetscErrorCode ierr;
380051e4cf2SJed Brown 
381051e4cf2SJed Brown   PetscFunctionBegin;
382051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
383a1601671SBarry 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);
384051e4cf2SJed 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);
385051e4cf2SJed Brown   PetscFunctionReturn(0);
386051e4cf2SJed Brown }
387e5c89e4eSSatish Balay 
3882d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3892d747510SLisandro Dalcin 
3902d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
3912d747510SLisandro Dalcin {
3922d747510SLisandro Dalcin   PetscErrorCode ierr;
3932d747510SLisandro Dalcin 
3942d747510SLisandro Dalcin   PetscFunctionBegin;
395589a23caSBarry Smith   ierr  = PetscStrncpy(programname,name,sizeof(programname));CHKERRQ(ierr);
3962d747510SLisandro Dalcin   PetscFunctionReturn(0);
3972d747510SLisandro Dalcin }
3982d747510SLisandro Dalcin 
3992d747510SLisandro Dalcin /*@C
4002d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4012d747510SLisandro Dalcin 
4022d747510SLisandro Dalcin     Not Collective
4032d747510SLisandro Dalcin 
4042d747510SLisandro Dalcin     Input Parameter:
4052d747510SLisandro Dalcin .   len - length of the string name
4062d747510SLisandro Dalcin 
4072d747510SLisandro Dalcin     Output Parameter:
4082d747510SLisandro Dalcin .   name - the name of the running program
4092d747510SLisandro Dalcin 
4102d747510SLisandro Dalcin    Level: advanced
4112d747510SLisandro Dalcin 
4122d747510SLisandro Dalcin     Notes:
4132d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4142d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4152d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4162d747510SLisandro Dalcin @*/
4172d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4182d747510SLisandro Dalcin {
4192d747510SLisandro Dalcin   PetscErrorCode ierr;
4202d747510SLisandro Dalcin 
4212d747510SLisandro Dalcin   PetscFunctionBegin;
4222d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4232d747510SLisandro Dalcin   PetscFunctionReturn(0);
4242d747510SLisandro Dalcin }
4252d747510SLisandro Dalcin 
426e5c89e4eSSatish Balay /*@C
427e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
428e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
429e5c89e4eSSatish Balay 
430e5c89e4eSSatish Balay    Not Collective
431e5c89e4eSSatish Balay 
432e5c89e4eSSatish Balay    Output Parameters:
433e5c89e4eSSatish Balay +  argc - count of number of command line arguments
434e5c89e4eSSatish Balay -  args - the command line arguments
435e5c89e4eSSatish Balay 
436e5c89e4eSSatish Balay    Level: intermediate
437e5c89e4eSSatish Balay 
438e5c89e4eSSatish Balay    Notes:
439e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
440e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
441e5c89e4eSSatish Balay 
442f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
443f177e3b1SBarry Smith 
444793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
445e5c89e4eSSatish Balay 
446e5c89e4eSSatish Balay @*/
4477087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
448e5c89e4eSSatish Balay {
449e5c89e4eSSatish Balay   PetscFunctionBegin;
45017186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
451e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
452e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
453e5c89e4eSSatish Balay   PetscFunctionReturn(0);
454e5c89e4eSSatish Balay }
455e5c89e4eSSatish Balay 
456793721a6SBarry Smith /*@C
457793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
458793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
459793721a6SBarry Smith 
460793721a6SBarry Smith    Not Collective
461793721a6SBarry Smith 
462793721a6SBarry Smith    Output Parameters:
463793721a6SBarry Smith .  args - the command line arguments
464793721a6SBarry Smith 
465793721a6SBarry Smith    Level: intermediate
466793721a6SBarry Smith 
467793721a6SBarry Smith    Notes:
468793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
469793721a6SBarry Smith 
470793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
471793721a6SBarry Smith 
472793721a6SBarry Smith @*/
4737087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
474793721a6SBarry Smith {
475793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
476793721a6SBarry Smith   PetscErrorCode ierr;
477793721a6SBarry Smith 
478793721a6SBarry Smith   PetscFunctionBegin;
47917186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4802d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
481785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
482793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
483793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
484793721a6SBarry Smith   }
4852d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
486793721a6SBarry Smith   PetscFunctionReturn(0);
487793721a6SBarry Smith }
488793721a6SBarry Smith 
489793721a6SBarry Smith /*@C
490793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
491793721a6SBarry Smith 
492793721a6SBarry Smith    Not Collective
493793721a6SBarry Smith 
494793721a6SBarry Smith    Output Parameters:
495793721a6SBarry Smith .  args - the command line arguments
496793721a6SBarry Smith 
497793721a6SBarry Smith    Level: intermediate
498793721a6SBarry Smith 
499793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
500793721a6SBarry Smith 
501793721a6SBarry Smith @*/
5027087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
503793721a6SBarry Smith {
504793721a6SBarry Smith   PetscInt       i = 0;
505793721a6SBarry Smith   PetscErrorCode ierr;
506793721a6SBarry Smith 
507793721a6SBarry Smith   PetscFunctionBegin;
508a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
509793721a6SBarry Smith   while (args[i]) {
510793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
511793721a6SBarry Smith     i++;
512793721a6SBarry Smith   }
513793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
514793721a6SBarry Smith   PetscFunctionReturn(0);
515793721a6SBarry Smith }
516793721a6SBarry Smith 
51711525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
51830befbd2SBarry Smith #include <petscconfiginfo.h>
51930befbd2SBarry Smith 
52095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
52111525c0dSBarry Smith {
52211525c0dSBarry Smith   if (!PetscGlobalRank) {
52330befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
52411525c0dSBarry Smith     int            port;
525ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
52611525c0dSBarry Smith     size_t         applinelen,introlen;
52711525c0dSBarry Smith     PetscErrorCode ierr;
528ffbd1cfbSBarry Smith     char           sawsurl[256];
52911525c0dSBarry Smith 
530c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
53111525c0dSBarry Smith     if (flg) {
53211525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
53311525c0dSBarry Smith 
534589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL);CHKERRQ(ierr);
53511525c0dSBarry Smith       if (sawslog[0]) {
53611525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
53711525c0dSBarry Smith       } else {
53811525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
53911525c0dSBarry Smith       }
54011525c0dSBarry Smith     }
541589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg);CHKERRQ(ierr);
54211525c0dSBarry Smith     if (flg) {
54311525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
54411525c0dSBarry Smith     }
545c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
546ffbd1cfbSBarry Smith     if (selectport) {
547ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
548ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
549ffbd1cfbSBarry Smith     } else {
550c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
55111525c0dSBarry Smith       if (flg) {
55211525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
55311525c0dSBarry Smith       }
554ffbd1cfbSBarry Smith     }
555589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg);CHKERRQ(ierr);
55611525c0dSBarry Smith     if (flg) {
55711525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
55811525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5599c1e0ce8SBarry Smith     } else {
560c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5619c1e0ce8SBarry Smith       if (flg) {
562589a23caSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root));CHKERRQ(ierr);
5639c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5649c1e0ce8SBarry Smith       }
56511525c0dSBarry Smith     }
566c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
56711525c0dSBarry Smith     if (flg2) {
56811525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
56911525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
570589a23caSBarry Smith       ierr = PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root);CHKERRQ(ierr);
57111525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
57211525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
57343da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
57411525c0dSBarry Smith     }
575589a23caSBarry Smith     ierr = PetscGetProgramName(programname,sizeof(programname));CHKERRQ(ierr);
57611525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
57711525c0dSBarry Smith     introlen   = 4096 + applinelen;
57830a8c9c0SSurtai Han     applinelen += 1024;
57911525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
58011525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
58111525c0dSBarry Smith 
58211525c0dSBarry Smith     if (rootlocal) {
58311525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
58411525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
58511525c0dSBarry Smith     }
58676a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
58711525c0dSBarry Smith     if (rootlocal && help) {
588928bb9adSStefano 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);
58911525c0dSBarry Smith     } else if (help) {
590928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
59111525c0dSBarry Smith     } else {
592928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
59311525c0dSBarry Smith     }
594b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
59530befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
59611525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
597a8d69d7bSBarry 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"
598df62c222SBarry 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"
599928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
60043da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
60111525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
60211525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
603ffbd1cfbSBarry Smith     if (selectport) {
604aa573868SBarry Smith       PetscBool silent;
6057d812c46SBarry Smith 
6067d812c46SBarry Smith       ierr = SAWs_Initialize();
6077d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6087d812c46SBarry Smith       while (ierr) {
6097d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6107d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6117d812c46SBarry Smith         ierr = SAWs_Initialize();
6127d812c46SBarry Smith       }
6137d812c46SBarry Smith 
614aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
615aa573868SBarry Smith       if (!silent) {
616ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
617ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
618ffbd1cfbSBarry Smith       }
6197d812c46SBarry Smith     } else {
6207d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
621aa573868SBarry Smith     }
6220af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6230af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6240af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6250af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6260af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
62711525c0dSBarry Smith   }
62811525c0dSBarry Smith   PetscFunctionReturn(0);
62911525c0dSBarry Smith }
63011525c0dSBarry Smith #endif
63111525c0dSBarry Smith 
6324dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6334dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6344dfee713SSatish Balay {
6354dfee713SSatish Balay   PetscFunctionBegin;
6364dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6374dfee713SSatish Balay     /* see MPI.py for details on this bug */
6384dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6394dfee713SSatish Balay #endif
6404dfee713SSatish Balay   PetscFunctionReturn(0);
6414dfee713SSatish Balay }
6424dfee713SSatish Balay 
643a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
644a56f64adSBarry Smith #include <adios.h>
64522580e64SBarry Smith #include <adios_read.h>
6467b56e58cSSatish Balay int64_t Petsc_adios_group;
647a56f64adSBarry Smith #endif
64855d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
64955d657eeSBarry Smith #include <adios2_c.h>
65055d657eeSBarry Smith #endif
651cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
652cd1b99a6SBarry Smith #include <omp.h>
653f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
654cd1b99a6SBarry Smith #endif
655a56f64adSBarry Smith 
656bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H)
657bc8a24ecSBarry Smith #include <dlfcn.h>
658bc8a24ecSBarry Smith #endif
659bc8a24ecSBarry Smith 
660e5c89e4eSSatish Balay /*@C
661e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
662e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
663e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
664e5c89e4eSSatish Balay    your program -- usually the very first line!
665e5c89e4eSSatish Balay 
666e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
667e5c89e4eSSatish Balay 
668e5c89e4eSSatish Balay    Input Parameters:
669e5c89e4eSSatish Balay +  argc - count of number of command line arguments
670e5c89e4eSSatish Balay .  args - the command line arguments
671be10d61cSLisandro Dalcin .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
672be10d61cSLisandro Dalcin           Use NULL or empty string to not check for code specific file.
673be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
674c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
6750298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
676e5c89e4eSSatish Balay 
67705827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
67805827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
67905827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
68005827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
68105827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
682e5c89e4eSSatish Balay 
683e5c89e4eSSatish Balay    Options Database Keys:
6847ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6857ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
686e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
6878a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
6888a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
6898a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
6908a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
6918a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
692bf4d2887SBarry Smith .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
693e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
694e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
695e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
69679dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
69779dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
69879dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
699aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
70092f119d6SBarry 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
70192f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
70292f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
703608c71bfSMatthew G. Knepley .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
70492f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
705e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
706e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
707e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
708e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
709e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7100841954dSBarry Smith -  -memory_view - Print memory usage at end of run
711e5c89e4eSSatish Balay 
712c5b5d8d5SVaclav Hapla    Options Database Keys for Option Database:
713c5b5d8d5SVaclav Hapla +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
714c5b5d8d5SVaclav Hapla .  -options_monitor - monitor all set options to standard output for the whole program run
715c5b5d8d5SVaclav Hapla -  -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet()
716c5b5d8d5SVaclav Hapla 
717c5b5d8d5SVaclav Hapla    Options -options_monitor_{all,cancel} are
718c5b5d8d5SVaclav Hapla    position-independent and apply to all options set since the PETSc start.
719c5b5d8d5SVaclav Hapla    They can be used also in option files.
720c5b5d8d5SVaclav Hapla 
721c5b5d8d5SVaclav Hapla    See PetscOptionsMonitorSet() to do monitoring programmatically.
722c5b5d8d5SVaclav Hapla 
723e5c89e4eSSatish Balay    Options Database Keys for Profiling:
724a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
725fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
726217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
727217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
728495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
729e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7309a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
73179dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
7329a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
733495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
734f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
735495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
736495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
737c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
73887c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
739c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
740495fc317SBarry Smith 
741609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
742e5c89e4eSSatish Balay 
743ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
744ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
745ffbd1cfbSBarry 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
746ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
747ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
748ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
749ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
750ffbd1cfbSBarry Smith 
751e5c89e4eSSatish Balay    Environmental Variables:
752e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
753e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
754e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
7554a971ea4SToby Isaac .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
7564a971ea4SToby Isaac .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
757e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
758e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
759e5c89e4eSSatish Balay 
760e5c89e4eSSatish Balay 
761e5c89e4eSSatish Balay    Level: beginner
762e5c89e4eSSatish Balay 
763e5c89e4eSSatish Balay    Notes:
764e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
765e5c89e4eSSatish Balay    it before PetscInitialize().
766e5c89e4eSSatish Balay 
767e5c89e4eSSatish Balay    Fortran Version:
768e5c89e4eSSatish Balay    In Fortran this routine has the format
769e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
770e5c89e4eSSatish Balay 
771e5c89e4eSSatish Balay +   ierr - error return code
772c5b5d8d5SVaclav Hapla -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
773c5b5d8d5SVaclav Hapla           Use PETSC_NULL_CHARACTER to not check for code specific file.
774c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
775e5c89e4eSSatish Balay 
776e5c89e4eSSatish Balay    Important Fortran Note:
7770eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7780298fd71SBarry Smith    null character string; you CANNOT just use NULL as
779a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
780e5c89e4eSSatish Balay 
78101cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
78201cb0274SBarry Smith    calling PetscInitialize().
783e5c89e4eSSatish Balay 
78401cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
785e5c89e4eSSatish Balay 
786e5c89e4eSSatish Balay @*/
7877087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
788e5c89e4eSSatish Balay {
789e5c89e4eSSatish Balay   PetscErrorCode ierr;
7904bb5149bSJed Brown   PetscMPIInt    flag, size;
79119c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
792e5c89e4eSSatish Balay   char           hostname[256];
793e5c89e4eSSatish Balay 
794e5c89e4eSSatish Balay   PetscFunctionBegin;
795e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7963d96e996SBarry Smith   /*
7973d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7983d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7993d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
8003d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
801110fc3b0SBarry Smith         Intel MPI Library v5.0 (2014)
8023d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
8033d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
8043d96e996SBarry Smith       listed above and since that time are compatible.
8053d96e996SBarry Smith 
8063d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
8073d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
8083d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
8093d96e996SBarry Smith       to perform the checking.
8103d96e996SBarry Smith 
8113d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
8123d96e996SBarry Smith 
8133d96e996SBarry Smith       Questions:
8143d96e996SBarry Smith 
8153d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
8163d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
8173d96e996SBarry Smith   */
8183d96e996SBarry Smith 
81919c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
82019c5658aSBarry Smith   {
82119c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
82219c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
82319c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8243d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
82519c5658aSBarry Smith #if defined(MPICH_VERSION)
8263d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
82719c5658aSBarry Smith     {
82819c5658aSBarry Smith       char *ver,*lf;
82919c5658aSBarry Smith       flg = PETSC_FALSE;
83019c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
83119c5658aSBarry Smith       if (ver) {
83219c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
83319c5658aSBarry Smith         if (lf) {
83419c5658aSBarry Smith           *lf = 0;
83519c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
83619c5658aSBarry Smith         }
83719c5658aSBarry Smith       }
83819c5658aSBarry Smith       if (!flg) {
83919c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8403d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
84119c5658aSBarry Smith       }
84219c5658aSBarry Smith     }
8433d96e996SBarry Smith #endif
8443d96e996SBarry 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?) */
84519c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
84619c5658aSBarry Smith     {
84716dc8964SSatish Balay       char *ver,bs[MPI_MAX_LIBRARY_VERSION_STRING],*bsf;
84819c5658aSBarry Smith       flg = PETSC_FALSE;
84916dc8964SSatish Balay #define PSTRSZ 2
85016dc8964SSatish Balay       char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI","FUJITSU MPI"};
85116dc8964SSatish Balay       char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v","Library "};
85216dc8964SSatish Balay       int i;
85316dc8964SSatish Balay       for (i=0; i<PSTRSZ; i++) {
85416dc8964SSatish Balay         ierr = PetscStrstr(mpilibraryversion,ompistr1[i],&ver);if (ierr) return ierr;
85519c5658aSBarry Smith         if (ver) {
85616dc8964SSatish Balay           PetscSNPrintf(bs,MPI_MAX_LIBRARY_VERSION_STRING,"%s%d.%d",ompistr2[i],OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
85719c5658aSBarry Smith           ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
85819c5658aSBarry Smith           if (bsf) flg = PETSC_TRUE;
85916dc8964SSatish Balay           break;
86016dc8964SSatish Balay         }
86119c5658aSBarry Smith       }
86219c5658aSBarry Smith       if (!flg) {
86319c5658aSBarry 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);
8643d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
86519c5658aSBarry Smith       }
86619c5658aSBarry Smith     }
86719c5658aSBarry Smith #endif
86819c5658aSBarry Smith   }
86919c5658aSBarry Smith #endif
87019c5658aSBarry Smith 
871bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM)
872bc8a24ecSBarry Smith   {
873bc8a24ecSBarry Smith     PetscInt cnt = 0;
874bc8a24ecSBarry 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 */
875bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++;
876bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++;
877bc8a24ecSBarry Smith     if (cnt > 1) {
878bc8a24ecSBarry Smith       fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
879bc8a24ecSBarry Smith       return PETSC_ERR_MPI_LIB_INCOMP;
880bc8a24ecSBarry Smith     }
881bc8a24ecSBarry Smith   }
882bc8a24ecSBarry Smith #endif
883bc8a24ecSBarry Smith 
884ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
885d89683f4Sbcordonn   PETSC_STDOUT = stdout;
886ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
887e5c89e4eSSatish Balay 
8888ad20175SVaclav Hapla   /*CHKERRQ can be used from now */
8898ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_TRUE;
8908ad20175SVaclav Hapla 
8910c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8920c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8930c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8940c30907bSSatish Balay #endif
8950c30907bSSatish Balay 
8964416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
897e5c89e4eSSatish Balay 
898e5c89e4eSSatish Balay   /*
899e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
900e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
901e5c89e4eSSatish Balay   */
902e5c89e4eSSatish Balay   if (argc && *argc) {
903e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
904e5c89e4eSSatish Balay   } else {
905e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
906e5c89e4eSSatish Balay   }
907e5c89e4eSSatish Balay 
908ffc4695bSBarry Smith   ierr = MPI_Initialized(&flag);CHKERRMPI(ierr);
909e5c89e4eSSatish Balay   if (!flag) {
910e32f2f54SBarry 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");
9114dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
9125e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
9135e765c61SJed Brown     {
9145e765c61SJed Brown       PetscMPIInt provided;
915ffc4695bSBarry Smith       ierr = MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided);CHKERRMPI(ierr);
9165e765c61SJed Brown     }
9175e765c61SJed Brown #else
918ffc4695bSBarry Smith     ierr = MPI_Init(argc,args);CHKERRMPI(ierr);
9195e765c61SJed Brown #endif
920e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
921e5c89e4eSSatish Balay   }
9224dfee713SSatish Balay 
923e5c89e4eSSatish Balay   if (argc && args) {
924e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
925e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
926e5c89e4eSSatish Balay   }
927e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
9285ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
9295ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
9305ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
931ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
932e5c89e4eSSatish Balay 
933a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
934ffc4695bSBarry Smith   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRMPI(ierr);
935e5c89e4eSSatish Balay 
9360f9be574SPeter Hill   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
937ffc4695bSBarry Smith     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRMPI(ierr);
938ffc4695bSBarry Smith     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRMPI(ierr);
9390f9be574SPeter Hill   }
94012801b39SBarry Smith 
941e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
942e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
943e5c89e4eSSatish Balay 
944ffc4695bSBarry Smith   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRMPI(ierr);
945ffc4695bSBarry Smith   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRMPI(ierr);
946e5c89e4eSSatish Balay 
9478ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9488ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9497cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
950e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
951e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
952e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
953e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
954e316c87fSJed Brown #endif
955e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9568ad47952SJed Brown 
957e5c89e4eSSatish Balay   /*
958e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
959e5c89e4eSSatish Balay      shared libraries the constructors for global variables
960e5c89e4eSSatish Balay      are not called; at least on IRIX.
961e5c89e4eSSatish Balay   */
962886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
963e5c89e4eSSatish Balay   {
9647a19d461SSatish Balay #if defined(PETSC_CLANGUAGE_CXX)
96550f81f78SJed Brown     PetscComplex ic(0.0,1.0);
966e5c89e4eSSatish Balay     PETSC_i = ic;
967252ecd31SSatish Balay #else
96850f81f78SJed Brown     PETSC_i = _Complex_I;
969b7940d39SSatish Balay #endif
970762437b8SSatish Balay   }
971886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
972e5c89e4eSSatish Balay 
973e5c89e4eSSatish Balay   /*
974e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
975e5c89e4eSSatish Balay      half of the entries and maxes the second half.
976e5c89e4eSSatish Balay   */
977ffc4695bSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRMPI(ierr);
978e5c89e4eSSatish Balay 
979ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
980ffc4695bSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRMPI(ierr);
981ffc4695bSBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRMPI(ierr);
9827c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
983ffc4695bSBarry Smith   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRMPI(ierr);
984ffc4695bSBarry Smith   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRMPI(ierr);
9858c764dc5SJose Roman #endif
986ffc4695bSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRMPI(ierr);
987ffc4695bSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRMPI(ierr);
988570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
989ffc4695bSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRMPI(ierr);
990ffc4695bSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRMPI(ierr);
991ffc4695bSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRMPI(ierr);
992ffc4695bSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRMPI(ierr);
993c90a1750SBarry Smith #endif
994c90a1750SBarry Smith 
995*de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
996ffc4695bSBarry Smith   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRMPI(ierr);
997cca4cb22SSatish Balay #endif
998cca4cb22SSatish Balay 
999ffc4695bSBarry Smith   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRMPI(ierr);
1000ffc4695bSBarry Smith   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRMPI(ierr);
1001e5c89e4eSSatish Balay 
100240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1003ffc4695bSBarry Smith   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRMPI(ierr);
1004ffc4695bSBarry Smith   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRMPI(ierr);
100544041f26SJed Brown #endif
1006e5c89e4eSSatish Balay 
1007e5c89e4eSSatish Balay   /*
1008480cf27aSJed Brown      Attributes to be set on PETSc communicators
1009480cf27aSJed Brown   */
1010ffc4695bSBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);CHKERRMPI(ierr);
1011ffc4695bSBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);CHKERRMPI(ierr);
1012ffc4695bSBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);CHKERRMPI(ierr);
1013ffc4695bSBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);CHKERRMPI(ierr);
1014480cf27aSJed Brown 
1015480cf27aSJed Brown   /*
1016e8fb0fc0SBarry Smith      Build the options database
1017e5c89e4eSSatish Balay   */
1018c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
1019e5c89e4eSSatish Balay 
1020703f3eceSBarry Smith   /* call a second time so it can look in the options database */
1021703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
10226dc8fec2Sbcordonn 
1023e5c89e4eSSatish Balay   /*
102457171095SVaclav Hapla      Check system options and print help
1025e5c89e4eSSatish Balay   */
102657171095SVaclav Hapla   ierr = PetscOptionsCheckInitial_Private(help);CHKERRQ(ierr);
1027e5c89e4eSSatish Balay 
1028d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
1029d45a07a7SBarry Smith 
1030e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
103111525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1032f4202a44SBarry Smith #endif
1033f4202a44SBarry Smith 
1034e5c89e4eSSatish Balay   /*
1035e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1036e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1037e5c89e4eSSatish Balay   */
1038e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1039e5c89e4eSSatish Balay 
1040ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
104102c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1042e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
104302c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1044cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1045101491d0SBarry Smith   {
1046101491d0SBarry Smith     PetscBool omp_view_flag;
1047cd1b99a6SBarry Smith     char      *threads = getenv("OMP_NUM_THREADS");
1048e5c89e4eSSatish Balay 
1049cd1b99a6SBarry Smith    if (threads) {
105002c9f0b5SLisandro Dalcin      ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
1051f90b075cSBarry Smith      (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
1052cd1b99a6SBarry Smith    } else {
1053cd1b99a6SBarry Smith #define NMAX  10000
1054101491d0SBarry Smith      int          i;
1055cd1b99a6SBarry Smith       PetscScalar *x;
1056cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1057cd1b99a6SBarry Smith #pragma omp parallel for
1058cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1059cd1b99a6SBarry Smith         x[i] = 0.0;
1060f90b075cSBarry Smith         PetscNumOMPThreads  = (PetscInt) omp_get_num_threads();
1061cd1b99a6SBarry Smith       }
1062cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
106302c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1064101491d0SBarry Smith     }
1065101491d0SBarry Smith     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
1066f90b075cSBarry 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);
1067101491d0SBarry Smith     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
1068101491d0SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1069101491d0SBarry Smith     if (flg) {
107002c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1071f90b075cSBarry Smith       omp_set_num_threads((int)PetscNumOMPThreads);
1072101491d0SBarry Smith     }
1073101491d0SBarry Smith     if (omp_view_flag) {
1074f90b075cSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
1075cd1b99a6SBarry Smith     }
1076cd1b99a6SBarry Smith   }
1077cd1b99a6SBarry Smith #endif
1078ef6c6fedSBoyana Norris 
1079951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1080e39fd77fSBarry Smith   /*
1081e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1082e39fd77fSBarry Smith 
1083e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1084e39fd77fSBarry Smith   */
108530815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
1086ffc4695bSBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRMPI(ierr);
108730815ce0SLisandro Dalcin   }
1088951e3c8eSBarry Smith #endif
1089e39fd77fSBarry Smith 
109041c0b4b3SShri Abhyankar   /*
109141c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
109241c0b4b3SShri Abhyankar   */
1093ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1094e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1095e1167bb9SShri Abhyankar #endif
1096e1167bb9SShri Abhyankar 
10972d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10982d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10992d53ad75SBarry Smith #endif
11002d53ad75SBarry Smith 
11015e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
110287c3beb6SLisandro Dalcin   {
110387c3beb6SLisandro Dalcin     PetscViewer viewer;
110422e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
11055e71baefSBarry Smith     if (flg) {
11065e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
110787c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
110887c3beb6SLisandro Dalcin     }
11095e71baefSBarry Smith   }
11105e71baefSBarry Smith #endif
1111dff31646SBarry Smith 
111287c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
111387c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
111487c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);}
111587c3beb6SLisandro Dalcin 
1116a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1117a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
11187b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1119a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
11209fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1121a56f64adSBarry Smith #endif
112255d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
112355d657eeSBarry Smith #endif
1124a56f64adSBarry Smith 
1125301d30feSBarry Smith   /*
112657171095SVaclav Hapla       Set flag that we are completely initialized
1127301d30feSBarry Smith   */
1128301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
11292db0e300SLisandro Dalcin 
11302db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
11312db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1132301d30feSBarry Smith   PetscFunctionReturn(0);
1133e5c89e4eSSatish Balay }
1134e5c89e4eSSatish Balay 
11354097062eSBarry Smith #if defined(PETSC_USE_LOG)
113695c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1137ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1138ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
113995c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11404097062eSBarry Smith #endif
1141e5c89e4eSSatish Balay 
1142008a6e76SBarry Smith /*
1143008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1144008a6e76SBarry Smith */
1145008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1146008a6e76SBarry Smith {
1147008a6e76SBarry Smith   PetscErrorCode ierr;
1148008a6e76SBarry Smith 
1149008a6e76SBarry Smith   PetscFunctionBegin;
1150008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1151ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRMPI(ierr);
1152008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1153ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRMPI(ierr);
1154008a6e76SBarry Smith #endif
1155ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRMPI(ierr);
1156ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRMPI(ierr);
1157008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1158ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRMPI(ierr);
1159ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRMPI(ierr);
1160ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRMPI(ierr);
1161008a6e76SBarry Smith #endif
1162008a6e76SBarry Smith 
1163*de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1164ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRMPI(ierr);
1165008a6e76SBarry Smith #endif
1166008a6e76SBarry Smith 
1167ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRMPI(ierr);
116840df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1169ffc4695bSBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRMPI(ierr);
1170008a6e76SBarry Smith #endif
1171ffc4695bSBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRMPI(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 
1216ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(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;
1321ffc4695bSBarry Smith     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRMPI(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);
14079245e749SBarry 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);
1449ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(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);
1453ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(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 
1509ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(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);
1513ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(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 
1527ffc4695bSBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRMPI(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);
1531ffc4695bSBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRMPI(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 
1547c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
1548c2b86a48SJunchao Zhang   if (PetscBeganKokkos) {
1549c2b86a48SJunchao Zhang     ierr = PetscKokkosFinalize_Private();CHKERRQ(ierr);
1550c2b86a48SJunchao Zhang     PetscBeganKokkos = PETSC_FALSE;
155145639126SStefano Zampini     PetscKokkosInitialized = PETSC_FALSE;
1552c2b86a48SJunchao Zhang   }
1553c2b86a48SJunchao Zhang #endif
1554c2b86a48SJunchao Zhang 
1555008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1556e5c89e4eSSatish Balay 
1557dbc8283eSBarry Smith   /*
1558efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1559efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1560efb80d3cSBarry Smith 
1561efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1562efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1563dbc8283eSBarry Smith  */
1564b770b1f6SSatish Balay   {
1565dbc8283eSBarry Smith     PetscCommCounter *counter;
1566dbc8283eSBarry Smith     PetscMPIInt      flg;
1567dbc8283eSBarry Smith     MPI_Comm         icomm;
1568265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
1569ffc4695bSBarry Smith     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRMPI(ierr);
1570dbc8283eSBarry Smith     if (flg) {
1571265f3f35SJed Brown       icomm = ucomm.comm;
1572ffc4695bSBarry Smith       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRMPI(ierr);
1573dbc8283eSBarry 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");
1574dbc8283eSBarry Smith 
1575ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1576ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRMPI(ierr);
1577ffc4695bSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRMPI(ierr);
1578dbc8283eSBarry Smith     }
1579ffc4695bSBarry Smith     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRMPI(ierr);
1580dbc8283eSBarry Smith     if (flg) {
1581265f3f35SJed Brown       icomm = ucomm.comm;
1582ffc4695bSBarry Smith       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRMPI(ierr);
1583dbc8283eSBarry 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");
1584dbc8283eSBarry Smith 
1585ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1586ffc4695bSBarry Smith       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRMPI(ierr);
1587ffc4695bSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRMPI(ierr);
1588dbc8283eSBarry Smith     }
1589b770b1f6SSatish Balay   }
1590dbc8283eSBarry Smith 
1591ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRMPI(ierr);
1592ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRMPI(ierr);
1593ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRMPI(ierr);
1594ffc4695bSBarry Smith   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRMPI(ierr);
1595480cf27aSJed Brown 
15965ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15975ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15985ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1599ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1600ef19f930SBarry Smith 
1601e5c89e4eSSatish Balay   if (PetscBeganMPI) {
160299608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
160399b1327fSBarry Smith     PetscMPIInt flag;
1604ffc4695bSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRMPI(ierr);
1605e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
160699608316SBarry Smith #endif
1607ffc4695bSBarry Smith     ierr = MPI_Finalize();CHKERRMPI(ierr);
1608e5c89e4eSSatish Balay   }
1609e5c89e4eSSatish Balay /*
1610e5c89e4eSSatish Balay 
1611e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1612e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1613e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1614e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1615e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16160e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1617e5c89e4eSSatish Balay    memory was not freed.
1618e5c89e4eSSatish Balay 
1619e5c89e4eSSatish Balay */
16201d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1621a297a907SKarl Rupp 
16228ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_FALSE;
1623e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1624e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
162556883f60SBarry Smith #if defined(PETSC_USE_GCOV)
162656883f60SBarry Smith   /*
162756883f60SBarry Smith      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
162856883f60SBarry Smith      gcov files are still being added to the directories as git tries to remove the directories.
162956883f60SBarry Smith    */
163056883f60SBarry Smith   __gcov_flush();
163156883f60SBarry Smith #endif
16323db9a53dSBarry Smith   PetscFunctionReturn(0);
1633e5c89e4eSSatish Balay }
1634e5c89e4eSSatish Balay 
163543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16368cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
163743db4dbbSBarry Smith {
163843db4dbbSBarry Smith   if (*a == *b) return 1;
163943db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
164043db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
164143db4dbbSBarry Smith   return 0;
164243db4dbbSBarry Smith }
1643a70650f6SBarry Smith #endif
164443db4dbbSBarry Smith 
164543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16468cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
164743db4dbbSBarry Smith {
164843db4dbbSBarry Smith   if (*a == *b) return 1;
164943db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
165043db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
165143db4dbbSBarry Smith   return 0;
165243db4dbbSBarry Smith }
165343db4dbbSBarry Smith #endif
1654