xref: /petsc/src/sys/objects/pinit.c (revision 2d7475101bf02d93e875ae8727f1067be42c09d1)
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;
3942218b76SBarry Smith PetscMPIInt Petsc_Shared_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 {
67*2d747510SLisandro Dalcin   PetscBool      flg;
68e5c89e4eSSatish Balay   PetscErrorCode ierr;
69e5c89e4eSSatish Balay 
70e5c89e4eSSatish Balay   PetscFunctionBegin;
71*2d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
72*2d747510SLisandro 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:
9295452b02SPatrick Sanan     this is called only by the PETSc MATLAB and 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
94945d1669SBarry Smith      be called multiple times from MATLAB and Julia without the problem of trying to initialize MPI more than once.
950f11a792SBarry Smith 
961ea5a559SBarry Smith      Turns off PETSc signal handling because that can interact with MATLAB's signal handling causing random crashes.
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 /*
114945d1669SBarry Smith       Used by MATLAB and 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 
37142218b76SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelShared(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 
400*2d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
401*2d747510SLisandro Dalcin 
402*2d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
403*2d747510SLisandro Dalcin {
404*2d747510SLisandro Dalcin   PetscErrorCode ierr;
405*2d747510SLisandro Dalcin 
406*2d747510SLisandro Dalcin   PetscFunctionBegin;
407*2d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
408*2d747510SLisandro Dalcin   PetscFunctionReturn(0);
409*2d747510SLisandro Dalcin }
410*2d747510SLisandro Dalcin 
411*2d747510SLisandro Dalcin /*@C
412*2d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
413*2d747510SLisandro Dalcin 
414*2d747510SLisandro Dalcin     Not Collective
415*2d747510SLisandro Dalcin 
416*2d747510SLisandro Dalcin     Input Parameter:
417*2d747510SLisandro Dalcin .   len - length of the string name
418*2d747510SLisandro Dalcin 
419*2d747510SLisandro Dalcin     Output Parameter:
420*2d747510SLisandro Dalcin .   name - the name of the running program
421*2d747510SLisandro Dalcin 
422*2d747510SLisandro Dalcin    Level: advanced
423*2d747510SLisandro Dalcin 
424*2d747510SLisandro Dalcin     Notes:
425*2d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
426*2d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
427*2d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
428*2d747510SLisandro Dalcin @*/
429*2d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
430*2d747510SLisandro Dalcin {
431*2d747510SLisandro Dalcin   PetscErrorCode ierr;
432*2d747510SLisandro Dalcin 
433*2d747510SLisandro Dalcin   PetscFunctionBegin;
434*2d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
435*2d747510SLisandro Dalcin   PetscFunctionReturn(0);
436*2d747510SLisandro Dalcin }
437*2d747510SLisandro 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()");
496*2d747510SLisandro 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   }
501*2d747510SLisandro 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) {
60630befbd2SBarry Smith       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);
60711525c0dSBarry Smith     } else if (help) {
60830a8c9c0SSurtai Han       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);
60911525c0dSBarry Smith     } else {
61030befbd2SBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);
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"
617df62c222SBarry Smith                                     "%s",version,petscconfigureoptions,appline);
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 
650e5c89e4eSSatish Balay /*@C
651e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
652e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
653e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
654e5c89e4eSSatish Balay    your program -- usually the very first line!
655e5c89e4eSSatish Balay 
656e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
657e5c89e4eSSatish Balay 
658e5c89e4eSSatish Balay    Input Parameters:
659e5c89e4eSSatish Balay +  argc - count of number of command line arguments
660e5c89e4eSSatish Balay .  args - the command line arguments
6610298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
662fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6630298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
664e5c89e4eSSatish Balay 
66505827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
66605827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
66705827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
66805827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
66905827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
670e5c89e4eSSatish Balay 
671e5c89e4eSSatish Balay    Options Database Keys:
6727ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6737ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
674e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
675e5c89e4eSSatish Balay .  -on_error_emacs <machinename> causes emacsclient to jump to error file
676b52f573bSBarry Smith .  -on_error_abort calls abort() when error detected (no traceback)
677e8fb0fc0SBarry Smith .  -on_error_mpiabort calls MPI_abort() when error detected
678e8fb0fc0SBarry Smith .  -error_output_stderr prints error messages to stderr instead of the default stdout
679e8fb0fc0SBarry Smith .  -error_output_none does not print the error messages (but handles errors in the same way as if this was not called)
680e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
681e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
682e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
683e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
6842fb0ec9aSBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries)
685e5c89e4eSSatish Balay .  -malloc no - Indicates not to use error-checking malloc
6862fb0ec9aSBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free
687aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
688dc92acbaSJed Brown .  -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds
689e5c89e4eSSatish Balay .  -fp_trap - Stops on floating point exceptions (Note that on the
690e5c89e4eSSatish Balay               IBM RS6000 this slows code by at least a factor of 10.)
691e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
692e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
693e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
694e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
695e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
6960841954dSBarry Smith -  -memory_view - Print memory usage at end of run
697e5c89e4eSSatish Balay 
698e5c89e4eSSatish Balay    Options Database Keys for Profiling:
699a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
700495fc317SBarry Smith +  -info <optional filename> - Prints verbose information to the screen
701495fc317SBarry Smith .  -info_exclude <null,vec,mat,pc,ksp,snes,ts> - Excludes some of the verbose messages
702217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
703217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
704495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
705e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7069a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
7079a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
708495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
709f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
710495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
711495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
712c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
71387c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
714c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
715495fc317SBarry Smith 
716609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
717e5c89e4eSSatish Balay 
718ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
719ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
720ffbd1cfbSBarry 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
721ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
722ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
723ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
724ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
725ffbd1cfbSBarry Smith 
726e5c89e4eSSatish Balay    Environmental Variables:
727e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
728e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
729e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
730e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
731e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
732e5c89e4eSSatish Balay 
733e5c89e4eSSatish Balay 
734e5c89e4eSSatish Balay    Level: beginner
735e5c89e4eSSatish Balay 
736e5c89e4eSSatish Balay    Notes:
737e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
738e5c89e4eSSatish Balay    it before PetscInitialize().
739e5c89e4eSSatish Balay 
740e5c89e4eSSatish Balay    Fortran Version:
741e5c89e4eSSatish Balay    In Fortran this routine has the format
742e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
743e5c89e4eSSatish Balay 
744e5c89e4eSSatish Balay +   ierr - error return code
7450eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
746fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
747e5c89e4eSSatish Balay 
748e5c89e4eSSatish Balay    Important Fortran Note:
7490eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7500298fd71SBarry Smith    null character string; you CANNOT just use NULL as
751a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
752e5c89e4eSSatish Balay 
75301cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
75401cb0274SBarry Smith    calling PetscInitialize().
755e5c89e4eSSatish Balay 
756e5c89e4eSSatish Balay    Concepts: initializing PETSc
757e5c89e4eSSatish Balay 
75801cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
759e5c89e4eSSatish Balay 
760e5c89e4eSSatish Balay @*/
7617087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
762e5c89e4eSSatish Balay {
763e5c89e4eSSatish Balay   PetscErrorCode ierr;
7644bb5149bSJed Brown   PetscMPIInt    flag, size;
76519c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
766e5c89e4eSSatish Balay   char           hostname[256];
767e5c89e4eSSatish Balay 
768e5c89e4eSSatish Balay   PetscFunctionBegin;
769e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7703d96e996SBarry Smith   /*
7713d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7723d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7733d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7743d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7753d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7763d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7773d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7783d96e996SBarry Smith       listed above and since that time are compatible.
7793d96e996SBarry Smith 
7803d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
7813d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
7823d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
7833d96e996SBarry Smith       to perform the checking.
7843d96e996SBarry Smith 
7853d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
7863d96e996SBarry Smith 
7873d96e996SBarry Smith       Questions:
7883d96e996SBarry Smith 
7893d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
7903d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
7913d96e996SBarry Smith   */
7923d96e996SBarry Smith 
79319c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
79419c5658aSBarry Smith   {
79519c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
79619c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
79719c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
7983d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
79919c5658aSBarry Smith #if defined(MPICH_VERSION)
8003d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
80119c5658aSBarry Smith     {
80219c5658aSBarry Smith       char *ver,*lf;
80319c5658aSBarry Smith       flg = PETSC_FALSE;
80419c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
80519c5658aSBarry Smith       if (ver) {
80619c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
80719c5658aSBarry Smith         if (lf) {
80819c5658aSBarry Smith           *lf = 0;
80919c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
81019c5658aSBarry Smith         }
81119c5658aSBarry Smith       }
81219c5658aSBarry Smith       if (!flg) {
81319c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8143d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
81519c5658aSBarry Smith       }
81619c5658aSBarry Smith     }
8173d96e996SBarry Smith #endif
8183d96e996SBarry 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?) */
81919c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
82019c5658aSBarry Smith     {
82119c5658aSBarry Smith       char *ver,bs[32],*bsf;
82219c5658aSBarry Smith       flg = PETSC_FALSE;
82319c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
82419c5658aSBarry Smith       if (ver) {
8252e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
82619c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
82719c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
82819c5658aSBarry Smith       }
82919c5658aSBarry Smith       if (!flg) {
83019c5658aSBarry 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);
8313d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
83219c5658aSBarry Smith       }
83319c5658aSBarry Smith     }
83419c5658aSBarry Smith #endif
83519c5658aSBarry Smith   }
83619c5658aSBarry Smith #endif
83719c5658aSBarry Smith 
838e5c89e4eSSatish Balay 
839ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
840d89683f4Sbcordonn   PETSC_STDOUT = stdout;
841ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
842e5c89e4eSSatish Balay 
8430c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8440c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8450c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8460c30907bSSatish Balay #endif
8470c30907bSSatish Balay 
8484416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
849e5c89e4eSSatish Balay 
850e5c89e4eSSatish Balay   /*
851e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
852e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
853e5c89e4eSSatish Balay   */
854e5c89e4eSSatish Balay   if (argc && *argc) {
855e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
856e5c89e4eSSatish Balay   } else {
857e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
858e5c89e4eSSatish Balay   }
859e5c89e4eSSatish Balay 
860e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
861e5c89e4eSSatish Balay   if (!flag) {
862e32f2f54SBarry 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");
8635e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8645e765c61SJed Brown     {
8655e765c61SJed Brown       PetscMPIInt provided;
8665e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
8675e765c61SJed Brown     }
8685e765c61SJed Brown #else
869e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
8705e765c61SJed Brown #endif
871e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
872e5c89e4eSSatish Balay   }
873e5c89e4eSSatish Balay   if (argc && args) {
874e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
875e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
876e5c89e4eSSatish Balay   }
877e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
8785ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
8795ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
8805ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
881ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
882e5c89e4eSSatish Balay 
883a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
884d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
885e5c89e4eSSatish Balay 
88612801b39SBarry Smith   ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
88712801b39SBarry Smith   ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
88812801b39SBarry Smith 
889e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
890e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
891e5c89e4eSSatish Balay 
892e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
893e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
894e5c89e4eSSatish Balay 
8958ad47952SJed Brown   MPIU_BOOL = MPI_INT;
8968ad47952SJed Brown   MPIU_ENUM = MPI_INT;
8978ad47952SJed Brown 
898e5c89e4eSSatish Balay   /*
899e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
900e5c89e4eSSatish Balay      shared libraries the constructors for global variables
901e5c89e4eSSatish Balay      are not called; at least on IRIX.
902e5c89e4eSSatish Balay   */
903886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
904e5c89e4eSSatish Balay   {
905252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
90650f81f78SJed Brown     PetscComplex ic(0.0,1.0);
907e5c89e4eSSatish Balay     PETSC_i = ic;
908252ecd31SSatish Balay #else
90950f81f78SJed Brown     PETSC_i = _Complex_I;
910b7940d39SSatish Balay #endif
911762437b8SSatish Balay   }
912762437b8SSatish Balay 
9132c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
914e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
915500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
916500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
917500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9182c876bd9SBarry Smith #endif
919886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
920e5c89e4eSSatish Balay 
921e5c89e4eSSatish Balay   /*
922e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
923e5c89e4eSSatish Balay      half of the entries and maxes the second half.
924e5c89e4eSSatish Balay   */
925367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
926e5c89e4eSSatish Balay 
927ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
928c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
929c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9307c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9318c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9328c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9338c764dc5SJose Roman #endif
934d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
935d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
936570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
937570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
938570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
939570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
940570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
941c90a1750SBarry Smith #endif
942c90a1750SBarry Smith 
943570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
944cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
945cca4cb22SSatish Balay #endif
946cca4cb22SSatish Balay 
947e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
948e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
949e5c89e4eSSatish Balay 
95044041f26SJed Brown #if defined(PETSC_USE_64BIT_INDICES) || !defined(MPI_2INT)
951e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
952e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
95344041f26SJed Brown #endif
954e5c89e4eSSatish Balay 
955ec957eceSBarry Smith 
956e5c89e4eSSatish Balay   /*
957480cf27aSJed Brown      Attributes to be set on PETSc communicators
958480cf27aSJed Brown   */
95912801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
96012801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
96112801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
96212801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelShared,&Petsc_Shared_keyval,(void*)0);CHKERRQ(ierr);
963480cf27aSJed Brown 
964480cf27aSJed Brown   /*
965e8fb0fc0SBarry Smith      Build the options database
966e5c89e4eSSatish Balay   */
967c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
968e5c89e4eSSatish Balay 
969703f3eceSBarry Smith   /* call a second time so it can look in the options database */
970703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
9716dc8fec2Sbcordonn 
972e5c89e4eSSatish Balay   /*
973e5c89e4eSSatish Balay      Print main application help message
974e5c89e4eSSatish Balay   */
975*2d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
976e5c89e4eSSatish Balay   if (help && flg) {
977e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
978e5c89e4eSSatish Balay   }
979e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
980e5c89e4eSSatish Balay 
981d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
982d45a07a7SBarry Smith 
983e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
98411525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
985f4202a44SBarry Smith #endif
986f4202a44SBarry Smith 
987896238b9SBarry Smith   /* Creates the logging data structures; this is enabled even if logging is not turned on */
988a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
989896238b9SBarry Smith   ierr = PetscLogInitialize();CHKERRQ(ierr);
990a9f03627SSatish Balay #endif
991e5c89e4eSSatish Balay 
992e5c89e4eSSatish Balay   /*
993e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
994e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
995e5c89e4eSSatish Balay   */
996e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
997e5c89e4eSSatish Balay 
998e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
999ae15b995SBarry Smith   ierr = PetscInfo1(0,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1000e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
1001ae15b995SBarry Smith   ierr = PetscInfo1(0,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1002e5c89e4eSSatish Balay 
1003e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr);
1004ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1005c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1006ef6c6fedSBoyana Norris 
1007951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1008e39fd77fSBarry Smith   /*
1009e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1010e39fd77fSBarry Smith 
1011e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1012e39fd77fSBarry Smith   */
1013e39fd77fSBarry Smith #if !defined(PETSC_WORDS_BIGENDIAN)
10140298fd71SBarry Smith   ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
1015e39fd77fSBarry Smith #endif
1016951e3c8eSBarry Smith #endif
1017e39fd77fSBarry Smith 
101841c0b4b3SShri Abhyankar   /*
101941c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
102041c0b4b3SShri Abhyankar   */
1021ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1022e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1023e1167bb9SShri Abhyankar #endif
1024e1167bb9SShri Abhyankar 
10252d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10262d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10272d53ad75SBarry Smith #endif
10282d53ad75SBarry Smith 
10295e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
103087c3beb6SLisandro Dalcin   {
103187c3beb6SLisandro Dalcin     PetscViewer viewer;
10325e71baefSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
10335e71baefSBarry Smith     if (flg) {
10345e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
103587c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
103687c3beb6SLisandro Dalcin     }
10375e71baefSBarry Smith   }
10385e71baefSBarry Smith #endif
1039dff31646SBarry Smith 
104087c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
104187c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
104287c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
104387c3beb6SLisandro Dalcin 
1044301d30feSBarry Smith   /*
1045301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1046301d30feSBarry Smith   */
1047301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
10482db0e300SLisandro Dalcin 
10492db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
10502db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1051301d30feSBarry Smith   PetscFunctionReturn(0);
1052e5c89e4eSSatish Balay }
1053e5c89e4eSSatish Balay 
10544097062eSBarry Smith #if defined(PETSC_USE_LOG)
105595c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1056ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1057ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
105895c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
10594097062eSBarry Smith #endif
1060e5c89e4eSSatish Balay 
1061008a6e76SBarry Smith /*
1062008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1063008a6e76SBarry Smith */
1064008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1065008a6e76SBarry Smith {
1066008a6e76SBarry Smith   PetscErrorCode ierr;
1067008a6e76SBarry Smith 
1068008a6e76SBarry Smith   PetscFunctionBegin;
1069008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1070008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1071008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1072008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1073008a6e76SBarry Smith #endif
1074008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1075008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1076008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1077008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1078008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1079008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1080008a6e76SBarry Smith #endif
1081008a6e76SBarry Smith 
1082008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1083008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1084008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1085008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1086008a6e76SBarry Smith #endif
1087008a6e76SBarry Smith #endif
1088008a6e76SBarry Smith 
1089008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1090008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1091008a6e76SBarry Smith #endif
1092008a6e76SBarry Smith 
1093008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
1094008a6e76SBarry Smith #if defined(PETSC_USE_64BIT_INDICES) || !defined(MPI_2INT)
1095008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1096008a6e76SBarry Smith #endif
1097008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1098008a6e76SBarry Smith   PetscFunctionReturn(0);
1099008a6e76SBarry Smith }
1100008a6e76SBarry Smith 
1101e5c89e4eSSatish Balay /*@C
1102e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1103e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1104e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1105e5c89e4eSSatish Balay 
1106e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1107e5c89e4eSSatish Balay 
1108e5c89e4eSSatish Balay    Options Database Keys:
110988c29154SBarry Smith +  -options_table - Calls PetscOptionsView()
1110e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11117eb1d149SBarry 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
1112e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1113e5c89e4eSSatish Balay .  -malloc_dump - Calls PetscMallocDump()
1114e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
1115e5c89e4eSSatish Balay -  -malloc_log - Prints summary of memory usage
1116e5c89e4eSSatish Balay 
1117e5c89e4eSSatish Balay    Level: beginner
1118e5c89e4eSSatish Balay 
1119e5c89e4eSSatish Balay    Note:
1120e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1121e5c89e4eSSatish Balay 
112288c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1123e5c89e4eSSatish Balay @*/
11247087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1125e5c89e4eSSatish Balay {
1126e5c89e4eSSatish Balay   PetscErrorCode ierr;
11274bb5149bSJed Brown   PetscMPIInt    rank;
1128a8d2bbe5SBarry Smith   PetscInt       nopt;
11292bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1130dff31646SBarry Smith   PetscBool      flg;
113110463e74SBarry Smith #if defined(PETSC_USE_LOG)
113210463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
113310463e74SBarry Smith #endif
1134e5c89e4eSSatish Balay 
1135e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
11364b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
11373cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1138e5c89e4eSSatish Balay   }
11393cbbc5ffSBarry Smith   PetscFunctionBegin;
11400298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1141b022a5c1SBarry Smith 
11421f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
11431f817a21SBarry Smith 
1144c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1145dff31646SBarry Smith   if (flg) {
11461f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
11471f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
11481f817a21SBarry Smith 
1149c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11501f817a21SBarry Smith     if (filename[0]) {
11511f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
11521f817a21SBarry Smith     }
1153dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1154dff31646SBarry Smith     cits[0] = 0;
1155dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
11561f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
11571f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
11581f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
11591f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
11601f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1161dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1162dff31646SBarry Smith   }
1163dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1164dff31646SBarry Smith 
1165c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
116604102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
116704102261SBarry Smith   {
116804102261SBarry Smith     PetscInt nmax = 2;
116904102261SBarry Smith     char     **buffs;
117004102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1171c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
117204102261SBarry Smith     if (flg1) {
117304102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
117404102261SBarry Smith       if (nmax == 1) {
117504102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
117604102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
117704102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
117804102261SBarry Smith       }
117904102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
118004102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
118104102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
118204102261SBarry Smith     }
118304102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
118404102261SBarry Smith   }
1185763ec7b1SBarry Smith   {
1186763ec7b1SBarry Smith     PetscInt nmax = 2;
1187763ec7b1SBarry Smith     char     **buffs;
1188763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1189763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1190763ec7b1SBarry Smith     if (flg1) {
1191763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1192763ec7b1SBarry Smith       if (nmax == 1) {
1193763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1194763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1195763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1196763ec7b1SBarry Smith       }
1197763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1198763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1199763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1200763ec7b1SBarry Smith     }
1201763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1202763ec7b1SBarry Smith   }
120304102261SBarry Smith #endif
120467234432SDmitry Karpeev   /*
120567234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
120667234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
120767234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
120867234432SDmitry Karpeev   */
120967234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
121004102261SBarry Smith 
12112d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12122d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12132d53ad75SBarry Smith #endif
12142d53ad75SBarry Smith 
12152d53ad75SBarry Smith 
1216e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1217dff31646SBarry Smith   flg = PETSC_FALSE;
1218c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1219d5649816SBarry Smith   if (flg) {
1220e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1221d5649816SBarry Smith   }
1222d5649816SBarry Smith #endif
1223d5649816SBarry Smith 
1224681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1225681455b2SBarry Smith   flg1 = PETSC_FALSE;
1226c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1227681455b2SBarry Smith   if (flg1) {
1228681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1229681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1230681455b2SBarry Smith   }
1231681455b2SBarry Smith #endif
1232681455b2SBarry Smith 
123367584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1234c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1235e5c89e4eSSatish Balay   if (!flg2) {
123690d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1237c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1238e5c89e4eSSatish Balay   }
1239e5c89e4eSSatish Balay   if (flg2) {
12400841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1241e5c89e4eSSatish Balay   }
124267584ceeSBarry Smith #endif
1243e5c89e4eSSatish Balay 
1244e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
124590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1246c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1247e5c89e4eSSatish Balay   if (flg1) {
1248e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1249205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1250e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1251e5c89e4eSSatish Balay   }
1252e5c89e4eSSatish Balay #endif
1253e5c89e4eSSatish Balay 
1254e5c89e4eSSatish Balay 
1255e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1256e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1257e5c89e4eSSatish Balay   mname[0] = 0;
1258c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1259e5c89e4eSSatish Balay   if (flg1) {
1260e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1261e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1262e5c89e4eSSatish Balay   }
1263e5c89e4eSSatish Balay #endif
1264356e5837SBarry Smith #endif
1265a297a907SKarl Rupp 
1266dd710f27SBarry Smith   /*
1267dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1268dd710f27SBarry Smith   */
1269dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1270dd710f27SBarry Smith 
1271356e5837SBarry Smith #if defined(PETSC_USE_LOG)
127287c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1273f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
127487c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
127587c3beb6SLisandro Dalcin 
1276356e5837SBarry Smith   mname[0] = 0;
1277c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1278e5c89e4eSSatish Balay   if (flg1) {
127991eabc43SBarry Smith     PetscViewer viewer;
128020a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
128191eabc43SBarry Smith     if (mname[0]) {
128291eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
128391eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
12846bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
128533f85c2fSBarry Smith     } else {
128633f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
12879a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
128833f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
12899a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
129033f85c2fSBarry Smith     }
1291e5c89e4eSSatish Balay   }
1292a297a907SKarl Rupp 
1293dd710f27SBarry Smith   /*
1294dd710f27SBarry Smith      Free any objects created by the last block of code.
1295dd710f27SBarry Smith   */
1296dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1297dd710f27SBarry Smith 
1298dd710f27SBarry Smith   mname[0] = 0;
1299c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1300c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
13017ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1302e5c89e4eSSatish Balay #endif
130310463e74SBarry Smith 
130433f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
130510463e74SBarry Smith 
130690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1307c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1308e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
130990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1310c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1311e5c89e4eSSatish Balay   if (flg1) {
1312e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1313e5c89e4eSSatish Balay   }
131490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
131590d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13168bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1317c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1318c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1319c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1320e4c476e2SSatish Balay 
1321e5c89e4eSSatish Balay   if (flg2) {
1322be56827dSJed Brown     PetscViewer viewer;
132302ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
132402ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1325c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1326be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1327e5c89e4eSSatish Balay   }
1328e5c89e4eSSatish Balay 
1329e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1330c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1331c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1332e5c89e4eSSatish Balay 
133333fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1334c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
1335c5929fdfSBarry Smith   ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
1336e5c89e4eSSatish Balay   if (flg3) {
1337e5c89e4eSSatish Balay     if (!flg2) { /* have not yet printed the options */
1338be56827dSJed Brown       PetscViewer viewer;
133902ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
134002ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1341c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1342be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1343e5c89e4eSSatish Balay     }
1344e5c89e4eSSatish Balay     if (!nopt) {
1345e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1346e5c89e4eSSatish Balay     } else if (nopt == 1) {
1347e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1348e5c89e4eSSatish Balay     } else {
13497582186dSLisandro Dalcin       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1350e5c89e4eSSatish Balay     }
1351df12ba86SBarry Smith   }
1352e5c89e4eSSatish Balay #if defined(PETSC_USE_DEBUG)
1353da8b8a77SBarry Smith   if (nopt && !flg3 && !flg1) {
1354e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
1355e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
1356c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1357e5c89e4eSSatish Balay   } else if (nopt && flg3) {
1358e5c89e4eSSatish Balay #else
1359e5c89e4eSSatish Balay   if (nopt && flg3) {
1360e5c89e4eSSatish Balay #endif
1361c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1362e5c89e4eSSatish Balay   }
1363e5c89e4eSSatish Balay 
1364e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1365d45a07a7SBarry Smith   if (!PetscGlobalRank) {
136687f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
136716ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1368d45a07a7SBarry Smith   }
1369ec957eceSBarry Smith #endif
1370ec957eceSBarry Smith 
13714097062eSBarry Smith #if defined(PETSC_USE_LOG)
137210463e74SBarry Smith   /*
1373dbc8283eSBarry Smith        List all objects the user may have forgot to free
13742eff7a51SBarry Smith   */
137505df10baSBarry Smith   if (PetscObjectsLog) {
1376c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1377a64a8e02SBarry Smith     if (flg1) {
1378a64a8e02SBarry Smith       MPI_Comm local_comm;
13797eb1d149SBarry Smith       char     string[64];
1380a64a8e02SBarry Smith 
1381c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1382a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1383a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
13847eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1385a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1386a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
13870a1571b3SBarry Smith     }
138805df10baSBarry Smith   }
13894097062eSBarry Smith #endif
13904097062eSBarry Smith 
13914097062eSBarry Smith #if defined(PETSC_USE_LOG)
1392dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1393dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1394a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
13954097062eSBarry Smith #endif
13962eff7a51SBarry Smith 
139733f85c2fSBarry Smith   /*
139833f85c2fSBarry Smith      Destroy any packages that registered a finalize
139933f85c2fSBarry Smith   */
140033f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
140133f85c2fSBarry Smith 
1402101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1403fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1404101409b8SToby Isaac #endif
1405101409b8SToby Isaac 
140633f85c2fSBarry Smith   /*
140748dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
140848dd1dffSBarry Smith 
140937e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
141048dd1dffSBarry Smith   */
141137e93019SBarry Smith 
14124028d114SSatish Balay   if (petsc_history) {
1413f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
1414e5c89e4eSSatish Balay     petsc_history = 0;
1415e5c89e4eSSatish Balay   }
14169de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1417e5c89e4eSSatish Balay 
14180298fd71SBarry Smith   ierr = PetscInfoAllow(PETSC_FALSE,NULL);CHKERRQ(ierr);
1419e5c89e4eSSatish Balay 
142067584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14218bb29257SSatish Balay   {
1422e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
1423e5c89e4eSSatish Balay     FILE *fd;
1424ed9cf6e9SBarry Smith     int  err;
1425e5c89e4eSSatish Balay 
1426e5c89e4eSSatish Balay     fname[0] = 0;
1427a297a907SKarl Rupp 
1428c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1429dc92acbaSJed Brown     flg2 = PETSC_FALSE;
1430c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);
14318bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
1432dc92acbaSJed Brown     if (PETSC_RUNNING_ON_VALGRIND) flg2 = PETSC_FALSE;
1433dc92acbaSJed Brown #else
1434dc92acbaSJed Brown     flg2 = PETSC_FALSE;         /* Skip reporting for optimized builds regardless of -malloc_test */
1435dc92acbaSJed Brown #endif
1436e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1437e5c89e4eSSatish Balay       char sname[PETSC_MAX_PATH_LEN];
1438e5c89e4eSSatish Balay 
14392e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1440e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1441e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1442ed9cf6e9SBarry Smith       err  = fclose(fd);
1443e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1444dc92acbaSJed Brown     } else if (flg1 || flg2) {
1445e5c89e4eSSatish Balay       MPI_Comm local_comm;
1446e5c89e4eSSatish Balay 
1447e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1448e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1449e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1450e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1451e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1452e5c89e4eSSatish Balay     }
1453e5c89e4eSSatish Balay   }
1454a64a8e02SBarry Smith 
14558bb29257SSatish Balay   {
1456e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
14570298fd71SBarry Smith     FILE *fd = NULL;
1458e5c89e4eSSatish Balay 
1459e5c89e4eSSatish Balay     fname[0] = 0;
1460a297a907SKarl Rupp 
1461c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_log",fname,250,&flg1);CHKERRQ(ierr);
1462c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-malloc_log_threshold",&flg2);CHKERRQ(ierr);
1463e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1464ed9cf6e9SBarry Smith       int err;
1465e5c89e4eSSatish Balay 
1466574034a9SJed Brown       if (!rank) {
1467574034a9SJed Brown         fd = fopen(fname,"w");
1468574034a9SJed Brown         if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",fname);
1469574034a9SJed Brown       }
1470e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(fd);CHKERRQ(ierr);
1471574034a9SJed Brown       if (fd) {
1472ed9cf6e9SBarry Smith         err = fclose(fd);
1473e32f2f54SBarry Smith         if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1474574034a9SJed Brown       }
1475574034a9SJed Brown     } else if (flg1 || flg2) {
1476e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(stdout);CHKERRQ(ierr);
1477e5c89e4eSSatish Balay     }
1478e5c89e4eSSatish Balay   }
147967584ceeSBarry Smith #endif
148020e2c332SMatthew G. Knepley 
14815486ca60SMatthew G. Knepley   /*
14825486ca60SMatthew G. Knepley      Close any open dynamic libraries
14835486ca60SMatthew G. Knepley   */
14845486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
14855486ca60SMatthew G. Knepley 
1486e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
14874416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1488e5c89e4eSSatish Balay 
1489e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
1490e5c89e4eSSatish Balay   PetscGlobalArgs = 0;
1491e5c89e4eSSatish Balay 
1492008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1493e5c89e4eSSatish Balay 
1494dbc8283eSBarry Smith   /*
1495efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1496efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1497efb80d3cSBarry Smith 
1498efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1499efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1500dbc8283eSBarry Smith  */
1501b770b1f6SSatish Balay   {
1502dbc8283eSBarry Smith     PetscCommCounter *counter;
1503dbc8283eSBarry Smith     PetscMPIInt      flg;
1504dbc8283eSBarry Smith     MPI_Comm         icomm;
1505265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
150647435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1507dbc8283eSBarry Smith     if (flg) {
1508265f3f35SJed Brown       icomm = ucomm.comm;
150947435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1510dbc8283eSBarry 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");
1511dbc8283eSBarry Smith 
151247435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
151347435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1514efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1515dbc8283eSBarry Smith     }
151647435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1517dbc8283eSBarry Smith     if (flg) {
1518265f3f35SJed Brown       icomm = ucomm.comm;
151947435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1520dbc8283eSBarry 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");
1521dbc8283eSBarry Smith 
152247435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
152347435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1524efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1525dbc8283eSBarry Smith     }
1526b770b1f6SSatish Balay   }
1527dbc8283eSBarry Smith 
152847435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
152947435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
153047435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
153147435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Shared_keyval);CHKERRQ(ierr);
1532480cf27aSJed Brown 
15335ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15345ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15355ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1536ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1537ef19f930SBarry Smith 
1538e5c89e4eSSatish Balay   if (PetscBeganMPI) {
153999608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
154099b1327fSBarry Smith     PetscMPIInt flag;
154199b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1542e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
154399608316SBarry Smith #endif
1544e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1545e5c89e4eSSatish Balay   }
1546e5c89e4eSSatish Balay /*
1547e5c89e4eSSatish Balay 
1548e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1549e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1550e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1551e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1552e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
15530e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1554e5c89e4eSSatish Balay    memory was not freed.
1555e5c89e4eSSatish Balay 
1556e5c89e4eSSatish Balay */
15571d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1558a297a907SKarl Rupp 
1559e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1560e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
15613db9a53dSBarry Smith   PetscFunctionReturn(0);
1562e5c89e4eSSatish Balay }
1563e5c89e4eSSatish Balay 
156443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
15658cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
156643db4dbbSBarry Smith {
156743db4dbbSBarry Smith   if (*a == *b) return 1;
156843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
156943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
157043db4dbbSBarry Smith   return 0;
157143db4dbbSBarry Smith }
1572a70650f6SBarry Smith #endif
157343db4dbbSBarry Smith 
157443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
15758cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
157643db4dbbSBarry Smith {
157743db4dbbSBarry Smith   if (*a == *b) return 1;
157843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
157943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
158043db4dbbSBarry Smith   return 0;
158143db4dbbSBarry Smith }
158243db4dbbSBarry Smith #endif
1583