xref: /petsc/src/sys/objects/pinit.c (revision 22e7e69cd8748a4e8659f05f7914590124c0fa91)
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};
466a6fc655SJed Brown const char *const PetscDataTypes[] = {"INT","DOUBLE","COMPLEX","LONG","SHORT","FLOAT",
472d53ad75SBarry Smith                                       "CHAR","LOGICAL","ENUM","BOOL","LONGDOUBLE","OBJECT","FUNCTION","PetscDataType","PETSC_",0};
48e5c89e4eSSatish Balay 
49ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
50ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
510f8e0872SSatish Balay 
52a2f94806SJed Brown PetscInt PetscHotRegionDepth;
53a2f94806SJed Brown 
54b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
55b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
56b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
57b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
58b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
59b22622e2STadeu Manoel #endif
60b22622e2STadeu Manoel 
61e5c89e4eSSatish Balay /*
62e5c89e4eSSatish Balay        Checks the options database for initializations related to the
63e5c89e4eSSatish Balay     PETSc components
64e5c89e4eSSatish Balay */
6595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode  PetscOptionsCheckInitial_Components(void)
66e5c89e4eSSatish Balay {
672d747510SLisandro Dalcin   PetscBool      flg;
68e5c89e4eSSatish Balay   PetscErrorCode ierr;
69e5c89e4eSSatish Balay 
70e5c89e4eSSatish Balay   PetscFunctionBegin;
712d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
722d747510SLisandro Dalcin   if (flg) {
73e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
74e8e7597cSSatish Balay     MPI_Comm comm = PETSC_COMM_WORLD;
75e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"------Additional PETSc component options--------\n");CHKERRQ(ierr);
768e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -log_exclude: <vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
778e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -info_exclude: <null,vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
78e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
79e5c89e4eSSatish Balay #endif
80e5c89e4eSSatish Balay   }
81e5c89e4eSSatish Balay   PetscFunctionReturn(0);
82e5c89e4eSSatish Balay }
83e5c89e4eSSatish Balay 
840f11a792SBarry Smith /*
85945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8672a42c3cSBarry Smith 
8772a42c3cSBarry Smith    Collective
8872a42c3cSBarry Smith 
8972a42c3cSBarry Smith    Level: advanced
9072a42c3cSBarry Smith 
9195452b02SPatrick Sanan     Notes:
92a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
930f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
94a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
950f11a792SBarry Smith 
96a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
971ea5a559SBarry Smith 
9872a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
990f11a792SBarry Smith */
100945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
10172a42c3cSBarry Smith {
10272a42c3cSBarry Smith   PetscErrorCode ierr;
10372a42c3cSBarry Smith   int            myargc   = argc;
10472a42c3cSBarry Smith   char           **myargs = args;
10572a42c3cSBarry Smith 
10672a42c3cSBarry Smith   PetscFunctionBegin;
1073bf036e2SBarry Smith   ierr = PetscInitialize(&myargc,&myargs,filename,help);CHKERRQ(ierr);
1081ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
109df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
11072a42c3cSBarry Smith   PetscFunctionReturn(ierr);
11172a42c3cSBarry Smith }
11272a42c3cSBarry Smith 
113f0865b08SBarry Smith /*
114a56f64adSBarry Smith       Used by Julia interface to get communicator
115f0865b08SBarry Smith */
116945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
117f0865b08SBarry Smith {
118f0865b08SBarry Smith   PetscFunctionBegin;
119f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
120f0865b08SBarry Smith   PetscFunctionReturn(0);
121f0865b08SBarry Smith }
122f0865b08SBarry Smith 
123e5c89e4eSSatish Balay /*@C
124e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
125e5c89e4eSSatish Balay         the command line arguments.
126e5c89e4eSSatish Balay 
127e5c89e4eSSatish Balay    Collective
128e5c89e4eSSatish Balay 
129e5c89e4eSSatish Balay    Level: advanced
130e5c89e4eSSatish Balay 
131e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
132e5c89e4eSSatish Balay @*/
1337087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
134e5c89e4eSSatish Balay {
135e5c89e4eSSatish Balay   PetscErrorCode ierr;
136e5c89e4eSSatish Balay   int            argc   = 0;
137e5c89e4eSSatish Balay   char           **args = 0;
138e5c89e4eSSatish Balay 
139e5c89e4eSSatish Balay   PetscFunctionBegin;
1400298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
141e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
142e5c89e4eSSatish Balay }
143e5c89e4eSSatish Balay 
144e5c89e4eSSatish Balay /*@
145e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
146e5c89e4eSSatish Balay 
14793b6d2d1SJed Brown    Level: beginner
148e5c89e4eSSatish Balay 
149e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
150e5c89e4eSSatish Balay @*/
1517087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool  *isInitialized)
152e5c89e4eSSatish Balay {
153e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
15493b6d2d1SJed Brown   return 0;
155e5c89e4eSSatish Balay }
156e5c89e4eSSatish Balay 
157e5c89e4eSSatish Balay /*@
158e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
159e5c89e4eSSatish Balay 
160e5c89e4eSSatish Balay    Level: developer
161e5c89e4eSSatish Balay 
162e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
163e5c89e4eSSatish Balay @*/
1647087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
165e5c89e4eSSatish Balay {
166e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
16793b6d2d1SJed Brown   return 0;
168e5c89e4eSSatish Balay }
169e5c89e4eSSatish Balay 
17095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void);
171e5c89e4eSSatish Balay 
172e5c89e4eSSatish Balay /*
173e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
174e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
175e5c89e4eSSatish Balay */
176367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
177e5c89e4eSSatish Balay 
178367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
179e5c89e4eSSatish Balay {
180e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
181e5c89e4eSSatish Balay 
182e5c89e4eSSatish Balay   PetscFunctionBegin;
183e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
184e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
185e5c89e4eSSatish Balay     MPI_Abort(MPI_COMM_WORLD,1);
186e5c89e4eSSatish Balay   }
187e5c89e4eSSatish Balay 
188e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
189e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
190e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
191e5c89e4eSSatish Balay   }
192812af9f3SBarry Smith   PetscFunctionReturnVoid();
193e5c89e4eSSatish Balay }
194e5c89e4eSSatish Balay 
195e5c89e4eSSatish Balay /*
196e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
197e5c89e4eSSatish Balay sum of the second entry.
198b693b147SBarry Smith 
19976f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
200367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
201b693b147SBarry Smith there would be no place to store the both needed results.
202e5c89e4eSSatish Balay */
20376ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
204e5c89e4eSSatish Balay {
205e5c89e4eSSatish Balay   PetscErrorCode ierr;
206e5c89e4eSSatish Balay 
207e5c89e4eSSatish Balay   PetscFunctionBegin;
208d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
209d6e4c47cSJed Brown   {
210d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
211367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
212d6e4c47cSJed Brown     *max = work.max;
213d6e4c47cSJed Brown     *sum = work.sum;
214d6e4c47cSJed Brown   }
215d6e4c47cSJed Brown #else
216d6e4c47cSJed Brown   {
217d6e4c47cSJed Brown     PetscMPIInt    size,rank;
218d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
219e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
220e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
221785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
222367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2236ac3741eSJed Brown     *max = work[rank].max;
2246ac3741eSJed Brown     *sum = work[rank].sum;
225e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
226d6e4c47cSJed Brown   }
227d6e4c47cSJed Brown #endif
228e5c89e4eSSatish Balay   PetscFunctionReturn(0);
229e5c89e4eSSatish Balay }
230e5c89e4eSSatish Balay 
231e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
232e5c89e4eSSatish Balay 
233570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
23406a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
235e5c89e4eSSatish Balay 
2368cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
237e5c89e4eSSatish Balay {
238e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
239e5c89e4eSSatish Balay 
240e5c89e4eSSatish Balay   PetscFunctionBegin;
2417c2de775SJed Brown   if (*datatype == MPIU_REAL) {
242e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
243a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2447c2de775SJed Brown   }
2457c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2467c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2477c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
248a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2497c2de775SJed Brown   }
2507c2de775SJed Brown #endif
2517c2de775SJed Brown   else {
2527c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
253e2e03761SBarry Smith     MPI_Abort(MPI_COMM_WORLD,1);
254e2e03761SBarry Smith   }
255812af9f3SBarry Smith   PetscFunctionReturnVoid();
256e5c89e4eSSatish Balay }
257e5c89e4eSSatish Balay #endif
258e5c89e4eSSatish Balay 
259570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
260d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
261d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
262d9822059SBarry Smith 
2638cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
264d9822059SBarry Smith {
265d9822059SBarry Smith   PetscInt i,count = *cnt;
266d9822059SBarry Smith 
267d9822059SBarry Smith   PetscFunctionBegin;
2687c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2698c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
270a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2717c2de775SJed Brown   }
2727c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2737c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2747c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2757c2de775SJed Brown     for (i=0; i<count; i++) {
2767c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2777c2de775SJed Brown     }
2787c2de775SJed Brown   }
2797c2de775SJed Brown #endif
2807c2de775SJed Brown   else {
2817c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
2828c764dc5SJose Roman     MPI_Abort(MPI_COMM_WORLD,1);
2838c764dc5SJose Roman   }
284d9822059SBarry Smith   PetscFunctionReturnVoid();
285d9822059SBarry Smith }
286d9822059SBarry Smith 
2878cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
288d9822059SBarry Smith {
289d9822059SBarry Smith   PetscInt    i,count = *cnt;
290d9822059SBarry Smith 
291d9822059SBarry Smith   PetscFunctionBegin;
2927c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2938c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
294a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2957c2de775SJed Brown   }
2967c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2977c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2987c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2997c2de775SJed Brown     for (i=0; i<count; i++) {
3007c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3017c2de775SJed Brown     }
3027c2de775SJed Brown   }
3037c2de775SJed Brown #endif
3047c2de775SJed Brown   else {
3058c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
3068c764dc5SJose Roman     MPI_Abort(MPI_COMM_WORLD,1);
3078c764dc5SJose Roman   }
308d9822059SBarry Smith   PetscFunctionReturnVoid();
309d9822059SBarry Smith }
310d9822059SBarry Smith #endif
311d9822059SBarry Smith 
312480cf27aSJed Brown /*
313480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
314480cf27aSJed Brown 
315ff0e51ddSBarry 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.
316480cf27aSJed Brown 
31712801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
318480cf27aSJed Brown 
319480cf27aSJed Brown */
3208cc058d9SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
321480cf27aSJed Brown {
322480cf27aSJed Brown   PetscErrorCode ierr;
323480cf27aSJed Brown 
324480cf27aSJed Brown   PetscFunctionBegin;
32512801b39SBarry Smith   ierr = PetscInfo1(0,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
32612801b39SBarry Smith   ierr = PetscFree(count_val);CHKERRMPI(ierr);
327480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
328480cf27aSJed Brown }
329480cf27aSJed Brown 
330480cf27aSJed Brown /*
33147435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
332da3039f7SJed Brown   calls MPI_Comm_free().
333da3039f7SJed Brown 
334da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
335480cf27aSJed Brown 
336ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
337480cf27aSJed Brown 
33812801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
339480cf27aSJed Brown 
340480cf27aSJed Brown */
341da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
342480cf27aSJed Brown {
343480cf27aSJed Brown   PetscErrorCode                    ierr;
344b89831e5SBarry Smith   PetscMPIInt                       flg;
345265f3f35SJed Brown   union {MPI_Comm comm; void *ptr;} icomm,ocomm;
346480cf27aSJed Brown 
347480cf27aSJed Brown   PetscFunctionBegin;
34812801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
349ec4fadc2SJed Brown   icomm.ptr = attr_val;
350da3039f7SJed Brown 
35147435625SJed Brown   ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
35212801b39SBarry Smith   if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm");
35312801b39SBarry Smith   if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm");
35447435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
35512801b39SBarry 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);
356da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
357b89831e5SBarry Smith }
358da3039f7SJed Brown 
359da3039f7SJed Brown /*
36047435625SJed 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.
361da3039f7SJed Brown  */
362da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
363da3039f7SJed Brown {
364da3039f7SJed Brown   PetscErrorCode ierr;
365da3039f7SJed Brown 
366da3039f7SJed Brown   PetscFunctionBegin;
36712801b39SBarry Smith   ierr = PetscInfo1(0,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
368480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
369480cf27aSJed Brown }
370480cf27aSJed Brown 
3715f7487a0SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Shm(MPI_Comm,PetscMPIInt,void *,void *);
37242218b76SBarry Smith 
373951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
374e39fd77fSBarry Smith #if !defined(PETSC_WORDS_BIGENDIAN)
3758cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3768cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3778cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
378e39fd77fSBarry Smith #endif
379951e3c8eSBarry Smith #endif
380e39fd77fSBarry Smith 
38112801b39SBarry Smith PetscMPIInt PETSC_MPI_ERROR_CLASS,PETSC_MPI_ERROR_CODE;
38212801b39SBarry Smith 
383eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
384eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3856ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
3866ae9a8a6SBarry Smith char **PetscGlobalArgs = 0;
387dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
388e5c89e4eSSatish Balay 
389dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
390051e4cf2SJed Brown {
391051e4cf2SJed Brown   PetscErrorCode ierr;
392051e4cf2SJed Brown 
393051e4cf2SJed Brown   PetscFunctionBegin;
394051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
395f191ac77SSatish 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);
396051e4cf2SJed 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);
397051e4cf2SJed Brown   PetscFunctionReturn(0);
398051e4cf2SJed Brown }
399e5c89e4eSSatish Balay 
4002d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4012d747510SLisandro Dalcin 
4022d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4032d747510SLisandro Dalcin {
4042d747510SLisandro Dalcin   PetscErrorCode ierr;
4052d747510SLisandro Dalcin 
4062d747510SLisandro Dalcin   PetscFunctionBegin;
4072d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
4082d747510SLisandro Dalcin   PetscFunctionReturn(0);
4092d747510SLisandro Dalcin }
4102d747510SLisandro Dalcin 
4112d747510SLisandro Dalcin /*@C
4122d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4132d747510SLisandro Dalcin 
4142d747510SLisandro Dalcin     Not Collective
4152d747510SLisandro Dalcin 
4162d747510SLisandro Dalcin     Input Parameter:
4172d747510SLisandro Dalcin .   len - length of the string name
4182d747510SLisandro Dalcin 
4192d747510SLisandro Dalcin     Output Parameter:
4202d747510SLisandro Dalcin .   name - the name of the running program
4212d747510SLisandro Dalcin 
4222d747510SLisandro Dalcin    Level: advanced
4232d747510SLisandro Dalcin 
4242d747510SLisandro Dalcin     Notes:
4252d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4262d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4272d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4282d747510SLisandro Dalcin @*/
4292d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4302d747510SLisandro Dalcin {
4312d747510SLisandro Dalcin   PetscErrorCode ierr;
4322d747510SLisandro Dalcin 
4332d747510SLisandro Dalcin   PetscFunctionBegin;
4342d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4352d747510SLisandro Dalcin   PetscFunctionReturn(0);
4362d747510SLisandro Dalcin }
4372d747510SLisandro Dalcin 
438e5c89e4eSSatish Balay /*@C
439e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
440e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
441e5c89e4eSSatish Balay 
442e5c89e4eSSatish Balay    Not Collective
443e5c89e4eSSatish Balay 
444e5c89e4eSSatish Balay    Output Parameters:
445e5c89e4eSSatish Balay +  argc - count of number of command line arguments
446e5c89e4eSSatish Balay -  args - the command line arguments
447e5c89e4eSSatish Balay 
448e5c89e4eSSatish Balay    Level: intermediate
449e5c89e4eSSatish Balay 
450e5c89e4eSSatish Balay    Notes:
451e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
452e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
453e5c89e4eSSatish Balay 
454f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
455f177e3b1SBarry Smith 
456e5c89e4eSSatish Balay    Concepts: command line arguments
457e5c89e4eSSatish Balay 
458793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
459e5c89e4eSSatish Balay 
460e5c89e4eSSatish Balay @*/
4617087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
462e5c89e4eSSatish Balay {
463e5c89e4eSSatish Balay   PetscFunctionBegin;
46417186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
465e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
466e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
467e5c89e4eSSatish Balay   PetscFunctionReturn(0);
468e5c89e4eSSatish Balay }
469e5c89e4eSSatish Balay 
470793721a6SBarry Smith /*@C
471793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
472793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
473793721a6SBarry Smith 
474793721a6SBarry Smith    Not Collective
475793721a6SBarry Smith 
476793721a6SBarry Smith    Output Parameters:
477793721a6SBarry Smith .  args - the command line arguments
478793721a6SBarry Smith 
479793721a6SBarry Smith    Level: intermediate
480793721a6SBarry Smith 
481793721a6SBarry Smith    Notes:
482793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
483793721a6SBarry Smith 
484793721a6SBarry Smith    Concepts: command line arguments
485793721a6SBarry Smith 
486793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
487793721a6SBarry Smith 
488793721a6SBarry Smith @*/
4897087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
490793721a6SBarry Smith {
491793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
492793721a6SBarry Smith   PetscErrorCode ierr;
493793721a6SBarry Smith 
494793721a6SBarry Smith   PetscFunctionBegin;
49517186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4962d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
497785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
498793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
499793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
500793721a6SBarry Smith   }
5012d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
502793721a6SBarry Smith   PetscFunctionReturn(0);
503793721a6SBarry Smith }
504793721a6SBarry Smith 
505793721a6SBarry Smith /*@C
506793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
507793721a6SBarry Smith 
508793721a6SBarry Smith    Not Collective
509793721a6SBarry Smith 
510793721a6SBarry Smith    Output Parameters:
511793721a6SBarry Smith .  args - the command line arguments
512793721a6SBarry Smith 
513793721a6SBarry Smith    Level: intermediate
514793721a6SBarry Smith 
515793721a6SBarry Smith    Concepts: command line arguments
516793721a6SBarry Smith 
517793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
518793721a6SBarry Smith 
519793721a6SBarry Smith @*/
5207087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
521793721a6SBarry Smith {
522793721a6SBarry Smith   PetscInt       i = 0;
523793721a6SBarry Smith   PetscErrorCode ierr;
524793721a6SBarry Smith 
525793721a6SBarry Smith   PetscFunctionBegin;
526a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
527793721a6SBarry Smith   while (args[i]) {
528793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
529793721a6SBarry Smith     i++;
530793721a6SBarry Smith   }
531793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
532793721a6SBarry Smith   PetscFunctionReturn(0);
533793721a6SBarry Smith }
534793721a6SBarry Smith 
53511525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
53630befbd2SBarry Smith #include <petscconfiginfo.h>
53730befbd2SBarry Smith 
53895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
53911525c0dSBarry Smith {
54011525c0dSBarry Smith   if (!PetscGlobalRank) {
54130befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
54211525c0dSBarry Smith     int            port;
543ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
54411525c0dSBarry Smith     size_t         applinelen,introlen;
54511525c0dSBarry Smith     PetscErrorCode ierr;
546ffbd1cfbSBarry Smith     char           sawsurl[256];
54711525c0dSBarry Smith 
548c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
54911525c0dSBarry Smith     if (flg) {
55011525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
55111525c0dSBarry Smith 
552c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
55311525c0dSBarry Smith       if (sawslog[0]) {
55411525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
55511525c0dSBarry Smith       } else {
55611525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
55711525c0dSBarry Smith       }
55811525c0dSBarry Smith     }
559c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
56011525c0dSBarry Smith     if (flg) {
56111525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
56211525c0dSBarry Smith     }
563c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
564ffbd1cfbSBarry Smith     if (selectport) {
565ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
566ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
567ffbd1cfbSBarry Smith     } else {
568c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
56911525c0dSBarry Smith       if (flg) {
57011525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
57111525c0dSBarry Smith       }
572ffbd1cfbSBarry Smith     }
573c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
57411525c0dSBarry Smith     if (flg) {
57511525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
57611525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5779c1e0ce8SBarry Smith     } else {
578c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5799c1e0ce8SBarry Smith       if (flg) {
5803c01dfcfSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
5819c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5829c1e0ce8SBarry Smith       }
58311525c0dSBarry Smith     }
584c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
58511525c0dSBarry Smith     if (flg2) {
58611525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
58711525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
58811525c0dSBarry Smith       ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr);
58911525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
59011525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
59143da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
59211525c0dSBarry Smith     }
59311525c0dSBarry Smith     ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr);
59411525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
59511525c0dSBarry Smith     introlen   = 4096 + applinelen;
59630a8c9c0SSurtai Han     applinelen += 1024;
59711525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
59811525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
59911525c0dSBarry Smith 
60011525c0dSBarry Smith     if (rootlocal) {
60111525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
60211525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
60311525c0dSBarry Smith     }
60476a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
60511525c0dSBarry Smith     if (rootlocal && help) {
606928bb9adSStefano 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);
60711525c0dSBarry Smith     } else if (help) {
608928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
60911525c0dSBarry Smith     } else {
610928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
61111525c0dSBarry Smith     }
612b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
61330befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
61411525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
61511525c0dSBarry 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"
616df62c222SBarry 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"
617928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
61843da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
61911525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
62011525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
621ffbd1cfbSBarry Smith     if (selectport) {
622aa573868SBarry Smith       PetscBool silent;
6237d812c46SBarry Smith 
6247d812c46SBarry Smith       ierr = SAWs_Initialize();
6257d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6267d812c46SBarry Smith       while (ierr) {
6277d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6287d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6297d812c46SBarry Smith         ierr = SAWs_Initialize();
6307d812c46SBarry Smith       }
6317d812c46SBarry Smith 
632aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
633aa573868SBarry Smith       if (!silent) {
634ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
635ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
636ffbd1cfbSBarry Smith       }
6377d812c46SBarry Smith     } else {
6387d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
639aa573868SBarry Smith     }
6400af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6410af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6420af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6430af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6440af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
64511525c0dSBarry Smith   }
64611525c0dSBarry Smith   PetscFunctionReturn(0);
64711525c0dSBarry Smith }
64811525c0dSBarry Smith #endif
64911525c0dSBarry Smith 
650a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
651a56f64adSBarry Smith #include <adios.h>
65222580e64SBarry Smith #include <adios_read.h>
6537b56e58cSSatish Balay int64_t Petsc_adios_group;
654a56f64adSBarry Smith #endif
65555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
65655d657eeSBarry Smith #include <adios2_c.h>
65755d657eeSBarry Smith #endif
658a56f64adSBarry Smith 
659e5c89e4eSSatish Balay /*@C
660e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
661e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
662e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
663e5c89e4eSSatish Balay    your program -- usually the very first line!
664e5c89e4eSSatish Balay 
665e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
666e5c89e4eSSatish Balay 
667e5c89e4eSSatish Balay    Input Parameters:
668e5c89e4eSSatish Balay +  argc - count of number of command line arguments
669e5c89e4eSSatish Balay .  args - the command line arguments
6700298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
671fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6720298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
673e5c89e4eSSatish Balay 
67405827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
67505827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
67605827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
67705827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
67805827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
679e5c89e4eSSatish Balay 
680e5c89e4eSSatish Balay    Options Database Keys:
6817ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6827ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
683e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
684e5c89e4eSSatish Balay .  -on_error_emacs <machinename> causes emacsclient to jump to error file
685b52f573bSBarry Smith .  -on_error_abort calls abort() when error detected (no traceback)
686e8fb0fc0SBarry Smith .  -on_error_mpiabort calls MPI_abort() when error detected
687e8fb0fc0SBarry Smith .  -error_output_stderr prints error messages to stderr instead of the default stdout
688e8fb0fc0SBarry Smith .  -error_output_none does not print the error messages (but handles errors in the same way as if this was not called)
689e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
690e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
691e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
692e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
6932fb0ec9aSBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries)
694e5c89e4eSSatish Balay .  -malloc no - Indicates not to use error-checking malloc
6952fb0ec9aSBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free
696aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
697dc92acbaSJed Brown .  -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds
698e5c89e4eSSatish Balay .  -fp_trap - Stops on floating point exceptions (Note that on the
699e5c89e4eSSatish Balay               IBM RS6000 this slows code by at least a factor of 10.)
700e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
701e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
702e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
703e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
704e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7050841954dSBarry Smith -  -memory_view - Print memory usage at end of run
706e5c89e4eSSatish Balay 
707e5c89e4eSSatish Balay    Options Database Keys for Profiling:
708a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
709495fc317SBarry Smith +  -info <optional filename> - Prints verbose information to the screen
710495fc317SBarry Smith .  -info_exclude <null,vec,mat,pc,ksp,snes,ts> - Excludes some of the verbose messages
711217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
712217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
713495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
714e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7159a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
7169a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
717495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
718f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
719495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
720495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
721c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
72287c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
723c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
724495fc317SBarry Smith 
725609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
726e5c89e4eSSatish Balay 
727ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
728ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
729ffbd1cfbSBarry 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
730ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
731ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
732ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
733ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
734ffbd1cfbSBarry Smith 
735e5c89e4eSSatish Balay    Environmental Variables:
736e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
737e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
738e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
739e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
740e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
741e5c89e4eSSatish Balay 
742e5c89e4eSSatish Balay 
743e5c89e4eSSatish Balay    Level: beginner
744e5c89e4eSSatish Balay 
745e5c89e4eSSatish Balay    Notes:
746e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
747e5c89e4eSSatish Balay    it before PetscInitialize().
748e5c89e4eSSatish Balay 
749e5c89e4eSSatish Balay    Fortran Version:
750e5c89e4eSSatish Balay    In Fortran this routine has the format
751e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
752e5c89e4eSSatish Balay 
753e5c89e4eSSatish Balay +   ierr - error return code
7540eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
755fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
756e5c89e4eSSatish Balay 
757e5c89e4eSSatish Balay    Important Fortran Note:
7580eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7590298fd71SBarry Smith    null character string; you CANNOT just use NULL as
760a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
761e5c89e4eSSatish Balay 
76201cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
76301cb0274SBarry Smith    calling PetscInitialize().
764e5c89e4eSSatish Balay 
765e5c89e4eSSatish Balay    Concepts: initializing PETSc
766e5c89e4eSSatish Balay 
76701cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
768e5c89e4eSSatish Balay 
769e5c89e4eSSatish Balay @*/
7707087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
771e5c89e4eSSatish Balay {
772e5c89e4eSSatish Balay   PetscErrorCode ierr;
7734bb5149bSJed Brown   PetscMPIInt    flag, size;
77419c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
775e5c89e4eSSatish Balay   char           hostname[256];
776e5c89e4eSSatish Balay 
777e5c89e4eSSatish Balay   PetscFunctionBegin;
778e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7793d96e996SBarry Smith   /*
7803d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7813d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7823d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7833d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7843d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7853d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7863d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7873d96e996SBarry Smith       listed above and since that time are compatible.
7883d96e996SBarry Smith 
7893d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
7903d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
7913d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
7923d96e996SBarry Smith       to perform the checking.
7933d96e996SBarry Smith 
7943d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
7953d96e996SBarry Smith 
7963d96e996SBarry Smith       Questions:
7973d96e996SBarry Smith 
7983d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
7993d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
8003d96e996SBarry Smith   */
8013d96e996SBarry Smith 
80219c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
80319c5658aSBarry Smith   {
80419c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
80519c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
80619c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8073d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
80819c5658aSBarry Smith #if defined(MPICH_VERSION)
8093d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
81019c5658aSBarry Smith     {
81119c5658aSBarry Smith       char *ver,*lf;
81219c5658aSBarry Smith       flg = PETSC_FALSE;
81319c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
81419c5658aSBarry Smith       if (ver) {
81519c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
81619c5658aSBarry Smith         if (lf) {
81719c5658aSBarry Smith           *lf = 0;
81819c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
81919c5658aSBarry Smith         }
82019c5658aSBarry Smith       }
82119c5658aSBarry Smith       if (!flg) {
82219c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8233d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
82419c5658aSBarry Smith       }
82519c5658aSBarry Smith     }
8263d96e996SBarry Smith #endif
8273d96e996SBarry 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?) */
82819c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
82919c5658aSBarry Smith     {
83019c5658aSBarry Smith       char *ver,bs[32],*bsf;
83119c5658aSBarry Smith       flg = PETSC_FALSE;
83219c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
83319c5658aSBarry Smith       if (ver) {
8342e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
83519c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
83619c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
83719c5658aSBarry Smith       }
83819c5658aSBarry Smith       if (!flg) {
83919c5658aSBarry 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);
8403d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
84119c5658aSBarry Smith       }
84219c5658aSBarry Smith     }
84319c5658aSBarry Smith #endif
84419c5658aSBarry Smith   }
84519c5658aSBarry Smith #endif
84619c5658aSBarry Smith 
847e5c89e4eSSatish Balay 
848ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
849d89683f4Sbcordonn   PETSC_STDOUT = stdout;
850ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
851e5c89e4eSSatish Balay 
8520c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8530c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8540c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8550c30907bSSatish Balay #endif
8560c30907bSSatish Balay 
8574416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
858e5c89e4eSSatish Balay 
859e5c89e4eSSatish Balay   /*
860e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
861e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
862e5c89e4eSSatish Balay   */
863e5c89e4eSSatish Balay   if (argc && *argc) {
864e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
865e5c89e4eSSatish Balay   } else {
866e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
867e5c89e4eSSatish Balay   }
868e5c89e4eSSatish Balay 
869e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
870e5c89e4eSSatish Balay   if (!flag) {
871e32f2f54SBarry 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");
8725e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8735e765c61SJed Brown     {
8745e765c61SJed Brown       PetscMPIInt provided;
8755e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
8765e765c61SJed Brown     }
8775e765c61SJed Brown #else
878e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
8795e765c61SJed Brown #endif
880e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
881e5c89e4eSSatish Balay   }
882e5c89e4eSSatish Balay   if (argc && args) {
883e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
884e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
885e5c89e4eSSatish Balay   }
886e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
8875ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
8885ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
8895ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
890ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
891e5c89e4eSSatish Balay 
892a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
893d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
894e5c89e4eSSatish Balay 
89512801b39SBarry Smith   ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
89612801b39SBarry Smith   ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
89712801b39SBarry Smith 
898e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
899e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
900e5c89e4eSSatish Balay 
901e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
902e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
903e5c89e4eSSatish Balay 
9048ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9058ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9068ad47952SJed Brown 
907e5c89e4eSSatish Balay   /*
908e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
909e5c89e4eSSatish Balay      shared libraries the constructors for global variables
910e5c89e4eSSatish Balay      are not called; at least on IRIX.
911e5c89e4eSSatish Balay   */
912886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
913e5c89e4eSSatish Balay   {
914252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
91550f81f78SJed Brown     PetscComplex ic(0.0,1.0);
916e5c89e4eSSatish Balay     PETSC_i = ic;
917252ecd31SSatish Balay #else
91850f81f78SJed Brown     PETSC_i = _Complex_I;
919b7940d39SSatish Balay #endif
920762437b8SSatish Balay   }
921762437b8SSatish Balay 
9222c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
923e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
924500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
925500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
926500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9272c876bd9SBarry Smith #endif
928886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
929e5c89e4eSSatish Balay 
930e5c89e4eSSatish Balay   /*
931e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
932e5c89e4eSSatish Balay      half of the entries and maxes the second half.
933e5c89e4eSSatish Balay   */
934367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
935e5c89e4eSSatish Balay 
936ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
937c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
938c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9397c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9408c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9418c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9428c764dc5SJose Roman #endif
943d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
944d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
945570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
946570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
947570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
948570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
949570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
950c90a1750SBarry Smith #endif
951c90a1750SBarry Smith 
952570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
953cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
954cca4cb22SSatish Balay #endif
955cca4cb22SSatish Balay 
956e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
957e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
958e5c89e4eSSatish Balay 
95944041f26SJed Brown #if defined(PETSC_USE_64BIT_INDICES) || !defined(MPI_2INT)
960e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
961e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
96244041f26SJed Brown #endif
963e5c89e4eSSatish Balay 
964ec957eceSBarry Smith 
965e5c89e4eSSatish Balay   /*
966480cf27aSJed Brown      Attributes to be set on PETSc communicators
967480cf27aSJed Brown   */
96812801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
96912801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
97012801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
9715f7487a0SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Shm,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
972480cf27aSJed Brown 
973480cf27aSJed Brown   /*
974e8fb0fc0SBarry Smith      Build the options database
975e5c89e4eSSatish Balay   */
976c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
977e5c89e4eSSatish Balay 
978703f3eceSBarry Smith   /* call a second time so it can look in the options database */
979703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
9806dc8fec2Sbcordonn 
981e5c89e4eSSatish Balay   /*
982e5c89e4eSSatish Balay      Print main application help message
983e5c89e4eSSatish Balay   */
9842d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
985e5c89e4eSSatish Balay   if (help && flg) {
986e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
987e5c89e4eSSatish Balay   }
988e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
989e5c89e4eSSatish Balay 
990d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
991d45a07a7SBarry Smith 
992e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
99311525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
994f4202a44SBarry Smith #endif
995f4202a44SBarry Smith 
996896238b9SBarry Smith   /* Creates the logging data structures; this is enabled even if logging is not turned on */
997a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
998896238b9SBarry Smith   ierr = PetscLogInitialize();CHKERRQ(ierr);
999a9f03627SSatish Balay #endif
1000e5c89e4eSSatish Balay 
1001e5c89e4eSSatish Balay   /*
1002e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1003e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1004e5c89e4eSSatish Balay   */
1005e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1006e5c89e4eSSatish Balay 
1007e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
1008ae15b995SBarry Smith   ierr = PetscInfo1(0,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1009e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
1010ae15b995SBarry Smith   ierr = PetscInfo1(0,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1011e5c89e4eSSatish Balay 
1012e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr);
1013ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1014c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1015ef6c6fedSBoyana Norris 
1016951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1017e39fd77fSBarry Smith   /*
1018e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1019e39fd77fSBarry Smith 
1020e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1021e39fd77fSBarry Smith   */
1022e39fd77fSBarry Smith #if !defined(PETSC_WORDS_BIGENDIAN)
10230298fd71SBarry Smith   ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
1024e39fd77fSBarry Smith #endif
1025951e3c8eSBarry Smith #endif
1026e39fd77fSBarry Smith 
102741c0b4b3SShri Abhyankar   /*
102841c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
102941c0b4b3SShri Abhyankar   */
1030ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1031e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1032e1167bb9SShri Abhyankar #endif
1033e1167bb9SShri Abhyankar 
10342d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10352d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10362d53ad75SBarry Smith #endif
10372d53ad75SBarry Smith 
10385e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
103987c3beb6SLisandro Dalcin   {
104087c3beb6SLisandro Dalcin     PetscViewer viewer;
1041*22e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
10425e71baefSBarry Smith     if (flg) {
10435e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
104487c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
104587c3beb6SLisandro Dalcin     }
10465e71baefSBarry Smith   }
10475e71baefSBarry Smith #endif
1048dff31646SBarry Smith 
104987c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
105087c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
105187c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
105287c3beb6SLisandro Dalcin 
1053a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1054a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
10557b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1056a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
10579fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1058a56f64adSBarry Smith #endif
105955d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
106055d657eeSBarry Smith #endif
1061a56f64adSBarry Smith 
1062301d30feSBarry Smith   /*
1063301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1064301d30feSBarry Smith   */
1065301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
10662db0e300SLisandro Dalcin 
10672db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
10682db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1069301d30feSBarry Smith   PetscFunctionReturn(0);
1070e5c89e4eSSatish Balay }
1071e5c89e4eSSatish Balay 
10724097062eSBarry Smith #if defined(PETSC_USE_LOG)
107395c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1074ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1075ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
107695c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
10774097062eSBarry Smith #endif
1078e5c89e4eSSatish Balay 
1079008a6e76SBarry Smith /*
1080008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1081008a6e76SBarry Smith */
1082008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1083008a6e76SBarry Smith {
1084008a6e76SBarry Smith   PetscErrorCode ierr;
1085008a6e76SBarry Smith 
1086008a6e76SBarry Smith   PetscFunctionBegin;
1087008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1088008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1089008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1090008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1091008a6e76SBarry Smith #endif
1092008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1093008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1094008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1095008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1096008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1097008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1098008a6e76SBarry Smith #endif
1099008a6e76SBarry Smith 
1100008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1101008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1102008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1103008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1104008a6e76SBarry Smith #endif
1105008a6e76SBarry Smith #endif
1106008a6e76SBarry Smith 
1107008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1108008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1109008a6e76SBarry Smith #endif
1110008a6e76SBarry Smith 
1111008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
1112008a6e76SBarry Smith #if defined(PETSC_USE_64BIT_INDICES) || !defined(MPI_2INT)
1113008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1114008a6e76SBarry Smith #endif
1115008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1116008a6e76SBarry Smith   PetscFunctionReturn(0);
1117008a6e76SBarry Smith }
1118008a6e76SBarry Smith 
1119e5c89e4eSSatish Balay /*@C
1120e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1121e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1122e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1123e5c89e4eSSatish Balay 
1124e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1125e5c89e4eSSatish Balay 
1126e5c89e4eSSatish Balay    Options Database Keys:
112726a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1128e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11297eb1d149SBarry 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
1130e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1131e5c89e4eSSatish Balay .  -malloc_dump - Calls PetscMallocDump()
1132e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
1133e5c89e4eSSatish Balay -  -malloc_log - Prints summary of memory usage
1134e5c89e4eSSatish Balay 
1135e5c89e4eSSatish Balay    Level: beginner
1136e5c89e4eSSatish Balay 
1137e5c89e4eSSatish Balay    Note:
1138e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1139e5c89e4eSSatish Balay 
114088c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1141e5c89e4eSSatish Balay @*/
11427087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1143e5c89e4eSSatish Balay {
1144e5c89e4eSSatish Balay   PetscErrorCode ierr;
11454bb5149bSJed Brown   PetscMPIInt    rank;
1146a8d2bbe5SBarry Smith   PetscInt       nopt;
11472bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1148dff31646SBarry Smith   PetscBool      flg;
114910463e74SBarry Smith #if defined(PETSC_USE_LOG)
115010463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
115110463e74SBarry Smith #endif
1152e5c89e4eSSatish Balay 
1153e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
11544b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
11553cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1156e5c89e4eSSatish Balay   }
11573cbbc5ffSBarry Smith   PetscFunctionBegin;
11580298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1159b022a5c1SBarry Smith 
11601f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1161a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
116222580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1163a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1164a56f64adSBarry Smith #endif
116555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
116655d657eeSBarry Smith #endif
1167c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1168dff31646SBarry Smith   if (flg) {
11691f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
11701f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
11711f817a21SBarry Smith 
1172c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11731f817a21SBarry Smith     if (filename[0]) {
11741f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
11751f817a21SBarry Smith     }
1176dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1177dff31646SBarry Smith     cits[0] = 0;
1178dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
11791f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
11801f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
11811f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
11821f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
11831f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1184dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1185dff31646SBarry Smith   }
1186dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1187dff31646SBarry Smith 
1188c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
118904102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
119004102261SBarry Smith   {
119104102261SBarry Smith     PetscInt nmax = 2;
119204102261SBarry Smith     char     **buffs;
119304102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1194c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
119504102261SBarry Smith     if (flg1) {
119604102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
119704102261SBarry Smith       if (nmax == 1) {
119804102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
119904102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
120004102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
120104102261SBarry Smith       }
120204102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
120304102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
120404102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
120504102261SBarry Smith     }
120604102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
120704102261SBarry Smith   }
1208763ec7b1SBarry Smith   {
1209763ec7b1SBarry Smith     PetscInt nmax = 2;
1210763ec7b1SBarry Smith     char     **buffs;
1211763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1212763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1213763ec7b1SBarry Smith     if (flg1) {
1214763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1215763ec7b1SBarry Smith       if (nmax == 1) {
1216763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1217763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1218763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1219763ec7b1SBarry Smith       }
1220763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1221763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1222763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1223763ec7b1SBarry Smith     }
1224763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1225763ec7b1SBarry Smith   }
122604102261SBarry Smith #endif
122767234432SDmitry Karpeev   /*
122867234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
122967234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
123067234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
123167234432SDmitry Karpeev   */
123267234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
123304102261SBarry Smith 
12342d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12352d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12362d53ad75SBarry Smith #endif
12372d53ad75SBarry Smith 
12382d53ad75SBarry Smith 
1239e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1240dff31646SBarry Smith   flg = PETSC_FALSE;
1241c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1242d5649816SBarry Smith   if (flg) {
1243e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1244d5649816SBarry Smith   }
1245d5649816SBarry Smith #endif
1246d5649816SBarry Smith 
1247681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1248681455b2SBarry Smith   flg1 = PETSC_FALSE;
1249c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1250681455b2SBarry Smith   if (flg1) {
1251681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1252681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1253681455b2SBarry Smith   }
1254681455b2SBarry Smith #endif
1255681455b2SBarry Smith 
125667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1257c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1258e5c89e4eSSatish Balay   if (!flg2) {
125990d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1260c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1261e5c89e4eSSatish Balay   }
1262e5c89e4eSSatish Balay   if (flg2) {
12630841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1264e5c89e4eSSatish Balay   }
126567584ceeSBarry Smith #endif
1266e5c89e4eSSatish Balay 
1267e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
126890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1269c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1270e5c89e4eSSatish Balay   if (flg1) {
1271e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1272205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1273e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1274e5c89e4eSSatish Balay   }
1275e5c89e4eSSatish Balay #endif
1276e5c89e4eSSatish Balay 
1277e5c89e4eSSatish Balay 
1278e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1279e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1280e5c89e4eSSatish Balay   mname[0] = 0;
1281c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1282e5c89e4eSSatish Balay   if (flg1) {
1283e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1284e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1285e5c89e4eSSatish Balay   }
1286e5c89e4eSSatish Balay #endif
1287356e5837SBarry Smith #endif
1288a297a907SKarl Rupp 
1289dd710f27SBarry Smith   /*
1290dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1291dd710f27SBarry Smith   */
1292dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1293dd710f27SBarry Smith 
1294356e5837SBarry Smith #if defined(PETSC_USE_LOG)
129587c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1296f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
129787c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
129887c3beb6SLisandro Dalcin 
1299356e5837SBarry Smith   mname[0] = 0;
1300c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1301e5c89e4eSSatish Balay   if (flg1) {
130291eabc43SBarry Smith     PetscViewer viewer;
130320a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
130491eabc43SBarry Smith     if (mname[0]) {
130591eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
130691eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13076bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
130833f85c2fSBarry Smith     } else {
130933f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13109a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
131133f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13129a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
131333f85c2fSBarry Smith     }
1314e5c89e4eSSatish Balay   }
1315a297a907SKarl Rupp 
1316dd710f27SBarry Smith   /*
1317dd710f27SBarry Smith      Free any objects created by the last block of code.
1318dd710f27SBarry Smith   */
1319dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1320dd710f27SBarry Smith 
1321dd710f27SBarry Smith   mname[0] = 0;
1322c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1323c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
13247ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1325e5c89e4eSSatish Balay #endif
132610463e74SBarry Smith 
132733f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
132810463e74SBarry Smith 
132990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1330c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1331e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
133290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1333c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1334e5c89e4eSSatish Balay   if (flg1) {
1335e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1336e5c89e4eSSatish Balay   }
133790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
133890d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13398bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1340c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1341c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1342c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1343e4c476e2SSatish Balay 
1344e5c89e4eSSatish Balay   if (flg2) {
1345be56827dSJed Brown     PetscViewer viewer;
134602ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
134702ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1348c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1349be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1350e5c89e4eSSatish Balay   }
1351e5c89e4eSSatish Balay 
1352e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1353c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1354c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1355e5c89e4eSSatish Balay 
135633fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1357c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
1358c5929fdfSBarry Smith   ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
1359e5c89e4eSSatish Balay   if (flg3) {
1360e5c89e4eSSatish Balay     if (!flg2) { /* have not yet printed the options */
1361be56827dSJed Brown       PetscViewer viewer;
136202ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
136302ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1364c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1365be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1366e5c89e4eSSatish Balay     }
1367e5c89e4eSSatish Balay     if (!nopt) {
1368e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1369e5c89e4eSSatish Balay     } else if (nopt == 1) {
1370e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1371e5c89e4eSSatish Balay     } else {
13727582186dSLisandro Dalcin       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1373e5c89e4eSSatish Balay     }
1374df12ba86SBarry Smith   }
1375e5c89e4eSSatish Balay #if defined(PETSC_USE_DEBUG)
1376da8b8a77SBarry Smith   if (nopt && !flg3 && !flg1) {
1377e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
1378e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
1379c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1380e5c89e4eSSatish Balay   } else if (nopt && flg3) {
1381e5c89e4eSSatish Balay #else
1382e5c89e4eSSatish Balay   if (nopt && flg3) {
1383e5c89e4eSSatish Balay #endif
1384c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1385e5c89e4eSSatish Balay   }
1386e5c89e4eSSatish Balay 
1387e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1388d45a07a7SBarry Smith   if (!PetscGlobalRank) {
138987f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
139016ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1391d45a07a7SBarry Smith   }
1392ec957eceSBarry Smith #endif
1393ec957eceSBarry Smith 
13944097062eSBarry Smith #if defined(PETSC_USE_LOG)
139510463e74SBarry Smith   /*
1396dbc8283eSBarry Smith        List all objects the user may have forgot to free
13972eff7a51SBarry Smith   */
139805df10baSBarry Smith   if (PetscObjectsLog) {
1399c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1400a64a8e02SBarry Smith     if (flg1) {
1401a64a8e02SBarry Smith       MPI_Comm local_comm;
14027eb1d149SBarry Smith       char     string[64];
1403a64a8e02SBarry Smith 
1404c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1405a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1406a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14077eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1408a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1409a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14100a1571b3SBarry Smith     }
141105df10baSBarry Smith   }
14124097062eSBarry Smith #endif
14134097062eSBarry Smith 
14144097062eSBarry Smith #if defined(PETSC_USE_LOG)
1415dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1416dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1417a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14184097062eSBarry Smith #endif
14192eff7a51SBarry Smith 
142033f85c2fSBarry Smith   /*
142133f85c2fSBarry Smith      Destroy any packages that registered a finalize
142233f85c2fSBarry Smith   */
142333f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
142433f85c2fSBarry Smith 
1425101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1426fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1427101409b8SToby Isaac #endif
1428101409b8SToby Isaac 
142933f85c2fSBarry Smith   /*
143048dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
143148dd1dffSBarry Smith 
143237e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
143348dd1dffSBarry Smith   */
143437e93019SBarry Smith 
14354028d114SSatish Balay   if (petsc_history) {
1436f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
1437e5c89e4eSSatish Balay     petsc_history = 0;
1438e5c89e4eSSatish Balay   }
14399de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1440e5c89e4eSSatish Balay 
14410298fd71SBarry Smith   ierr = PetscInfoAllow(PETSC_FALSE,NULL);CHKERRQ(ierr);
1442e5c89e4eSSatish Balay 
144367584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14448bb29257SSatish Balay   {
1445e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
1446e5c89e4eSSatish Balay     FILE *fd;
1447ed9cf6e9SBarry Smith     int  err;
1448e5c89e4eSSatish Balay 
1449e5c89e4eSSatish Balay     fname[0] = 0;
1450a297a907SKarl Rupp 
1451c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1452dc92acbaSJed Brown     flg2 = PETSC_FALSE;
1453c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);
14548bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
1455dc92acbaSJed Brown     if (PETSC_RUNNING_ON_VALGRIND) flg2 = PETSC_FALSE;
1456dc92acbaSJed Brown #else
1457dc92acbaSJed Brown     flg2 = PETSC_FALSE;         /* Skip reporting for optimized builds regardless of -malloc_test */
1458dc92acbaSJed Brown #endif
1459e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1460e5c89e4eSSatish Balay       char sname[PETSC_MAX_PATH_LEN];
1461e5c89e4eSSatish Balay 
14622e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1463e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1464e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1465ed9cf6e9SBarry Smith       err  = fclose(fd);
1466e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1467dc92acbaSJed Brown     } else if (flg1 || flg2) {
1468e5c89e4eSSatish Balay       MPI_Comm local_comm;
1469e5c89e4eSSatish Balay 
1470e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1471e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1472e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1473e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1474e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1475e5c89e4eSSatish Balay     }
1476e5c89e4eSSatish Balay   }
1477a64a8e02SBarry Smith 
14788bb29257SSatish Balay   {
1479e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
14800298fd71SBarry Smith     FILE *fd = NULL;
1481e5c89e4eSSatish Balay 
1482e5c89e4eSSatish Balay     fname[0] = 0;
1483a297a907SKarl Rupp 
1484c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_log",fname,250,&flg1);CHKERRQ(ierr);
1485c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-malloc_log_threshold",&flg2);CHKERRQ(ierr);
1486e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1487ed9cf6e9SBarry Smith       int err;
1488e5c89e4eSSatish Balay 
1489574034a9SJed Brown       if (!rank) {
1490574034a9SJed Brown         fd = fopen(fname,"w");
1491574034a9SJed Brown         if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",fname);
1492574034a9SJed Brown       }
1493e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(fd);CHKERRQ(ierr);
1494574034a9SJed Brown       if (fd) {
1495ed9cf6e9SBarry Smith         err = fclose(fd);
1496e32f2f54SBarry Smith         if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1497574034a9SJed Brown       }
1498574034a9SJed Brown     } else if (flg1 || flg2) {
1499e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(stdout);CHKERRQ(ierr);
1500e5c89e4eSSatish Balay     }
1501e5c89e4eSSatish Balay   }
150267584ceeSBarry Smith #endif
150320e2c332SMatthew G. Knepley 
15045486ca60SMatthew G. Knepley   /*
15055486ca60SMatthew G. Knepley      Close any open dynamic libraries
15065486ca60SMatthew G. Knepley   */
15075486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15085486ca60SMatthew G. Knepley 
1509e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15104416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1511e5c89e4eSSatish Balay 
1512e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
1513e5c89e4eSSatish Balay   PetscGlobalArgs = 0;
1514e5c89e4eSSatish Balay 
1515008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1516e5c89e4eSSatish Balay 
1517dbc8283eSBarry Smith   /*
1518efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1519efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1520efb80d3cSBarry Smith 
1521efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1522efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1523dbc8283eSBarry Smith  */
1524b770b1f6SSatish Balay   {
1525dbc8283eSBarry Smith     PetscCommCounter *counter;
1526dbc8283eSBarry Smith     PetscMPIInt      flg;
1527dbc8283eSBarry Smith     MPI_Comm         icomm;
1528265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
152947435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1530dbc8283eSBarry Smith     if (flg) {
1531265f3f35SJed Brown       icomm = ucomm.comm;
153247435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1533dbc8283eSBarry 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");
1534dbc8283eSBarry Smith 
153547435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
153647435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1537efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1538dbc8283eSBarry Smith     }
153947435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1540dbc8283eSBarry Smith     if (flg) {
1541265f3f35SJed Brown       icomm = ucomm.comm;
154247435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1543dbc8283eSBarry 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");
1544dbc8283eSBarry Smith 
154547435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
154647435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1547efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1548dbc8283eSBarry Smith     }
1549b770b1f6SSatish Balay   }
1550dbc8283eSBarry Smith 
155147435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
155247435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
155347435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
15545f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1555480cf27aSJed Brown 
15565ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15575ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15585ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1559ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1560ef19f930SBarry Smith 
1561e5c89e4eSSatish Balay   if (PetscBeganMPI) {
156299608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
156399b1327fSBarry Smith     PetscMPIInt flag;
156499b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1565e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
156699608316SBarry Smith #endif
1567e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1568e5c89e4eSSatish Balay   }
1569e5c89e4eSSatish Balay /*
1570e5c89e4eSSatish Balay 
1571e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1572e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1573e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1574e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1575e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
15760e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1577e5c89e4eSSatish Balay    memory was not freed.
1578e5c89e4eSSatish Balay 
1579e5c89e4eSSatish Balay */
15801d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1581a297a907SKarl Rupp 
1582e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1583e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
15843db9a53dSBarry Smith   PetscFunctionReturn(0);
1585e5c89e4eSSatish Balay }
1586e5c89e4eSSatish Balay 
158743db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
15888cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
158943db4dbbSBarry Smith {
159043db4dbbSBarry Smith   if (*a == *b) return 1;
159143db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
159243db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
159343db4dbbSBarry Smith   return 0;
159443db4dbbSBarry Smith }
1595a70650f6SBarry Smith #endif
159643db4dbbSBarry Smith 
159743db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
15988cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
159943db4dbbSBarry Smith {
160043db4dbbSBarry Smith   if (*a == *b) return 1;
160143db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
160243db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
160343db4dbbSBarry Smith   return 0;
160443db4dbbSBarry Smith }
160543db4dbbSBarry Smith #endif
1606