xref: /petsc/src/sys/objects/pinit.c (revision 40df0d723ea94fee6b28ca5943ca31376c9d8193)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>        /*I  "petscsys.h"   I*/
6022afb99SBarry Smith #include <petscvalgrind.h>
7665c2dedSJed Brown #include <petscviewer.h>
88101f56cSMatthew Knepley 
9a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
1195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
12a9f03627SSatish Balay #endif
13f2d66bcaSShri Abhyankar 
142d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1595c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
162d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
172d53ad75SBarry Smith #endif
182d53ad75SBarry Smith 
19a6790183SBarry Smith #if defined(PETSC_HAVE_SAWS)
2016ad0300SBarry Smith #include <petscviewersaws.h>
21a6790183SBarry Smith #endif
22e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
23e5c89e4eSSatish Balay 
2495c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
25e5c89e4eSSatish Balay 
2695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
2795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
2895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
2995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
3195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
320069ddf5SShri Abhyankar 
33e5c89e4eSSatish Balay /* user may set this BEFORE calling PetscInitialize() */
34e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
35e5c89e4eSSatish Balay 
36480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
37480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
38480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
395f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
40480cf27aSJed Brown 
41e5c89e4eSSatish Balay /*
42e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
43e5c89e4eSSatish Balay */
446a6fc655SJed Brown const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",0};
456a6fc655SJed Brown const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",0};
46e5c89e4eSSatish Balay 
47ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
48ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
490f8e0872SSatish Balay 
50a2f94806SJed Brown PetscInt PetscHotRegionDepth;
51a2f94806SJed Brown 
52b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
53b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
54b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
55b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
56b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
57b22622e2STadeu Manoel #endif
58b22622e2STadeu Manoel 
59e5c89e4eSSatish Balay /*
60e5c89e4eSSatish Balay        Checks the options database for initializations related to the
61e5c89e4eSSatish Balay     PETSc components
62e5c89e4eSSatish Balay */
6395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode  PetscOptionsCheckInitial_Components(void)
64e5c89e4eSSatish Balay {
652d747510SLisandro Dalcin   PetscBool      flg;
66e5c89e4eSSatish Balay   PetscErrorCode ierr;
67e5c89e4eSSatish Balay 
68e5c89e4eSSatish Balay   PetscFunctionBegin;
692d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
702d747510SLisandro Dalcin   if (flg) {
71e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
72e8e7597cSSatish Balay     MPI_Comm comm = PETSC_COMM_WORLD;
73e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"------Additional PETSc component options--------\n");CHKERRQ(ierr);
748e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -log_exclude: <vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
758e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -info_exclude: <null,vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
76e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
77e5c89e4eSSatish Balay #endif
78e5c89e4eSSatish Balay   }
79e5c89e4eSSatish Balay   PetscFunctionReturn(0);
80e5c89e4eSSatish Balay }
81e5c89e4eSSatish Balay 
820f11a792SBarry Smith /*
83945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8472a42c3cSBarry Smith 
8572a42c3cSBarry Smith    Collective
8672a42c3cSBarry Smith 
8772a42c3cSBarry Smith    Level: advanced
8872a42c3cSBarry Smith 
8995452b02SPatrick Sanan     Notes:
90a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
910f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
92a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
930f11a792SBarry Smith 
94a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
951ea5a559SBarry Smith 
9672a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
970f11a792SBarry Smith */
98945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
9972a42c3cSBarry Smith {
10072a42c3cSBarry Smith   PetscErrorCode ierr;
10172a42c3cSBarry Smith   int            myargc   = argc;
10272a42c3cSBarry Smith   char           **myargs = args;
10372a42c3cSBarry Smith 
10472a42c3cSBarry Smith   PetscFunctionBegin;
1053bf036e2SBarry Smith   ierr = PetscInitialize(&myargc,&myargs,filename,help);CHKERRQ(ierr);
1061ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
107df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
10872a42c3cSBarry Smith   PetscFunctionReturn(ierr);
10972a42c3cSBarry Smith }
11072a42c3cSBarry Smith 
111f0865b08SBarry Smith /*
112a56f64adSBarry Smith       Used by Julia interface to get communicator
113f0865b08SBarry Smith */
114945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
115f0865b08SBarry Smith {
116f0865b08SBarry Smith   PetscFunctionBegin;
117f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
118f0865b08SBarry Smith   PetscFunctionReturn(0);
119f0865b08SBarry Smith }
120f0865b08SBarry Smith 
121e5c89e4eSSatish Balay /*@C
122e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
123e5c89e4eSSatish Balay         the command line arguments.
124e5c89e4eSSatish Balay 
125e5c89e4eSSatish Balay    Collective
126e5c89e4eSSatish Balay 
127e5c89e4eSSatish Balay    Level: advanced
128e5c89e4eSSatish Balay 
129e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
130e5c89e4eSSatish Balay @*/
1317087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
132e5c89e4eSSatish Balay {
133e5c89e4eSSatish Balay   PetscErrorCode ierr;
134e5c89e4eSSatish Balay   int            argc   = 0;
135e5c89e4eSSatish Balay   char           **args = 0;
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay   PetscFunctionBegin;
1380298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
139e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
140e5c89e4eSSatish Balay }
141e5c89e4eSSatish Balay 
142e5c89e4eSSatish Balay /*@
143e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
144e5c89e4eSSatish Balay 
14593b6d2d1SJed Brown    Level: beginner
146e5c89e4eSSatish Balay 
147e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
148e5c89e4eSSatish Balay @*/
1497087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool  *isInitialized)
150e5c89e4eSSatish Balay {
151e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
15293b6d2d1SJed Brown   return 0;
153e5c89e4eSSatish Balay }
154e5c89e4eSSatish Balay 
155e5c89e4eSSatish Balay /*@
156e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
157e5c89e4eSSatish Balay 
158e5c89e4eSSatish Balay    Level: developer
159e5c89e4eSSatish Balay 
160e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
161e5c89e4eSSatish Balay @*/
1627087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
163e5c89e4eSSatish Balay {
164e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
16593b6d2d1SJed Brown   return 0;
166e5c89e4eSSatish Balay }
167e5c89e4eSSatish Balay 
16895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void);
169e5c89e4eSSatish Balay 
170e5c89e4eSSatish Balay /*
171e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
172e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
173e5c89e4eSSatish Balay */
174367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
175e5c89e4eSSatish Balay 
176367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
177e5c89e4eSSatish Balay {
178e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
179e5c89e4eSSatish Balay 
180e5c89e4eSSatish Balay   PetscFunctionBegin;
181e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
182e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
183e5c89e4eSSatish Balay     MPI_Abort(MPI_COMM_WORLD,1);
184e5c89e4eSSatish Balay   }
185e5c89e4eSSatish Balay 
186e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
187e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
188e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
189e5c89e4eSSatish Balay   }
190812af9f3SBarry Smith   PetscFunctionReturnVoid();
191e5c89e4eSSatish Balay }
192e5c89e4eSSatish Balay 
193e5c89e4eSSatish Balay /*
194e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
195e5c89e4eSSatish Balay sum of the second entry.
196b693b147SBarry Smith 
19776f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
198367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
199b693b147SBarry Smith there would be no place to store the both needed results.
200e5c89e4eSSatish Balay */
20176ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
202e5c89e4eSSatish Balay {
203e5c89e4eSSatish Balay   PetscErrorCode ierr;
204e5c89e4eSSatish Balay 
205e5c89e4eSSatish Balay   PetscFunctionBegin;
206d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
207d6e4c47cSJed Brown   {
208d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
209367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
210d6e4c47cSJed Brown     *max = work.max;
211d6e4c47cSJed Brown     *sum = work.sum;
212d6e4c47cSJed Brown   }
213d6e4c47cSJed Brown #else
214d6e4c47cSJed Brown   {
215d6e4c47cSJed Brown     PetscMPIInt    size,rank;
216d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
217e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
218e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
219785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
220367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2216ac3741eSJed Brown     *max = work[rank].max;
2226ac3741eSJed Brown     *sum = work[rank].sum;
223e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
224d6e4c47cSJed Brown   }
225d6e4c47cSJed Brown #endif
226e5c89e4eSSatish Balay   PetscFunctionReturn(0);
227e5c89e4eSSatish Balay }
228e5c89e4eSSatish Balay 
229e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
230e5c89e4eSSatish Balay 
231570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
23206a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
233e5c89e4eSSatish Balay 
2348cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
235e5c89e4eSSatish Balay {
236e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
237e5c89e4eSSatish Balay 
238e5c89e4eSSatish Balay   PetscFunctionBegin;
2397c2de775SJed Brown   if (*datatype == MPIU_REAL) {
240e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
241a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2427c2de775SJed Brown   }
2437c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2447c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2457c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
246a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2477c2de775SJed Brown   }
2487c2de775SJed Brown #endif
2497c2de775SJed Brown   else {
2507c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
251e2e03761SBarry Smith     MPI_Abort(MPI_COMM_WORLD,1);
252e2e03761SBarry Smith   }
253812af9f3SBarry Smith   PetscFunctionReturnVoid();
254e5c89e4eSSatish Balay }
255e5c89e4eSSatish Balay #endif
256e5c89e4eSSatish Balay 
257570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
258d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
259d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
260d9822059SBarry Smith 
2618cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
262d9822059SBarry Smith {
263d9822059SBarry Smith   PetscInt i,count = *cnt;
264d9822059SBarry Smith 
265d9822059SBarry Smith   PetscFunctionBegin;
2667c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2678c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
268a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2697c2de775SJed Brown   }
2707c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2717c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2727c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2737c2de775SJed Brown     for (i=0; i<count; i++) {
2747c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2757c2de775SJed Brown     }
2767c2de775SJed Brown   }
2777c2de775SJed Brown #endif
2787c2de775SJed Brown   else {
2797c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
2808c764dc5SJose Roman     MPI_Abort(MPI_COMM_WORLD,1);
2818c764dc5SJose Roman   }
282d9822059SBarry Smith   PetscFunctionReturnVoid();
283d9822059SBarry Smith }
284d9822059SBarry Smith 
2858cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
286d9822059SBarry Smith {
287d9822059SBarry Smith   PetscInt    i,count = *cnt;
288d9822059SBarry Smith 
289d9822059SBarry Smith   PetscFunctionBegin;
2907c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2918c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
292a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2937c2de775SJed Brown   }
2947c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2957c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2967c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2977c2de775SJed Brown     for (i=0; i<count; i++) {
2987c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2997c2de775SJed Brown     }
3007c2de775SJed Brown   }
3017c2de775SJed Brown #endif
3027c2de775SJed Brown   else {
3038c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
3048c764dc5SJose Roman     MPI_Abort(MPI_COMM_WORLD,1);
3058c764dc5SJose Roman   }
306d9822059SBarry Smith   PetscFunctionReturnVoid();
307d9822059SBarry Smith }
308d9822059SBarry Smith #endif
309d9822059SBarry Smith 
310480cf27aSJed Brown /*
311480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
312480cf27aSJed Brown 
313ff0e51ddSBarry 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.
314480cf27aSJed Brown 
31512801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
316480cf27aSJed Brown 
317480cf27aSJed Brown */
3188cc058d9SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
319480cf27aSJed Brown {
320480cf27aSJed Brown   PetscErrorCode ierr;
321480cf27aSJed Brown 
322480cf27aSJed Brown   PetscFunctionBegin;
32312801b39SBarry Smith   ierr = PetscInfo1(0,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
32412801b39SBarry Smith   ierr = PetscFree(count_val);CHKERRMPI(ierr);
325480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
326480cf27aSJed Brown }
327480cf27aSJed Brown 
328480cf27aSJed Brown /*
32947435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
330da3039f7SJed Brown   calls MPI_Comm_free().
331da3039f7SJed Brown 
332da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
333480cf27aSJed Brown 
334ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
335480cf27aSJed Brown 
33612801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
337480cf27aSJed Brown 
338480cf27aSJed Brown */
339da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
340480cf27aSJed Brown {
341480cf27aSJed Brown   PetscErrorCode                    ierr;
342b89831e5SBarry Smith   PetscMPIInt                       flg;
343265f3f35SJed Brown   union {MPI_Comm comm; void *ptr;} icomm,ocomm;
344480cf27aSJed Brown 
345480cf27aSJed Brown   PetscFunctionBegin;
34612801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
347ec4fadc2SJed Brown   icomm.ptr = attr_val;
348da3039f7SJed Brown 
34947435625SJed Brown   ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
35012801b39SBarry Smith   if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm");
35112801b39SBarry Smith   if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm");
35247435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
35312801b39SBarry Smith   ierr = PetscInfo1(0,"User MPI_Comm %ld is being freed after removing reference from inner PETSc comm to this outer comm\n",(long)comm);CHKERRMPI(ierr);
354da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
355b89831e5SBarry Smith }
356da3039f7SJed Brown 
357da3039f7SJed Brown /*
35847435625SJed Brown  * This is invoked on the inner comm when Petsc_DelComm_Outer calls MPI_Comm_delete_attr.  It should not be reached any other way.
359da3039f7SJed Brown  */
360da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
361da3039f7SJed Brown {
362da3039f7SJed Brown   PetscErrorCode ierr;
363da3039f7SJed Brown 
364da3039f7SJed Brown   PetscFunctionBegin;
36512801b39SBarry Smith   ierr = PetscInfo1(0,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
366480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
367480cf27aSJed Brown }
368480cf27aSJed Brown 
3695f7487a0SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Shm(MPI_Comm,PetscMPIInt,void *,void *);
37042218b76SBarry Smith 
371951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
372e39fd77fSBarry Smith #if !defined(PETSC_WORDS_BIGENDIAN)
3738cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3748cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3758cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
376e39fd77fSBarry Smith #endif
377951e3c8eSBarry Smith #endif
378e39fd77fSBarry Smith 
37912801b39SBarry Smith PetscMPIInt PETSC_MPI_ERROR_CLASS,PETSC_MPI_ERROR_CODE;
38012801b39SBarry Smith 
381eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
382eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3836ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
3846ae9a8a6SBarry Smith char **PetscGlobalArgs = 0;
385dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
386e5c89e4eSSatish Balay 
387dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
388051e4cf2SJed Brown {
389051e4cf2SJed Brown   PetscErrorCode ierr;
390051e4cf2SJed Brown 
391051e4cf2SJed Brown   PetscFunctionBegin;
392051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
393f191ac77SSatish Balay   ierr = PetscCitationsRegister("@TechReport{petsc-user-ref,\n  Author = {Satish Balay and Shrirang Abhyankar and Mark F. Adams and Jed Brown and Peter Brune\n            and Kris Buschelman and Lisandro Dalcin and Victor Eijkhout and William D. Gropp\n            and 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.9},\n  Institution = {Argonne National Laboratory},\n  Year = {2018}\n}\n",NULL);CHKERRQ(ierr);
394051e4cf2SJed 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);
395051e4cf2SJed Brown   PetscFunctionReturn(0);
396051e4cf2SJed Brown }
397e5c89e4eSSatish Balay 
3982d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3992d747510SLisandro Dalcin 
4002d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4012d747510SLisandro Dalcin {
4022d747510SLisandro Dalcin   PetscErrorCode ierr;
4032d747510SLisandro Dalcin 
4042d747510SLisandro Dalcin   PetscFunctionBegin;
4052d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
4062d747510SLisandro Dalcin   PetscFunctionReturn(0);
4072d747510SLisandro Dalcin }
4082d747510SLisandro Dalcin 
4092d747510SLisandro Dalcin /*@C
4102d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4112d747510SLisandro Dalcin 
4122d747510SLisandro Dalcin     Not Collective
4132d747510SLisandro Dalcin 
4142d747510SLisandro Dalcin     Input Parameter:
4152d747510SLisandro Dalcin .   len - length of the string name
4162d747510SLisandro Dalcin 
4172d747510SLisandro Dalcin     Output Parameter:
4182d747510SLisandro Dalcin .   name - the name of the running program
4192d747510SLisandro Dalcin 
4202d747510SLisandro Dalcin    Level: advanced
4212d747510SLisandro Dalcin 
4222d747510SLisandro Dalcin     Notes:
4232d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4242d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4252d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4262d747510SLisandro Dalcin @*/
4272d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4282d747510SLisandro Dalcin {
4292d747510SLisandro Dalcin   PetscErrorCode ierr;
4302d747510SLisandro Dalcin 
4312d747510SLisandro Dalcin   PetscFunctionBegin;
4322d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4332d747510SLisandro Dalcin   PetscFunctionReturn(0);
4342d747510SLisandro Dalcin }
4352d747510SLisandro Dalcin 
436e5c89e4eSSatish Balay /*@C
437e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
438e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
439e5c89e4eSSatish Balay 
440e5c89e4eSSatish Balay    Not Collective
441e5c89e4eSSatish Balay 
442e5c89e4eSSatish Balay    Output Parameters:
443e5c89e4eSSatish Balay +  argc - count of number of command line arguments
444e5c89e4eSSatish Balay -  args - the command line arguments
445e5c89e4eSSatish Balay 
446e5c89e4eSSatish Balay    Level: intermediate
447e5c89e4eSSatish Balay 
448e5c89e4eSSatish Balay    Notes:
449e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
450e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
451e5c89e4eSSatish Balay 
452f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
453f177e3b1SBarry Smith 
454e5c89e4eSSatish Balay    Concepts: command line arguments
455e5c89e4eSSatish Balay 
456793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
457e5c89e4eSSatish Balay 
458e5c89e4eSSatish Balay @*/
4597087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
460e5c89e4eSSatish Balay {
461e5c89e4eSSatish Balay   PetscFunctionBegin;
46217186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
463e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
464e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
465e5c89e4eSSatish Balay   PetscFunctionReturn(0);
466e5c89e4eSSatish Balay }
467e5c89e4eSSatish Balay 
468793721a6SBarry Smith /*@C
469793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
470793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
471793721a6SBarry Smith 
472793721a6SBarry Smith    Not Collective
473793721a6SBarry Smith 
474793721a6SBarry Smith    Output Parameters:
475793721a6SBarry Smith .  args - the command line arguments
476793721a6SBarry Smith 
477793721a6SBarry Smith    Level: intermediate
478793721a6SBarry Smith 
479793721a6SBarry Smith    Notes:
480793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
481793721a6SBarry Smith 
482793721a6SBarry Smith    Concepts: command line arguments
483793721a6SBarry Smith 
484793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
485793721a6SBarry Smith 
486793721a6SBarry Smith @*/
4877087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
488793721a6SBarry Smith {
489793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
490793721a6SBarry Smith   PetscErrorCode ierr;
491793721a6SBarry Smith 
492793721a6SBarry Smith   PetscFunctionBegin;
49317186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4942d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
495785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
496793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
497793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
498793721a6SBarry Smith   }
4992d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
500793721a6SBarry Smith   PetscFunctionReturn(0);
501793721a6SBarry Smith }
502793721a6SBarry Smith 
503793721a6SBarry Smith /*@C
504793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
505793721a6SBarry Smith 
506793721a6SBarry Smith    Not Collective
507793721a6SBarry Smith 
508793721a6SBarry Smith    Output Parameters:
509793721a6SBarry Smith .  args - the command line arguments
510793721a6SBarry Smith 
511793721a6SBarry Smith    Level: intermediate
512793721a6SBarry Smith 
513793721a6SBarry Smith    Concepts: command line arguments
514793721a6SBarry Smith 
515793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
516793721a6SBarry Smith 
517793721a6SBarry Smith @*/
5187087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
519793721a6SBarry Smith {
520793721a6SBarry Smith   PetscInt       i = 0;
521793721a6SBarry Smith   PetscErrorCode ierr;
522793721a6SBarry Smith 
523793721a6SBarry Smith   PetscFunctionBegin;
524a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
525793721a6SBarry Smith   while (args[i]) {
526793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
527793721a6SBarry Smith     i++;
528793721a6SBarry Smith   }
529793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
530793721a6SBarry Smith   PetscFunctionReturn(0);
531793721a6SBarry Smith }
532793721a6SBarry Smith 
53311525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
53430befbd2SBarry Smith #include <petscconfiginfo.h>
53530befbd2SBarry Smith 
53695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
53711525c0dSBarry Smith {
53811525c0dSBarry Smith   if (!PetscGlobalRank) {
53930befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
54011525c0dSBarry Smith     int            port;
541ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
54211525c0dSBarry Smith     size_t         applinelen,introlen;
54311525c0dSBarry Smith     PetscErrorCode ierr;
544ffbd1cfbSBarry Smith     char           sawsurl[256];
54511525c0dSBarry Smith 
546c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
54711525c0dSBarry Smith     if (flg) {
54811525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
54911525c0dSBarry Smith 
550c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
55111525c0dSBarry Smith       if (sawslog[0]) {
55211525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
55311525c0dSBarry Smith       } else {
55411525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
55511525c0dSBarry Smith       }
55611525c0dSBarry Smith     }
557c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
55811525c0dSBarry Smith     if (flg) {
55911525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
56011525c0dSBarry Smith     }
561c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
562ffbd1cfbSBarry Smith     if (selectport) {
563ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
564ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
565ffbd1cfbSBarry Smith     } else {
566c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
56711525c0dSBarry Smith       if (flg) {
56811525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
56911525c0dSBarry Smith       }
570ffbd1cfbSBarry Smith     }
571c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
57211525c0dSBarry Smith     if (flg) {
57311525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
57411525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5759c1e0ce8SBarry Smith     } else {
576c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5779c1e0ce8SBarry Smith       if (flg) {
5783c01dfcfSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
5799c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5809c1e0ce8SBarry Smith       }
58111525c0dSBarry Smith     }
582c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
58311525c0dSBarry Smith     if (flg2) {
58411525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
58511525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
58611525c0dSBarry Smith       ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr);
58711525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
58811525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
58943da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
59011525c0dSBarry Smith     }
59111525c0dSBarry Smith     ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr);
59211525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
59311525c0dSBarry Smith     introlen   = 4096 + applinelen;
59430a8c9c0SSurtai Han     applinelen += 1024;
59511525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
59611525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
59711525c0dSBarry Smith 
59811525c0dSBarry Smith     if (rootlocal) {
59911525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
60011525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
60111525c0dSBarry Smith     }
60276a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
60311525c0dSBarry Smith     if (rootlocal && help) {
604928bb9adSStefano 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);
60511525c0dSBarry Smith     } else if (help) {
606928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
60711525c0dSBarry Smith     } else {
608928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
60911525c0dSBarry Smith     }
610b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
61130befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
61211525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
61311525c0dSBarry Smith                                     "<center><h2> <a href=\"http://www.mcs.anl.gov/petsc\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n"
614df62c222SBarry 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"
615928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
61643da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
61711525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
61811525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
619ffbd1cfbSBarry Smith     if (selectport) {
620aa573868SBarry Smith       PetscBool silent;
6217d812c46SBarry Smith 
6227d812c46SBarry Smith       ierr = SAWs_Initialize();
6237d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6247d812c46SBarry Smith       while (ierr) {
6257d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6267d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6277d812c46SBarry Smith         ierr = SAWs_Initialize();
6287d812c46SBarry Smith       }
6297d812c46SBarry Smith 
630aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
631aa573868SBarry Smith       if (!silent) {
632ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
633ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
634ffbd1cfbSBarry Smith       }
6357d812c46SBarry Smith     } else {
6367d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
637aa573868SBarry Smith     }
6380af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6390af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6400af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6410af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6420af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
64311525c0dSBarry Smith   }
64411525c0dSBarry Smith   PetscFunctionReturn(0);
64511525c0dSBarry Smith }
64611525c0dSBarry Smith #endif
64711525c0dSBarry Smith 
648a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
649a56f64adSBarry Smith #include <adios.h>
65022580e64SBarry Smith #include <adios_read.h>
6517b56e58cSSatish Balay int64_t Petsc_adios_group;
652a56f64adSBarry Smith #endif
65355d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
65455d657eeSBarry Smith #include <adios2_c.h>
65555d657eeSBarry Smith #endif
656a56f64adSBarry Smith 
657e5c89e4eSSatish Balay /*@C
658e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
659e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
660e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
661e5c89e4eSSatish Balay    your program -- usually the very first line!
662e5c89e4eSSatish Balay 
663e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
664e5c89e4eSSatish Balay 
665e5c89e4eSSatish Balay    Input Parameters:
666e5c89e4eSSatish Balay +  argc - count of number of command line arguments
667e5c89e4eSSatish Balay .  args - the command line arguments
6680298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
669fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6700298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
671e5c89e4eSSatish Balay 
67205827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
67305827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
67405827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
67505827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
67605827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
677e5c89e4eSSatish Balay 
678e5c89e4eSSatish Balay    Options Database Keys:
6797ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6807ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
681e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
6828a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
6838a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
6848a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
6858a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
6868a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
687e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
688e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
689e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
690e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
6912fb0ec9aSBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries)
692e5c89e4eSSatish Balay .  -malloc no - Indicates not to use error-checking malloc
6932fb0ec9aSBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free
694aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
695dc92acbaSJed Brown .  -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds
696e5c89e4eSSatish Balay .  -fp_trap - Stops on floating point exceptions (Note that on the
697e5c89e4eSSatish Balay               IBM RS6000 this slows code by at least a factor of 10.)
698e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
699e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
700e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
701e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
702e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7030841954dSBarry Smith -  -memory_view - Print memory usage at end of run
704e5c89e4eSSatish Balay 
705e5c89e4eSSatish Balay    Options Database Keys for Profiling:
706a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
707495fc317SBarry Smith +  -info <optional filename> - Prints verbose information to the screen
708495fc317SBarry Smith .  -info_exclude <null,vec,mat,pc,ksp,snes,ts> - Excludes some of the verbose messages
709217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
710217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
711495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
712e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7139a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
7149a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
715495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
716f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
717495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
718495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
719c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
72087c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
721c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
722495fc317SBarry Smith 
723609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
724e5c89e4eSSatish Balay 
725ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
726ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
727ffbd1cfbSBarry 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
728ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
729ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
730ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
731ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
732ffbd1cfbSBarry Smith 
733e5c89e4eSSatish Balay    Environmental Variables:
734e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
735e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
736e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
737e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
738e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
739e5c89e4eSSatish Balay 
740e5c89e4eSSatish Balay 
741e5c89e4eSSatish Balay    Level: beginner
742e5c89e4eSSatish Balay 
743e5c89e4eSSatish Balay    Notes:
744e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
745e5c89e4eSSatish Balay    it before PetscInitialize().
746e5c89e4eSSatish Balay 
747e5c89e4eSSatish Balay    Fortran Version:
748e5c89e4eSSatish Balay    In Fortran this routine has the format
749e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
750e5c89e4eSSatish Balay 
751e5c89e4eSSatish Balay +   ierr - error return code
7520eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
753fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
754e5c89e4eSSatish Balay 
755e5c89e4eSSatish Balay    Important Fortran Note:
7560eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7570298fd71SBarry Smith    null character string; you CANNOT just use NULL as
758a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
759e5c89e4eSSatish Balay 
76001cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
76101cb0274SBarry Smith    calling PetscInitialize().
762e5c89e4eSSatish Balay 
763e5c89e4eSSatish Balay    Concepts: initializing PETSc
764e5c89e4eSSatish Balay 
76501cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
766e5c89e4eSSatish Balay 
767e5c89e4eSSatish Balay @*/
7687087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
769e5c89e4eSSatish Balay {
770e5c89e4eSSatish Balay   PetscErrorCode ierr;
7714bb5149bSJed Brown   PetscMPIInt    flag, size;
77219c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
773e5c89e4eSSatish Balay   char           hostname[256];
774e5c89e4eSSatish Balay 
775e5c89e4eSSatish Balay   PetscFunctionBegin;
776e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7773d96e996SBarry Smith   /*
7783d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7793d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7803d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7813d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7823d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7833d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7843d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7853d96e996SBarry Smith       listed above and since that time are compatible.
7863d96e996SBarry Smith 
7873d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
7883d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
7893d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
7903d96e996SBarry Smith       to perform the checking.
7913d96e996SBarry Smith 
7923d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
7933d96e996SBarry Smith 
7943d96e996SBarry Smith       Questions:
7953d96e996SBarry Smith 
7963d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
7973d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
7983d96e996SBarry Smith   */
7993d96e996SBarry Smith 
80019c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
80119c5658aSBarry Smith   {
80219c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
80319c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
80419c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8053d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
80619c5658aSBarry Smith #if defined(MPICH_VERSION)
8073d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
80819c5658aSBarry Smith     {
80919c5658aSBarry Smith       char *ver,*lf;
81019c5658aSBarry Smith       flg = PETSC_FALSE;
81119c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
81219c5658aSBarry Smith       if (ver) {
81319c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
81419c5658aSBarry Smith         if (lf) {
81519c5658aSBarry Smith           *lf = 0;
81619c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
81719c5658aSBarry Smith         }
81819c5658aSBarry Smith       }
81919c5658aSBarry Smith       if (!flg) {
82019c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8213d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
82219c5658aSBarry Smith       }
82319c5658aSBarry Smith     }
8243d96e996SBarry Smith #endif
8253d96e996SBarry 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?) */
82619c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
82719c5658aSBarry Smith     {
82819c5658aSBarry Smith       char *ver,bs[32],*bsf;
82919c5658aSBarry Smith       flg = PETSC_FALSE;
83019c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
83119c5658aSBarry Smith       if (ver) {
8322e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
83319c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
83419c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
83519c5658aSBarry Smith       }
83619c5658aSBarry Smith       if (!flg) {
83719c5658aSBarry 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);
8383d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
83919c5658aSBarry Smith       }
84019c5658aSBarry Smith     }
84119c5658aSBarry Smith #endif
84219c5658aSBarry Smith   }
84319c5658aSBarry Smith #endif
84419c5658aSBarry Smith 
845e5c89e4eSSatish Balay 
846ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
847d89683f4Sbcordonn   PETSC_STDOUT = stdout;
848ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
849e5c89e4eSSatish Balay 
8500c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8510c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8520c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8530c30907bSSatish Balay #endif
8540c30907bSSatish Balay 
8554416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
856e5c89e4eSSatish Balay 
857e5c89e4eSSatish Balay   /*
858e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
859e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
860e5c89e4eSSatish Balay   */
861e5c89e4eSSatish Balay   if (argc && *argc) {
862e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
863e5c89e4eSSatish Balay   } else {
864e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
865e5c89e4eSSatish Balay   }
866e5c89e4eSSatish Balay 
867e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
868e5c89e4eSSatish Balay   if (!flag) {
869e32f2f54SBarry 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");
8705e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8715e765c61SJed Brown     {
8725e765c61SJed Brown       PetscMPIInt provided;
8735e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
8745e765c61SJed Brown     }
8755e765c61SJed Brown #else
876e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
8775e765c61SJed Brown #endif
878e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
879e5c89e4eSSatish Balay   }
880e5c89e4eSSatish Balay   if (argc && args) {
881e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
882e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
883e5c89e4eSSatish Balay   }
884e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
8855ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
8865ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
8875ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
888ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
889e5c89e4eSSatish Balay 
890a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
891d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
892e5c89e4eSSatish Balay 
89312801b39SBarry Smith   ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
89412801b39SBarry Smith   ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
89512801b39SBarry Smith 
896e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
897e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
898e5c89e4eSSatish Balay 
899e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
900e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
901e5c89e4eSSatish Balay 
9028ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9038ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9048ad47952SJed Brown 
905e5c89e4eSSatish Balay   /*
906e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
907e5c89e4eSSatish Balay      shared libraries the constructors for global variables
908e5c89e4eSSatish Balay      are not called; at least on IRIX.
909e5c89e4eSSatish Balay   */
910886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
911e5c89e4eSSatish Balay   {
912252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
91350f81f78SJed Brown     PetscComplex ic(0.0,1.0);
914e5c89e4eSSatish Balay     PETSC_i = ic;
915252ecd31SSatish Balay #else
91650f81f78SJed Brown     PETSC_i = _Complex_I;
917b7940d39SSatish Balay #endif
918762437b8SSatish Balay   }
919762437b8SSatish Balay 
9202c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
921e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
922500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
923500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
924500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9252c876bd9SBarry Smith #endif
926886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
927e5c89e4eSSatish Balay 
928e5c89e4eSSatish Balay   /*
929e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
930e5c89e4eSSatish Balay      half of the entries and maxes the second half.
931e5c89e4eSSatish Balay   */
932367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
933e5c89e4eSSatish Balay 
934ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
935c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
936c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9377c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9388c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9398c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9408c764dc5SJose Roman #endif
941d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
942d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
943570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
944570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
945570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
946570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
947570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
948c90a1750SBarry Smith #endif
949c90a1750SBarry Smith 
950570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
951cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
952cca4cb22SSatish Balay #endif
953cca4cb22SSatish Balay 
954e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
955e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
956e5c89e4eSSatish Balay 
957*40df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
958e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
959e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
96044041f26SJed Brown #endif
961e5c89e4eSSatish Balay 
962ec957eceSBarry Smith 
963e5c89e4eSSatish Balay   /*
964480cf27aSJed Brown      Attributes to be set on PETSc communicators
965480cf27aSJed Brown   */
96612801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
96712801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
96812801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
9695f7487a0SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Shm,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
970480cf27aSJed Brown 
971480cf27aSJed Brown   /*
972e8fb0fc0SBarry Smith      Build the options database
973e5c89e4eSSatish Balay   */
974c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
975e5c89e4eSSatish Balay 
976703f3eceSBarry Smith   /* call a second time so it can look in the options database */
977703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
9786dc8fec2Sbcordonn 
979e5c89e4eSSatish Balay   /*
980e5c89e4eSSatish Balay      Print main application help message
981e5c89e4eSSatish Balay   */
9822d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
983e5c89e4eSSatish Balay   if (help && flg) {
984e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
985e5c89e4eSSatish Balay   }
986e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
987e5c89e4eSSatish Balay 
988d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
989d45a07a7SBarry Smith 
990e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
99111525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
992f4202a44SBarry Smith #endif
993f4202a44SBarry Smith 
994896238b9SBarry Smith   /* Creates the logging data structures; this is enabled even if logging is not turned on */
995a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
996896238b9SBarry Smith   ierr = PetscLogInitialize();CHKERRQ(ierr);
997a9f03627SSatish Balay #endif
998e5c89e4eSSatish Balay 
999e5c89e4eSSatish Balay   /*
1000e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1001e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1002e5c89e4eSSatish Balay   */
1003e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1004e5c89e4eSSatish Balay 
1005e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
1006ae15b995SBarry Smith   ierr = PetscInfo1(0,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1007e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
1008ae15b995SBarry Smith   ierr = PetscInfo1(0,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1009e5c89e4eSSatish Balay 
1010e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr);
1011ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1012c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1013ef6c6fedSBoyana Norris 
1014951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1015e39fd77fSBarry Smith   /*
1016e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1017e39fd77fSBarry Smith 
1018e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1019e39fd77fSBarry Smith   */
1020e39fd77fSBarry Smith #if !defined(PETSC_WORDS_BIGENDIAN)
10210298fd71SBarry Smith   ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
1022e39fd77fSBarry Smith #endif
1023951e3c8eSBarry Smith #endif
1024e39fd77fSBarry Smith 
102541c0b4b3SShri Abhyankar   /*
102641c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
102741c0b4b3SShri Abhyankar   */
1028ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1029e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1030e1167bb9SShri Abhyankar #endif
1031e1167bb9SShri Abhyankar 
10322d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10332d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10342d53ad75SBarry Smith #endif
10352d53ad75SBarry Smith 
10365e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
103787c3beb6SLisandro Dalcin   {
103887c3beb6SLisandro Dalcin     PetscViewer viewer;
103922e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
10405e71baefSBarry Smith     if (flg) {
10415e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
104287c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
104387c3beb6SLisandro Dalcin     }
10445e71baefSBarry Smith   }
10455e71baefSBarry Smith #endif
1046dff31646SBarry Smith 
104787c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
104887c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
104987c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
105087c3beb6SLisandro Dalcin 
1051a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1052a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
10537b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1054a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
10559fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1056a56f64adSBarry Smith #endif
105755d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
105855d657eeSBarry Smith #endif
1059a56f64adSBarry Smith 
1060301d30feSBarry Smith   /*
1061301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1062301d30feSBarry Smith   */
1063301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
10642db0e300SLisandro Dalcin 
10652db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
10662db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1067301d30feSBarry Smith   PetscFunctionReturn(0);
1068e5c89e4eSSatish Balay }
1069e5c89e4eSSatish Balay 
10704097062eSBarry Smith #if defined(PETSC_USE_LOG)
107195c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1072ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1073ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
107495c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
10754097062eSBarry Smith #endif
1076e5c89e4eSSatish Balay 
1077008a6e76SBarry Smith /*
1078008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1079008a6e76SBarry Smith */
1080008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1081008a6e76SBarry Smith {
1082008a6e76SBarry Smith   PetscErrorCode ierr;
1083008a6e76SBarry Smith 
1084008a6e76SBarry Smith   PetscFunctionBegin;
1085008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1086008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1087008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1088008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1089008a6e76SBarry Smith #endif
1090008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1091008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1092008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1093008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1094008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1095008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1096008a6e76SBarry Smith #endif
1097008a6e76SBarry Smith 
1098008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1099008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1100008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1101008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1102008a6e76SBarry Smith #endif
1103008a6e76SBarry Smith #endif
1104008a6e76SBarry Smith 
1105008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1106008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1107008a6e76SBarry Smith #endif
1108008a6e76SBarry Smith 
1109008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
1110*40df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1111008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1112008a6e76SBarry Smith #endif
1113008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1114008a6e76SBarry Smith   PetscFunctionReturn(0);
1115008a6e76SBarry Smith }
1116008a6e76SBarry Smith 
1117e5c89e4eSSatish Balay /*@C
1118e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1119e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1120e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1121e5c89e4eSSatish Balay 
1122e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1123e5c89e4eSSatish Balay 
1124e5c89e4eSSatish Balay    Options Database Keys:
112526a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1126e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11277eb1d149SBarry 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
1128e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1129e5c89e4eSSatish Balay .  -malloc_dump - Calls PetscMallocDump()
1130e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
1131e5c89e4eSSatish Balay -  -malloc_log - Prints summary of memory usage
1132e5c89e4eSSatish Balay 
1133e5c89e4eSSatish Balay    Level: beginner
1134e5c89e4eSSatish Balay 
1135e5c89e4eSSatish Balay    Note:
1136e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1137e5c89e4eSSatish Balay 
113888c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1139e5c89e4eSSatish Balay @*/
11407087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1141e5c89e4eSSatish Balay {
1142e5c89e4eSSatish Balay   PetscErrorCode ierr;
11434bb5149bSJed Brown   PetscMPIInt    rank;
1144a8d2bbe5SBarry Smith   PetscInt       nopt;
11452bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1146dff31646SBarry Smith   PetscBool      flg;
114710463e74SBarry Smith #if defined(PETSC_USE_LOG)
114810463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
114910463e74SBarry Smith #endif
1150e5c89e4eSSatish Balay 
1151e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
11524b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
11533cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1154e5c89e4eSSatish Balay   }
11553cbbc5ffSBarry Smith   PetscFunctionBegin;
11560298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1157b022a5c1SBarry Smith 
11581f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1159a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
116022580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1161a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1162a56f64adSBarry Smith #endif
116355d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
116455d657eeSBarry Smith #endif
1165c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1166dff31646SBarry Smith   if (flg) {
11671f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
11681f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
11691f817a21SBarry Smith 
1170c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11711f817a21SBarry Smith     if (filename[0]) {
11721f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
11731f817a21SBarry Smith     }
1174dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1175dff31646SBarry Smith     cits[0] = 0;
1176dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
11771f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
11781f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
11791f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
11801f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
11811f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1182dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1183dff31646SBarry Smith   }
1184dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1185dff31646SBarry Smith 
1186c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
118704102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
118804102261SBarry Smith   {
118904102261SBarry Smith     PetscInt nmax = 2;
119004102261SBarry Smith     char     **buffs;
119104102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1192c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
119304102261SBarry Smith     if (flg1) {
119404102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
119504102261SBarry Smith       if (nmax == 1) {
119604102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
119704102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
119804102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
119904102261SBarry Smith       }
120004102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
120104102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
120204102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
120304102261SBarry Smith     }
120404102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
120504102261SBarry Smith   }
1206763ec7b1SBarry Smith   {
1207763ec7b1SBarry Smith     PetscInt nmax = 2;
1208763ec7b1SBarry Smith     char     **buffs;
1209763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1210763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1211763ec7b1SBarry Smith     if (flg1) {
1212763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1213763ec7b1SBarry Smith       if (nmax == 1) {
1214763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1215763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1216763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1217763ec7b1SBarry Smith       }
1218763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1219763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1220763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1221763ec7b1SBarry Smith     }
1222763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1223763ec7b1SBarry Smith   }
122404102261SBarry Smith #endif
122567234432SDmitry Karpeev   /*
122667234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
122767234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
122867234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
122967234432SDmitry Karpeev   */
123067234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
123104102261SBarry Smith 
12322d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12332d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12342d53ad75SBarry Smith #endif
12352d53ad75SBarry Smith 
12362d53ad75SBarry Smith 
1237e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1238dff31646SBarry Smith   flg = PETSC_FALSE;
1239c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1240d5649816SBarry Smith   if (flg) {
1241e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1242d5649816SBarry Smith   }
1243d5649816SBarry Smith #endif
1244d5649816SBarry Smith 
1245681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1246681455b2SBarry Smith   flg1 = PETSC_FALSE;
1247c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1248681455b2SBarry Smith   if (flg1) {
1249681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1250681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1251681455b2SBarry Smith   }
1252681455b2SBarry Smith #endif
1253681455b2SBarry Smith 
125467584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1255c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1256e5c89e4eSSatish Balay   if (!flg2) {
125790d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1258c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1259e5c89e4eSSatish Balay   }
1260e5c89e4eSSatish Balay   if (flg2) {
12610841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1262e5c89e4eSSatish Balay   }
126367584ceeSBarry Smith #endif
1264e5c89e4eSSatish Balay 
1265e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
126690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1267c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1268e5c89e4eSSatish Balay   if (flg1) {
1269e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1270205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1271e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1272e5c89e4eSSatish Balay   }
1273e5c89e4eSSatish Balay #endif
1274e5c89e4eSSatish Balay 
1275e5c89e4eSSatish Balay 
1276e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1277e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1278e5c89e4eSSatish Balay   mname[0] = 0;
1279c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1280e5c89e4eSSatish Balay   if (flg1) {
1281e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1282e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1283e5c89e4eSSatish Balay   }
1284e5c89e4eSSatish Balay #endif
1285356e5837SBarry Smith #endif
1286a297a907SKarl Rupp 
1287dd710f27SBarry Smith   /*
1288dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1289dd710f27SBarry Smith   */
1290dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1291dd710f27SBarry Smith 
1292356e5837SBarry Smith #if defined(PETSC_USE_LOG)
129387c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1294f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
129587c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
129687c3beb6SLisandro Dalcin 
1297356e5837SBarry Smith   mname[0] = 0;
1298c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1299e5c89e4eSSatish Balay   if (flg1) {
130091eabc43SBarry Smith     PetscViewer viewer;
130120a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
130291eabc43SBarry Smith     if (mname[0]) {
130391eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
130491eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13056bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
130633f85c2fSBarry Smith     } else {
130733f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13089a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
130933f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13109a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
131133f85c2fSBarry Smith     }
1312e5c89e4eSSatish Balay   }
1313a297a907SKarl Rupp 
1314dd710f27SBarry Smith   /*
1315dd710f27SBarry Smith      Free any objects created by the last block of code.
1316dd710f27SBarry Smith   */
1317dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1318dd710f27SBarry Smith 
1319dd710f27SBarry Smith   mname[0] = 0;
1320c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1321c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
13227ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1323e5c89e4eSSatish Balay #endif
132410463e74SBarry Smith 
132533f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
132610463e74SBarry Smith 
132790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1328c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1329e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
133090d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1331c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1332e5c89e4eSSatish Balay   if (flg1) {
1333e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1334e5c89e4eSSatish Balay   }
133590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
133690d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13378bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1338c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1339c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1340c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1341e4c476e2SSatish Balay 
1342e5c89e4eSSatish Balay   if (flg2) {
1343be56827dSJed Brown     PetscViewer viewer;
134402ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
134502ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1346c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1347be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1348e5c89e4eSSatish Balay   }
1349e5c89e4eSSatish Balay 
1350e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1351c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1352c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1353e5c89e4eSSatish Balay 
135433fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1355c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
1356c5929fdfSBarry Smith   ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
1357e5c89e4eSSatish Balay   if (flg3) {
1358e5c89e4eSSatish Balay     if (!flg2) { /* have not yet printed the options */
1359be56827dSJed Brown       PetscViewer viewer;
136002ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
136102ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1362c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1363be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1364e5c89e4eSSatish Balay     }
1365e5c89e4eSSatish Balay     if (!nopt) {
1366e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1367e5c89e4eSSatish Balay     } else if (nopt == 1) {
1368e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1369e5c89e4eSSatish Balay     } else {
13707582186dSLisandro Dalcin       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1371e5c89e4eSSatish Balay     }
1372df12ba86SBarry Smith   }
1373e5c89e4eSSatish Balay #if defined(PETSC_USE_DEBUG)
1374da8b8a77SBarry Smith   if (nopt && !flg3 && !flg1) {
1375e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
1376e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
1377c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1378e5c89e4eSSatish Balay   } else if (nopt && flg3) {
1379e5c89e4eSSatish Balay #else
1380e5c89e4eSSatish Balay   if (nopt && flg3) {
1381e5c89e4eSSatish Balay #endif
1382c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1383e5c89e4eSSatish Balay   }
1384e5c89e4eSSatish Balay 
1385e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1386d45a07a7SBarry Smith   if (!PetscGlobalRank) {
138787f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
138816ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1389d45a07a7SBarry Smith   }
1390ec957eceSBarry Smith #endif
1391ec957eceSBarry Smith 
13924097062eSBarry Smith #if defined(PETSC_USE_LOG)
139310463e74SBarry Smith   /*
1394dbc8283eSBarry Smith        List all objects the user may have forgot to free
13952eff7a51SBarry Smith   */
139605df10baSBarry Smith   if (PetscObjectsLog) {
1397c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1398a64a8e02SBarry Smith     if (flg1) {
1399a64a8e02SBarry Smith       MPI_Comm local_comm;
14007eb1d149SBarry Smith       char     string[64];
1401a64a8e02SBarry Smith 
1402c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1403a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1404a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14057eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1406a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1407a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14080a1571b3SBarry Smith     }
140905df10baSBarry Smith   }
14104097062eSBarry Smith #endif
14114097062eSBarry Smith 
14124097062eSBarry Smith #if defined(PETSC_USE_LOG)
1413dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1414dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1415a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14164097062eSBarry Smith #endif
14172eff7a51SBarry Smith 
141833f85c2fSBarry Smith   /*
141933f85c2fSBarry Smith      Destroy any packages that registered a finalize
142033f85c2fSBarry Smith   */
142133f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
142233f85c2fSBarry Smith 
1423101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1424fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1425101409b8SToby Isaac #endif
1426101409b8SToby Isaac 
142733f85c2fSBarry Smith   /*
142848dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
142948dd1dffSBarry Smith 
143037e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
143148dd1dffSBarry Smith   */
143237e93019SBarry Smith 
14334028d114SSatish Balay   if (petsc_history) {
1434f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
1435e5c89e4eSSatish Balay     petsc_history = 0;
1436e5c89e4eSSatish Balay   }
14379de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1438e5c89e4eSSatish Balay 
14390298fd71SBarry Smith   ierr = PetscInfoAllow(PETSC_FALSE,NULL);CHKERRQ(ierr);
1440e5c89e4eSSatish Balay 
144167584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14428bb29257SSatish Balay   {
1443e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
1444e5c89e4eSSatish Balay     FILE *fd;
1445ed9cf6e9SBarry Smith     int  err;
1446e5c89e4eSSatish Balay 
1447e5c89e4eSSatish Balay     fname[0] = 0;
1448a297a907SKarl Rupp 
1449c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1450dc92acbaSJed Brown     flg2 = PETSC_FALSE;
1451c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);
14528bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
1453dc92acbaSJed Brown     if (PETSC_RUNNING_ON_VALGRIND) flg2 = PETSC_FALSE;
1454dc92acbaSJed Brown #else
1455dc92acbaSJed Brown     flg2 = PETSC_FALSE;         /* Skip reporting for optimized builds regardless of -malloc_test */
1456dc92acbaSJed Brown #endif
1457e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1458e5c89e4eSSatish Balay       char sname[PETSC_MAX_PATH_LEN];
1459e5c89e4eSSatish Balay 
14602e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1461e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1462e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1463ed9cf6e9SBarry Smith       err  = fclose(fd);
1464e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1465dc92acbaSJed Brown     } else if (flg1 || flg2) {
1466e5c89e4eSSatish Balay       MPI_Comm local_comm;
1467e5c89e4eSSatish Balay 
1468e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1469e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1470e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1471e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1472e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1473e5c89e4eSSatish Balay     }
1474e5c89e4eSSatish Balay   }
1475a64a8e02SBarry Smith 
14768bb29257SSatish Balay   {
1477e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
14780298fd71SBarry Smith     FILE *fd = NULL;
1479e5c89e4eSSatish Balay 
1480e5c89e4eSSatish Balay     fname[0] = 0;
1481a297a907SKarl Rupp 
1482c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_log",fname,250,&flg1);CHKERRQ(ierr);
1483c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-malloc_log_threshold",&flg2);CHKERRQ(ierr);
1484e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1485ed9cf6e9SBarry Smith       int err;
1486e5c89e4eSSatish Balay 
1487574034a9SJed Brown       if (!rank) {
1488574034a9SJed Brown         fd = fopen(fname,"w");
1489574034a9SJed Brown         if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",fname);
1490574034a9SJed Brown       }
1491e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(fd);CHKERRQ(ierr);
1492574034a9SJed Brown       if (fd) {
1493ed9cf6e9SBarry Smith         err = fclose(fd);
1494e32f2f54SBarry Smith         if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1495574034a9SJed Brown       }
1496574034a9SJed Brown     } else if (flg1 || flg2) {
1497e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(stdout);CHKERRQ(ierr);
1498e5c89e4eSSatish Balay     }
1499e5c89e4eSSatish Balay   }
150067584ceeSBarry Smith #endif
150120e2c332SMatthew G. Knepley 
15025486ca60SMatthew G. Knepley   /*
15035486ca60SMatthew G. Knepley      Close any open dynamic libraries
15045486ca60SMatthew G. Knepley   */
15055486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15065486ca60SMatthew G. Knepley 
1507e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15084416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1509e5c89e4eSSatish Balay 
1510e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
1511e5c89e4eSSatish Balay   PetscGlobalArgs = 0;
1512e5c89e4eSSatish Balay 
1513008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1514e5c89e4eSSatish Balay 
1515dbc8283eSBarry Smith   /*
1516efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1517efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1518efb80d3cSBarry Smith 
1519efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1520efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1521dbc8283eSBarry Smith  */
1522b770b1f6SSatish Balay   {
1523dbc8283eSBarry Smith     PetscCommCounter *counter;
1524dbc8283eSBarry Smith     PetscMPIInt      flg;
1525dbc8283eSBarry Smith     MPI_Comm         icomm;
1526265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
152747435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1528dbc8283eSBarry Smith     if (flg) {
1529265f3f35SJed Brown       icomm = ucomm.comm;
153047435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1531dbc8283eSBarry 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");
1532dbc8283eSBarry Smith 
153347435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
153447435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1535efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1536dbc8283eSBarry Smith     }
153747435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1538dbc8283eSBarry Smith     if (flg) {
1539265f3f35SJed Brown       icomm = ucomm.comm;
154047435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1541dbc8283eSBarry 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");
1542dbc8283eSBarry Smith 
154347435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
154447435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1545efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1546dbc8283eSBarry Smith     }
1547b770b1f6SSatish Balay   }
1548dbc8283eSBarry Smith 
154947435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
155047435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
155147435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
15525f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1553480cf27aSJed Brown 
15545ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15555ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15565ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1557ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1558ef19f930SBarry Smith 
1559e5c89e4eSSatish Balay   if (PetscBeganMPI) {
156099608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
156199b1327fSBarry Smith     PetscMPIInt flag;
156299b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1563e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
156499608316SBarry Smith #endif
1565e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1566e5c89e4eSSatish Balay   }
1567e5c89e4eSSatish Balay /*
1568e5c89e4eSSatish Balay 
1569e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1570e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1571e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1572e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1573e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
15740e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1575e5c89e4eSSatish Balay    memory was not freed.
1576e5c89e4eSSatish Balay 
1577e5c89e4eSSatish Balay */
15781d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1579a297a907SKarl Rupp 
1580e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1581e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
15823db9a53dSBarry Smith   PetscFunctionReturn(0);
1583e5c89e4eSSatish Balay }
1584e5c89e4eSSatish Balay 
158543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
15868cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
158743db4dbbSBarry Smith {
158843db4dbbSBarry Smith   if (*a == *b) return 1;
158943db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
159043db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
159143db4dbbSBarry Smith   return 0;
159243db4dbbSBarry Smith }
1593a70650f6SBarry Smith #endif
159443db4dbbSBarry Smith 
159543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
15968cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
159743db4dbbSBarry Smith {
159843db4dbbSBarry Smith   if (*a == *b) return 1;
159943db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
160043db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
160143db4dbbSBarry Smith   return 0;
160243db4dbbSBarry Smith }
160343db4dbbSBarry Smith #endif
1604