xref: /petsc/src/sys/objects/pinit.c (revision cd1b99a63362bb6d9ae5109b3aca94abc60d19f9)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>        /*I  "petscsys.h"   I*/
6022afb99SBarry Smith #include <petscvalgrind.h>
7665c2dedSJed Brown #include <petscviewer.h>
88101f56cSMatthew Knepley 
9a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogInitialize(void);
1195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
12a9f03627SSatish Balay #endif
13f2d66bcaSShri Abhyankar 
142d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1595c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
162d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
172d53ad75SBarry Smith #endif
182d53ad75SBarry Smith 
19a6790183SBarry Smith #if defined(PETSC_HAVE_SAWS)
2016ad0300SBarry Smith #include <petscviewersaws.h>
21a6790183SBarry Smith #endif
22e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
23e5c89e4eSSatish Balay 
2495c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
25e5c89e4eSSatish Balay 
2695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
2795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
2895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
2995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
3195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
320069ddf5SShri Abhyankar 
33e5c89e4eSSatish Balay /* user may set this BEFORE calling PetscInitialize() */
34e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
35e5c89e4eSSatish Balay 
36480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
37480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
38480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
395f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
40480cf27aSJed Brown 
41e5c89e4eSSatish Balay /*
42e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
43e5c89e4eSSatish Balay */
446a6fc655SJed Brown const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",0};
456a6fc655SJed Brown const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",0};
46e5c89e4eSSatish Balay 
47ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
48ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
490f8e0872SSatish Balay 
50a2f94806SJed Brown PetscInt PetscHotRegionDepth;
51a2f94806SJed Brown 
52b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
53b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
54b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
55b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
56b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
57b22622e2STadeu Manoel #endif
58b22622e2STadeu Manoel 
59e5c89e4eSSatish Balay /*
60e5c89e4eSSatish Balay        Checks the options database for initializations related to the
61e5c89e4eSSatish Balay     PETSc components
62e5c89e4eSSatish Balay */
6395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode  PetscOptionsCheckInitial_Components(void)
64e5c89e4eSSatish Balay {
652d747510SLisandro Dalcin   PetscBool      flg;
66e5c89e4eSSatish Balay   PetscErrorCode ierr;
67e5c89e4eSSatish Balay 
68e5c89e4eSSatish Balay   PetscFunctionBegin;
692d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
702d747510SLisandro Dalcin   if (flg) {
71e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
72e8e7597cSSatish Balay     MPI_Comm comm = PETSC_COMM_WORLD;
73e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"------Additional PETSc component options--------\n");CHKERRQ(ierr);
748e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -log_exclude: <vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
758e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -info_exclude: <null,vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
76e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
77e5c89e4eSSatish Balay #endif
78e5c89e4eSSatish Balay   }
79e5c89e4eSSatish Balay   PetscFunctionReturn(0);
80e5c89e4eSSatish Balay }
81e5c89e4eSSatish Balay 
820f11a792SBarry Smith /*
83945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8472a42c3cSBarry Smith 
8572a42c3cSBarry Smith    Collective
8672a42c3cSBarry Smith 
8772a42c3cSBarry Smith    Level: advanced
8872a42c3cSBarry Smith 
8995452b02SPatrick Sanan     Notes:
90a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
910f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
92a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
930f11a792SBarry Smith 
94a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
951ea5a559SBarry Smith 
9672a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
970f11a792SBarry Smith */
98945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
9972a42c3cSBarry Smith {
10072a42c3cSBarry Smith   PetscErrorCode ierr;
10172a42c3cSBarry Smith   int            myargc   = argc;
10272a42c3cSBarry Smith   char           **myargs = args;
10372a42c3cSBarry Smith 
10472a42c3cSBarry Smith   PetscFunctionBegin;
1053bf036e2SBarry Smith   ierr = PetscInitialize(&myargc,&myargs,filename,help);CHKERRQ(ierr);
1061ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
107df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
10872a42c3cSBarry Smith   PetscFunctionReturn(ierr);
10972a42c3cSBarry Smith }
11072a42c3cSBarry Smith 
111f0865b08SBarry Smith /*
112a56f64adSBarry Smith       Used by Julia interface to get communicator
113f0865b08SBarry Smith */
114945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
115f0865b08SBarry Smith {
116f0865b08SBarry Smith   PetscFunctionBegin;
117f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
118f0865b08SBarry Smith   PetscFunctionReturn(0);
119f0865b08SBarry Smith }
120f0865b08SBarry Smith 
121e5c89e4eSSatish Balay /*@C
122e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
123e5c89e4eSSatish Balay         the command line arguments.
124e5c89e4eSSatish Balay 
125e5c89e4eSSatish Balay    Collective
126e5c89e4eSSatish Balay 
127e5c89e4eSSatish Balay    Level: advanced
128e5c89e4eSSatish Balay 
129e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
130e5c89e4eSSatish Balay @*/
1317087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
132e5c89e4eSSatish Balay {
133e5c89e4eSSatish Balay   PetscErrorCode ierr;
134e5c89e4eSSatish Balay   int            argc   = 0;
135e5c89e4eSSatish Balay   char           **args = 0;
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay   PetscFunctionBegin;
1380298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
139e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
140e5c89e4eSSatish Balay }
141e5c89e4eSSatish Balay 
142e5c89e4eSSatish Balay /*@
143e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
144e5c89e4eSSatish Balay 
14593b6d2d1SJed Brown    Level: beginner
146e5c89e4eSSatish Balay 
147e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
148e5c89e4eSSatish Balay @*/
1497087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool  *isInitialized)
150e5c89e4eSSatish Balay {
151e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
15293b6d2d1SJed Brown   return 0;
153e5c89e4eSSatish Balay }
154e5c89e4eSSatish Balay 
155e5c89e4eSSatish Balay /*@
156e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
157e5c89e4eSSatish Balay 
158e5c89e4eSSatish Balay    Level: developer
159e5c89e4eSSatish Balay 
160e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
161e5c89e4eSSatish Balay @*/
1627087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
163e5c89e4eSSatish Balay {
164e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
16593b6d2d1SJed Brown   return 0;
166e5c89e4eSSatish Balay }
167e5c89e4eSSatish Balay 
16895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void);
169e5c89e4eSSatish Balay 
170e5c89e4eSSatish Balay /*
171e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
172e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
173e5c89e4eSSatish Balay */
174367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
175e5c89e4eSSatish Balay 
176367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
177e5c89e4eSSatish Balay {
178e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
179e5c89e4eSSatish Balay 
180e5c89e4eSSatish Balay   PetscFunctionBegin;
181e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
182e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
183e5c89e4eSSatish Balay     MPI_Abort(MPI_COMM_WORLD,1);
184e5c89e4eSSatish Balay   }
185e5c89e4eSSatish Balay 
186e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
187e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
188e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
189e5c89e4eSSatish Balay   }
190812af9f3SBarry Smith   PetscFunctionReturnVoid();
191e5c89e4eSSatish Balay }
192e5c89e4eSSatish Balay 
193e5c89e4eSSatish Balay /*
194e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
195e5c89e4eSSatish Balay sum of the second entry.
196b693b147SBarry Smith 
19776f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
198367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
199b693b147SBarry Smith there would be no place to store the both needed results.
200e5c89e4eSSatish Balay */
20176ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
202e5c89e4eSSatish Balay {
203e5c89e4eSSatish Balay   PetscErrorCode ierr;
204e5c89e4eSSatish Balay 
205e5c89e4eSSatish Balay   PetscFunctionBegin;
206d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
207d6e4c47cSJed Brown   {
208d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
209367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
210d6e4c47cSJed Brown     *max = work.max;
211d6e4c47cSJed Brown     *sum = work.sum;
212d6e4c47cSJed Brown   }
213d6e4c47cSJed Brown #else
214d6e4c47cSJed Brown   {
215d6e4c47cSJed Brown     PetscMPIInt    size,rank;
216d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
217e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
218e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
219785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
220367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2216ac3741eSJed Brown     *max = work[rank].max;
2226ac3741eSJed Brown     *sum = work[rank].sum;
223e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
224d6e4c47cSJed Brown   }
225d6e4c47cSJed Brown #endif
226e5c89e4eSSatish Balay   PetscFunctionReturn(0);
227e5c89e4eSSatish Balay }
228e5c89e4eSSatish Balay 
229e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
230e5c89e4eSSatish Balay 
231570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
23206a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
233e5c89e4eSSatish Balay 
2348cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
235e5c89e4eSSatish Balay {
236e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
237e5c89e4eSSatish Balay 
238e5c89e4eSSatish Balay   PetscFunctionBegin;
2397c2de775SJed Brown   if (*datatype == MPIU_REAL) {
240e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
241a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2427c2de775SJed Brown   }
2437c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2447c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2457c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
246a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2477c2de775SJed Brown   }
2487c2de775SJed Brown #endif
2497c2de775SJed Brown   else {
2507c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
251e2e03761SBarry Smith     MPI_Abort(MPI_COMM_WORLD,1);
252e2e03761SBarry Smith   }
253812af9f3SBarry Smith   PetscFunctionReturnVoid();
254e5c89e4eSSatish Balay }
255e5c89e4eSSatish Balay #endif
256e5c89e4eSSatish Balay 
257570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
258d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
259d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
260d9822059SBarry Smith 
2618cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
262d9822059SBarry Smith {
263d9822059SBarry Smith   PetscInt i,count = *cnt;
264d9822059SBarry Smith 
265d9822059SBarry Smith   PetscFunctionBegin;
2667c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2678c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
268a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2697c2de775SJed Brown   }
2707c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2717c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2727c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2737c2de775SJed Brown     for (i=0; i<count; i++) {
2747c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2757c2de775SJed Brown     }
2767c2de775SJed Brown   }
2777c2de775SJed Brown #endif
2787c2de775SJed Brown   else {
2797c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
2808c764dc5SJose Roman     MPI_Abort(MPI_COMM_WORLD,1);
2818c764dc5SJose Roman   }
282d9822059SBarry Smith   PetscFunctionReturnVoid();
283d9822059SBarry Smith }
284d9822059SBarry Smith 
2858cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
286d9822059SBarry Smith {
287d9822059SBarry Smith   PetscInt    i,count = *cnt;
288d9822059SBarry Smith 
289d9822059SBarry Smith   PetscFunctionBegin;
2907c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2918c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
292a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2937c2de775SJed Brown   }
2947c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2957c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2967c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2977c2de775SJed Brown     for (i=0; i<count; i++) {
2987c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2997c2de775SJed Brown     }
3007c2de775SJed Brown   }
3017c2de775SJed Brown #endif
3027c2de775SJed Brown   else {
3038c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
3048c764dc5SJose Roman     MPI_Abort(MPI_COMM_WORLD,1);
3058c764dc5SJose Roman   }
306d9822059SBarry Smith   PetscFunctionReturnVoid();
307d9822059SBarry Smith }
308d9822059SBarry Smith #endif
309d9822059SBarry Smith 
310480cf27aSJed Brown /*
311480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
312480cf27aSJed Brown 
313ff0e51ddSBarry Smith    This is called by MPI, not by users. This is called by MPI_Comm_free() when the communicator that has this  data as an attribute is freed.
314480cf27aSJed Brown 
31512801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
316480cf27aSJed Brown 
317480cf27aSJed Brown */
3188cc058d9SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
319480cf27aSJed Brown {
320480cf27aSJed Brown   PetscErrorCode ierr;
321480cf27aSJed Brown 
322480cf27aSJed Brown   PetscFunctionBegin;
32312801b39SBarry Smith   ierr = PetscInfo1(0,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
32412801b39SBarry Smith   ierr = PetscFree(count_val);CHKERRMPI(ierr);
325480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
326480cf27aSJed Brown }
327480cf27aSJed Brown 
328480cf27aSJed Brown /*
32947435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
330da3039f7SJed Brown   calls MPI_Comm_free().
331da3039f7SJed Brown 
332da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
333480cf27aSJed Brown 
334ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
335480cf27aSJed Brown 
33612801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
337480cf27aSJed Brown 
338480cf27aSJed Brown */
339da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
340480cf27aSJed Brown {
341480cf27aSJed Brown   PetscErrorCode                    ierr;
342b89831e5SBarry Smith   PetscMPIInt                       flg;
343265f3f35SJed Brown   union {MPI_Comm comm; void *ptr;} icomm,ocomm;
344480cf27aSJed Brown 
345480cf27aSJed Brown   PetscFunctionBegin;
34612801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
347ec4fadc2SJed Brown   icomm.ptr = attr_val;
348da3039f7SJed Brown 
34947435625SJed Brown   ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
35012801b39SBarry Smith   if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm");
35112801b39SBarry Smith   if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm");
35247435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
35312801b39SBarry Smith   ierr = PetscInfo1(0,"User MPI_Comm %ld is being freed after removing reference from inner PETSc comm to this outer comm\n",(long)comm);CHKERRMPI(ierr);
354da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
355b89831e5SBarry Smith }
356da3039f7SJed Brown 
357da3039f7SJed Brown /*
35847435625SJed Brown  * This is invoked on the inner comm when Petsc_DelComm_Outer calls MPI_Comm_delete_attr.  It should not be reached any other way.
359da3039f7SJed Brown  */
360da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
361da3039f7SJed Brown {
362da3039f7SJed Brown   PetscErrorCode ierr;
363da3039f7SJed Brown 
364da3039f7SJed Brown   PetscFunctionBegin;
36512801b39SBarry Smith   ierr = PetscInfo1(0,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
366480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
367480cf27aSJed Brown }
368480cf27aSJed Brown 
3695f7487a0SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Shm(MPI_Comm,PetscMPIInt,void *,void *);
37042218b76SBarry Smith 
371951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3728cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3738cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3748cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
375e39fd77fSBarry Smith #endif
376e39fd77fSBarry Smith 
37712801b39SBarry Smith PetscMPIInt PETSC_MPI_ERROR_CLASS,PETSC_MPI_ERROR_CODE;
37812801b39SBarry Smith 
379eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
380eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3816ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
3826ae9a8a6SBarry Smith char **PetscGlobalArgs = 0;
383dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
384e5c89e4eSSatish Balay 
385dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
386051e4cf2SJed Brown {
387051e4cf2SJed Brown   PetscErrorCode ierr;
388051e4cf2SJed Brown 
389051e4cf2SJed Brown   PetscFunctionBegin;
390051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
391a1601671SBarry Smith   ierr = PetscCitationsRegister("@TechReport{petsc-user-ref,\n  Author = {Satish Balay and Shrirang Abhyankar and Mark F. Adams and Jed Brown \n            and Peter Brune and Kris Buschelman and Lisandro Dalcin and\n            Victor Eijkhout and William D. Gropp and Dmitry Karpeyev and\n            Dinesh Kaushik and Matthew G. Knepley and Dave A. May and Lois Curfman McInnes\n            and Richard Tran Mills and Todd Munson and Karl Rupp and Patrick Sanan\n            and Barry F. Smith and Stefano Zampini and Hong Zhang and Hong Zhang},\n  Title = {{PETS}c Users Manual},\n  Number = {ANL-95/11 - Revision 3.11},\n  Institution = {Argonne National Laboratory},\n  Year = {2019}\n}\n",NULL);CHKERRQ(ierr);
392051e4cf2SJed 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);
393051e4cf2SJed Brown   PetscFunctionReturn(0);
394051e4cf2SJed Brown }
395e5c89e4eSSatish Balay 
3962d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3972d747510SLisandro Dalcin 
3982d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
3992d747510SLisandro Dalcin {
4002d747510SLisandro Dalcin   PetscErrorCode ierr;
4012d747510SLisandro Dalcin 
4022d747510SLisandro Dalcin   PetscFunctionBegin;
4032d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
4042d747510SLisandro Dalcin   PetscFunctionReturn(0);
4052d747510SLisandro Dalcin }
4062d747510SLisandro Dalcin 
4072d747510SLisandro Dalcin /*@C
4082d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4092d747510SLisandro Dalcin 
4102d747510SLisandro Dalcin     Not Collective
4112d747510SLisandro Dalcin 
4122d747510SLisandro Dalcin     Input Parameter:
4132d747510SLisandro Dalcin .   len - length of the string name
4142d747510SLisandro Dalcin 
4152d747510SLisandro Dalcin     Output Parameter:
4162d747510SLisandro Dalcin .   name - the name of the running program
4172d747510SLisandro Dalcin 
4182d747510SLisandro Dalcin    Level: advanced
4192d747510SLisandro Dalcin 
4202d747510SLisandro Dalcin     Notes:
4212d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4222d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4232d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4242d747510SLisandro Dalcin @*/
4252d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4262d747510SLisandro Dalcin {
4272d747510SLisandro Dalcin   PetscErrorCode ierr;
4282d747510SLisandro Dalcin 
4292d747510SLisandro Dalcin   PetscFunctionBegin;
4302d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4312d747510SLisandro Dalcin   PetscFunctionReturn(0);
4322d747510SLisandro Dalcin }
4332d747510SLisandro Dalcin 
434e5c89e4eSSatish Balay /*@C
435e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
436e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
437e5c89e4eSSatish Balay 
438e5c89e4eSSatish Balay    Not Collective
439e5c89e4eSSatish Balay 
440e5c89e4eSSatish Balay    Output Parameters:
441e5c89e4eSSatish Balay +  argc - count of number of command line arguments
442e5c89e4eSSatish Balay -  args - the command line arguments
443e5c89e4eSSatish Balay 
444e5c89e4eSSatish Balay    Level: intermediate
445e5c89e4eSSatish Balay 
446e5c89e4eSSatish Balay    Notes:
447e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
448e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
449e5c89e4eSSatish Balay 
450f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
451f177e3b1SBarry Smith 
452e5c89e4eSSatish Balay    Concepts: command line arguments
453e5c89e4eSSatish Balay 
454793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
455e5c89e4eSSatish Balay 
456e5c89e4eSSatish Balay @*/
4577087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
458e5c89e4eSSatish Balay {
459e5c89e4eSSatish Balay   PetscFunctionBegin;
46017186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
461e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
462e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
463e5c89e4eSSatish Balay   PetscFunctionReturn(0);
464e5c89e4eSSatish Balay }
465e5c89e4eSSatish Balay 
466793721a6SBarry Smith /*@C
467793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
468793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
469793721a6SBarry Smith 
470793721a6SBarry Smith    Not Collective
471793721a6SBarry Smith 
472793721a6SBarry Smith    Output Parameters:
473793721a6SBarry Smith .  args - the command line arguments
474793721a6SBarry Smith 
475793721a6SBarry Smith    Level: intermediate
476793721a6SBarry Smith 
477793721a6SBarry Smith    Notes:
478793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
479793721a6SBarry Smith 
480793721a6SBarry Smith    Concepts: command line arguments
481793721a6SBarry Smith 
482793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
483793721a6SBarry Smith 
484793721a6SBarry Smith @*/
4857087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
486793721a6SBarry Smith {
487793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
488793721a6SBarry Smith   PetscErrorCode ierr;
489793721a6SBarry Smith 
490793721a6SBarry Smith   PetscFunctionBegin;
49117186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4922d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
493785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
494793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
495793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
496793721a6SBarry Smith   }
4972d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
498793721a6SBarry Smith   PetscFunctionReturn(0);
499793721a6SBarry Smith }
500793721a6SBarry Smith 
501793721a6SBarry Smith /*@C
502793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
503793721a6SBarry Smith 
504793721a6SBarry Smith    Not Collective
505793721a6SBarry Smith 
506793721a6SBarry Smith    Output Parameters:
507793721a6SBarry Smith .  args - the command line arguments
508793721a6SBarry Smith 
509793721a6SBarry Smith    Level: intermediate
510793721a6SBarry Smith 
511793721a6SBarry Smith    Concepts: command line arguments
512793721a6SBarry Smith 
513793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
514793721a6SBarry Smith 
515793721a6SBarry Smith @*/
5167087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
517793721a6SBarry Smith {
518793721a6SBarry Smith   PetscInt       i = 0;
519793721a6SBarry Smith   PetscErrorCode ierr;
520793721a6SBarry Smith 
521793721a6SBarry Smith   PetscFunctionBegin;
522a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
523793721a6SBarry Smith   while (args[i]) {
524793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
525793721a6SBarry Smith     i++;
526793721a6SBarry Smith   }
527793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
528793721a6SBarry Smith   PetscFunctionReturn(0);
529793721a6SBarry Smith }
530793721a6SBarry Smith 
53111525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
53230befbd2SBarry Smith #include <petscconfiginfo.h>
53330befbd2SBarry Smith 
53495c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
53511525c0dSBarry Smith {
53611525c0dSBarry Smith   if (!PetscGlobalRank) {
53730befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
53811525c0dSBarry Smith     int            port;
539ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
54011525c0dSBarry Smith     size_t         applinelen,introlen;
54111525c0dSBarry Smith     PetscErrorCode ierr;
542ffbd1cfbSBarry Smith     char           sawsurl[256];
54311525c0dSBarry Smith 
544c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
54511525c0dSBarry Smith     if (flg) {
54611525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
54711525c0dSBarry Smith 
548c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
54911525c0dSBarry Smith       if (sawslog[0]) {
55011525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
55111525c0dSBarry Smith       } else {
55211525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
55311525c0dSBarry Smith       }
55411525c0dSBarry Smith     }
555c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
55611525c0dSBarry Smith     if (flg) {
55711525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
55811525c0dSBarry Smith     }
559c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
560ffbd1cfbSBarry Smith     if (selectport) {
561ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
562ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
563ffbd1cfbSBarry Smith     } else {
564c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
56511525c0dSBarry Smith       if (flg) {
56611525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
56711525c0dSBarry Smith       }
568ffbd1cfbSBarry Smith     }
569c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
57011525c0dSBarry Smith     if (flg) {
57111525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
57211525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5739c1e0ce8SBarry Smith     } else {
574c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5759c1e0ce8SBarry Smith       if (flg) {
5763c01dfcfSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
5779c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5789c1e0ce8SBarry Smith       }
57911525c0dSBarry Smith     }
580c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
58111525c0dSBarry Smith     if (flg2) {
58211525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
58311525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
58411525c0dSBarry Smith       ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr);
58511525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
58611525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
58743da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
58811525c0dSBarry Smith     }
58911525c0dSBarry Smith     ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr);
59011525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
59111525c0dSBarry Smith     introlen   = 4096 + applinelen;
59230a8c9c0SSurtai Han     applinelen += 1024;
59311525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
59411525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
59511525c0dSBarry Smith 
59611525c0dSBarry Smith     if (rootlocal) {
59711525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
59811525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
59911525c0dSBarry Smith     }
60076a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
60111525c0dSBarry Smith     if (rootlocal && help) {
602928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running <a href=\"%s.c.html\">%s</a> %s</center><br><center><pre>%s</pre></center><br>\n",programname,programname,options,help);CHKERRQ(ierr);
60311525c0dSBarry Smith     } else if (help) {
604928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
60511525c0dSBarry Smith     } else {
606928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
60711525c0dSBarry Smith     }
608b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
60930befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
61011525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
611a8d69d7bSBarry Smith                                     "<center><h2> <a href=\"https://www.mcs.anl.gov/petsc\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n"
612df62c222SBarry 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"
613928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
61443da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
61511525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
61611525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
617ffbd1cfbSBarry Smith     if (selectport) {
618aa573868SBarry Smith       PetscBool silent;
6197d812c46SBarry Smith 
6207d812c46SBarry Smith       ierr = SAWs_Initialize();
6217d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6227d812c46SBarry Smith       while (ierr) {
6237d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6247d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6257d812c46SBarry Smith         ierr = SAWs_Initialize();
6267d812c46SBarry Smith       }
6277d812c46SBarry Smith 
628aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
629aa573868SBarry Smith       if (!silent) {
630ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
631ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
632ffbd1cfbSBarry Smith       }
6337d812c46SBarry Smith     } else {
6347d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
635aa573868SBarry Smith     }
6360af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6370af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6380af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6390af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6400af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
64111525c0dSBarry Smith   }
64211525c0dSBarry Smith   PetscFunctionReturn(0);
64311525c0dSBarry Smith }
64411525c0dSBarry Smith #endif
64511525c0dSBarry Smith 
646a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
647a56f64adSBarry Smith #include <adios.h>
64822580e64SBarry Smith #include <adios_read.h>
6497b56e58cSSatish Balay int64_t Petsc_adios_group;
650a56f64adSBarry Smith #endif
65155d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
65255d657eeSBarry Smith #include <adios2_c.h>
65355d657eeSBarry Smith #endif
654*cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
655*cd1b99a6SBarry Smith #include <omp.h>
656*cd1b99a6SBarry Smith #endif
657a56f64adSBarry Smith 
658e5c89e4eSSatish Balay /*@C
659e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
660e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
661e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
662e5c89e4eSSatish Balay    your program -- usually the very first line!
663e5c89e4eSSatish Balay 
664e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
665e5c89e4eSSatish Balay 
666e5c89e4eSSatish Balay    Input Parameters:
667e5c89e4eSSatish Balay +  argc - count of number of command line arguments
668e5c89e4eSSatish Balay .  args - the command line arguments
6690298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
670fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6710298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
672e5c89e4eSSatish Balay 
67305827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
67405827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
67505827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
67605827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
67705827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
678e5c89e4eSSatish Balay 
679e5c89e4eSSatish Balay    Options Database Keys:
6807ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6817ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
682e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
6838a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
6848a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
6858a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
6868a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
6878a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
688e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
689e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
690e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
691e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
6922fb0ec9aSBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries)
693e5c89e4eSSatish Balay .  -malloc no - Indicates not to use error-checking malloc
6942fb0ec9aSBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free
695aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
696dc92acbaSJed Brown .  -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds
697e5c89e4eSSatish Balay .  -fp_trap - Stops on floating point exceptions (Note that on the
698e5c89e4eSSatish Balay               IBM RS6000 this slows code by at least a factor of 10.)
699e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
700e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
701e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
702e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
703e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7040841954dSBarry Smith -  -memory_view - Print memory usage at end of run
705e5c89e4eSSatish Balay 
706e5c89e4eSSatish Balay    Options Database Keys for Profiling:
707a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
708495fc317SBarry Smith +  -info <optional filename> - Prints verbose information to the screen
709495fc317SBarry Smith .  -info_exclude <null,vec,mat,pc,ksp,snes,ts> - Excludes some of the verbose messages
710217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
711217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
712495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
713e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7149a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
7159a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
716495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
717f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
718495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
719495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
720c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
72187c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
722c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
723495fc317SBarry Smith 
724609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
725e5c89e4eSSatish Balay 
726ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
727ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
728ffbd1cfbSBarry 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
729ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
730ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
731ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
732ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
733ffbd1cfbSBarry Smith 
734e5c89e4eSSatish Balay    Environmental Variables:
735e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
736e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
737e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
738e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
739e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
740e5c89e4eSSatish Balay 
741e5c89e4eSSatish Balay 
742e5c89e4eSSatish Balay    Level: beginner
743e5c89e4eSSatish Balay 
744e5c89e4eSSatish Balay    Notes:
745e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
746e5c89e4eSSatish Balay    it before PetscInitialize().
747e5c89e4eSSatish Balay 
748e5c89e4eSSatish Balay    Fortran Version:
749e5c89e4eSSatish Balay    In Fortran this routine has the format
750e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
751e5c89e4eSSatish Balay 
752e5c89e4eSSatish Balay +   ierr - error return code
7530eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
754fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
755e5c89e4eSSatish Balay 
756e5c89e4eSSatish Balay    Important Fortran Note:
7570eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7580298fd71SBarry Smith    null character string; you CANNOT just use NULL as
759a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
760e5c89e4eSSatish Balay 
76101cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
76201cb0274SBarry Smith    calling PetscInitialize().
763e5c89e4eSSatish Balay 
764e5c89e4eSSatish Balay    Concepts: initializing PETSc
765e5c89e4eSSatish Balay 
76601cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
767e5c89e4eSSatish Balay 
768e5c89e4eSSatish Balay @*/
7697087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
770e5c89e4eSSatish Balay {
771e5c89e4eSSatish Balay   PetscErrorCode ierr;
7724bb5149bSJed Brown   PetscMPIInt    flag, size;
77319c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
774e5c89e4eSSatish Balay   char           hostname[256];
775e5c89e4eSSatish Balay 
776e5c89e4eSSatish Balay   PetscFunctionBegin;
777e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7783d96e996SBarry Smith   /*
7793d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7803d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7813d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7823d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7833d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7843d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7853d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7863d96e996SBarry Smith       listed above and since that time are compatible.
7873d96e996SBarry Smith 
7883d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
7893d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
7903d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
7913d96e996SBarry Smith       to perform the checking.
7923d96e996SBarry Smith 
7933d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
7943d96e996SBarry Smith 
7953d96e996SBarry Smith       Questions:
7963d96e996SBarry Smith 
7973d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
7983d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
7993d96e996SBarry Smith   */
8003d96e996SBarry Smith 
80119c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
80219c5658aSBarry Smith   {
80319c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
80419c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
80519c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8063d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
80719c5658aSBarry Smith #if defined(MPICH_VERSION)
8083d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
80919c5658aSBarry Smith     {
81019c5658aSBarry Smith       char *ver,*lf;
81119c5658aSBarry Smith       flg = PETSC_FALSE;
81219c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
81319c5658aSBarry Smith       if (ver) {
81419c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
81519c5658aSBarry Smith         if (lf) {
81619c5658aSBarry Smith           *lf = 0;
81719c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
81819c5658aSBarry Smith         }
81919c5658aSBarry Smith       }
82019c5658aSBarry Smith       if (!flg) {
82119c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8223d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
82319c5658aSBarry Smith       }
82419c5658aSBarry Smith     }
8253d96e996SBarry Smith #endif
8263d96e996SBarry 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?) */
82719c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
82819c5658aSBarry Smith     {
82919c5658aSBarry Smith       char *ver,bs[32],*bsf;
83019c5658aSBarry Smith       flg = PETSC_FALSE;
83119c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
83219c5658aSBarry Smith       if (ver) {
8332e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
83419c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
83519c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
83619c5658aSBarry Smith       }
83719c5658aSBarry Smith       if (!flg) {
83819c5658aSBarry 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);
8393d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
84019c5658aSBarry Smith       }
84119c5658aSBarry Smith     }
84219c5658aSBarry Smith #endif
84319c5658aSBarry Smith   }
84419c5658aSBarry Smith #endif
84519c5658aSBarry Smith 
846e5c89e4eSSatish Balay 
847ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
848d89683f4Sbcordonn   PETSC_STDOUT = stdout;
849ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
850e5c89e4eSSatish Balay 
8510c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8520c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8530c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8540c30907bSSatish Balay #endif
8550c30907bSSatish Balay 
8564416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
857e5c89e4eSSatish Balay 
858e5c89e4eSSatish Balay   /*
859e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
860e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
861e5c89e4eSSatish Balay   */
862e5c89e4eSSatish Balay   if (argc && *argc) {
863e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
864e5c89e4eSSatish Balay   } else {
865e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
866e5c89e4eSSatish Balay   }
867e5c89e4eSSatish Balay 
868e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
869e5c89e4eSSatish Balay   if (!flag) {
870e32f2f54SBarry 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");
8715e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8725e765c61SJed Brown     {
8735e765c61SJed Brown       PetscMPIInt provided;
8745e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
8755e765c61SJed Brown     }
8765e765c61SJed Brown #else
877e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
8785e765c61SJed Brown #endif
879e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
880e5c89e4eSSatish Balay   }
881e5c89e4eSSatish Balay   if (argc && args) {
882e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
883e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
884e5c89e4eSSatish Balay   }
885e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
8865ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
8875ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
8885ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
889ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
890e5c89e4eSSatish Balay 
891a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
892d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
893e5c89e4eSSatish Balay 
89412801b39SBarry Smith   ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
89512801b39SBarry Smith   ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
89612801b39SBarry Smith 
897e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
898e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
899e5c89e4eSSatish Balay 
900e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
901e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
902e5c89e4eSSatish Balay 
9038ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9048ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9057cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
906e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
907e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
908e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
909e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
910e316c87fSJed Brown #endif
911e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9128ad47952SJed Brown 
913e5c89e4eSSatish Balay   /*
914e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
915e5c89e4eSSatish Balay      shared libraries the constructors for global variables
916e5c89e4eSSatish Balay      are not called; at least on IRIX.
917e5c89e4eSSatish Balay   */
918886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
919e5c89e4eSSatish Balay   {
920252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
92150f81f78SJed Brown     PetscComplex ic(0.0,1.0);
922e5c89e4eSSatish Balay     PETSC_i = ic;
923252ecd31SSatish Balay #else
92450f81f78SJed Brown     PETSC_i = _Complex_I;
925b7940d39SSatish Balay #endif
926762437b8SSatish Balay   }
927762437b8SSatish Balay 
9282c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
929e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
930500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
931500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
932500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9332c876bd9SBarry Smith #endif
934886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
935e5c89e4eSSatish Balay 
936e5c89e4eSSatish Balay   /*
937e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
938e5c89e4eSSatish Balay      half of the entries and maxes the second half.
939e5c89e4eSSatish Balay   */
940367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
941e5c89e4eSSatish Balay 
942ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
943c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
944c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9457c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9468c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9478c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9488c764dc5SJose Roman #endif
949d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
950d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
951570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
952570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
953570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
954570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
955570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
956c90a1750SBarry Smith #endif
957c90a1750SBarry Smith 
958570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
959cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
960cca4cb22SSatish Balay #endif
961cca4cb22SSatish Balay 
962e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
963e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
964e5c89e4eSSatish Balay 
96540df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
966e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
967e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
96844041f26SJed Brown #endif
969e5c89e4eSSatish Balay 
970ec957eceSBarry Smith 
971e5c89e4eSSatish Balay   /*
972480cf27aSJed Brown      Attributes to be set on PETSc communicators
973480cf27aSJed Brown   */
97412801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
97512801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
97612801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
9775f7487a0SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Shm,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
978480cf27aSJed Brown 
979480cf27aSJed Brown   /*
980e8fb0fc0SBarry Smith      Build the options database
981e5c89e4eSSatish Balay   */
982c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
983e5c89e4eSSatish Balay 
984703f3eceSBarry Smith   /* call a second time so it can look in the options database */
985703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
9866dc8fec2Sbcordonn 
987e5c89e4eSSatish Balay   /*
988e5c89e4eSSatish Balay      Print main application help message
989e5c89e4eSSatish Balay   */
9902d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
991e5c89e4eSSatish Balay   if (help && flg) {
992e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
993e5c89e4eSSatish Balay   }
994e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
995e5c89e4eSSatish Balay 
996d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
997d45a07a7SBarry Smith 
998e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
99911525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1000f4202a44SBarry Smith #endif
1001f4202a44SBarry Smith 
1002896238b9SBarry Smith   /* Creates the logging data structures; this is enabled even if logging is not turned on */
1003a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1004896238b9SBarry Smith   ierr = PetscLogInitialize();CHKERRQ(ierr);
1005a9f03627SSatish Balay #endif
1006e5c89e4eSSatish Balay 
1007e5c89e4eSSatish Balay   /*
1008e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1009e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1010e5c89e4eSSatish Balay   */
1011e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1012e5c89e4eSSatish Balay 
1013e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
1014ae15b995SBarry Smith   ierr = PetscInfo1(0,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1015e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
1016ae15b995SBarry Smith   ierr = PetscInfo1(0,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1017*cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1018*cd1b99a6SBarry Smith   if (PetscLogPrintInfo) {
1019*cd1b99a6SBarry Smith     char *threads = getenv("OMP_NUM_THREADS");
1020e5c89e4eSSatish Balay 
1021*cd1b99a6SBarry Smith     if (threads) {
1022*cd1b99a6SBarry Smith       ierr = PetscInfo1(0,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);
1023*cd1b99a6SBarry Smith     } else {
1024*cd1b99a6SBarry Smith #define NMAX  10000
1025*cd1b99a6SBarry Smith       int         i,nth;
1026*cd1b99a6SBarry Smith       PetscScalar *x;
1027*cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1028*cd1b99a6SBarry Smith #pragma omp parallel for
1029*cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1030*cd1b99a6SBarry Smith         x[i] = 0.0;
1031*cd1b99a6SBarry Smith         nth  = omp_get_num_threads();
1032*cd1b99a6SBarry Smith       }
1033*cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
1034*cd1b99a6SBarry Smith       ierr = PetscInfo1(0,"Number of OpenMP threads %d (number not set with OMP_NUM_THREADS)\n",nth);
1035*cd1b99a6SBarry Smith     }
1036*cd1b99a6SBarry Smith   }
1037*cd1b99a6SBarry Smith #endif
1038e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr);
1039ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1040c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1041ef6c6fedSBoyana Norris 
1042951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1043e39fd77fSBarry Smith   /*
1044e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1045e39fd77fSBarry Smith 
1046e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1047e39fd77fSBarry Smith   */
104830815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
10490298fd71SBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
105030815ce0SLisandro Dalcin   }
1051951e3c8eSBarry Smith #endif
1052e39fd77fSBarry Smith 
105341c0b4b3SShri Abhyankar   /*
105441c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
105541c0b4b3SShri Abhyankar   */
1056ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1057e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1058e1167bb9SShri Abhyankar #endif
1059e1167bb9SShri Abhyankar 
10602d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10612d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10622d53ad75SBarry Smith #endif
10632d53ad75SBarry Smith 
10645e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
106587c3beb6SLisandro Dalcin   {
106687c3beb6SLisandro Dalcin     PetscViewer viewer;
106722e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
10685e71baefSBarry Smith     if (flg) {
10695e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
107087c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
107187c3beb6SLisandro Dalcin     }
10725e71baefSBarry Smith   }
10735e71baefSBarry Smith #endif
1074dff31646SBarry Smith 
107587c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
107687c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
107787c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
107887c3beb6SLisandro Dalcin 
1079a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1080a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
10817b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1082a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
10839fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1084a56f64adSBarry Smith #endif
108555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
108655d657eeSBarry Smith #endif
1087a56f64adSBarry Smith 
1088301d30feSBarry Smith   /*
1089301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1090301d30feSBarry Smith   */
1091301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
10922db0e300SLisandro Dalcin 
10932db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
10942db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1095301d30feSBarry Smith   PetscFunctionReturn(0);
1096e5c89e4eSSatish Balay }
1097e5c89e4eSSatish Balay 
10984097062eSBarry Smith #if defined(PETSC_USE_LOG)
109995c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1100ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1101ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
110295c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11034097062eSBarry Smith #endif
1104e5c89e4eSSatish Balay 
1105008a6e76SBarry Smith /*
1106008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1107008a6e76SBarry Smith */
1108008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1109008a6e76SBarry Smith {
1110008a6e76SBarry Smith   PetscErrorCode ierr;
1111008a6e76SBarry Smith 
1112008a6e76SBarry Smith   PetscFunctionBegin;
1113008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1114008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1115008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1116008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1117008a6e76SBarry Smith #endif
1118008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1119008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1120008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1121008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1122008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1123008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1124008a6e76SBarry Smith #endif
1125008a6e76SBarry Smith 
1126008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1127008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1128008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1129008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1130008a6e76SBarry Smith #endif
1131008a6e76SBarry Smith #endif
1132008a6e76SBarry Smith 
1133008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1134008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1135008a6e76SBarry Smith #endif
1136008a6e76SBarry Smith 
1137008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
113840df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1139008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1140008a6e76SBarry Smith #endif
1141008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1142008a6e76SBarry Smith   PetscFunctionReturn(0);
1143008a6e76SBarry Smith }
1144008a6e76SBarry Smith 
1145e5c89e4eSSatish Balay /*@C
1146e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1147e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1148e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1149e5c89e4eSSatish Balay 
1150e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1151e5c89e4eSSatish Balay 
1152e5c89e4eSSatish Balay    Options Database Keys:
115326a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1154e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11557eb1d149SBarry 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
1156e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
1157e5c89e4eSSatish Balay .  -malloc_dump - Calls PetscMallocDump()
1158e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
1159e5c89e4eSSatish Balay -  -malloc_log - Prints summary of memory usage
1160e5c89e4eSSatish Balay 
1161e5c89e4eSSatish Balay    Level: beginner
1162e5c89e4eSSatish Balay 
1163e5c89e4eSSatish Balay    Note:
1164e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1165e5c89e4eSSatish Balay 
116688c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1167e5c89e4eSSatish Balay @*/
11687087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1169e5c89e4eSSatish Balay {
1170e5c89e4eSSatish Balay   PetscErrorCode ierr;
11714bb5149bSJed Brown   PetscMPIInt    rank;
1172a8d2bbe5SBarry Smith   PetscInt       nopt;
11732bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1174dff31646SBarry Smith   PetscBool      flg;
117510463e74SBarry Smith #if defined(PETSC_USE_LOG)
117610463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
117710463e74SBarry Smith #endif
1178e5c89e4eSSatish Balay 
1179e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
11804b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
11813cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1182e5c89e4eSSatish Balay   }
11833cbbc5ffSBarry Smith   PetscFunctionBegin;
11840298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1185b022a5c1SBarry Smith 
11861f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1187a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
118822580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1189a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1190a56f64adSBarry Smith #endif
119155d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
119255d657eeSBarry Smith #endif
1193c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1194dff31646SBarry Smith   if (flg) {
11951f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
11961f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
11971f817a21SBarry Smith 
1198c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11991f817a21SBarry Smith     if (filename[0]) {
12001f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
12011f817a21SBarry Smith     }
1202dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1203dff31646SBarry Smith     cits[0] = 0;
1204dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
12051f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
12061f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12071f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
12081f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12091f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1210dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1211dff31646SBarry Smith   }
1212dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1213dff31646SBarry Smith 
1214c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
121504102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
121604102261SBarry Smith   {
121704102261SBarry Smith     PetscInt nmax = 2;
121804102261SBarry Smith     char     **buffs;
121904102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1220c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
122104102261SBarry Smith     if (flg1) {
122204102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
122304102261SBarry Smith       if (nmax == 1) {
122404102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
122504102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
122604102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
122704102261SBarry Smith       }
122804102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
122904102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
123004102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
123104102261SBarry Smith     }
123204102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
123304102261SBarry Smith   }
1234763ec7b1SBarry Smith   {
1235763ec7b1SBarry Smith     PetscInt nmax = 2;
1236763ec7b1SBarry Smith     char     **buffs;
1237763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1238763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1239763ec7b1SBarry Smith     if (flg1) {
1240763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1241763ec7b1SBarry Smith       if (nmax == 1) {
1242763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1243763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1244763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1245763ec7b1SBarry Smith       }
1246763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1247763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1248763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1249763ec7b1SBarry Smith     }
1250763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1251763ec7b1SBarry Smith   }
125204102261SBarry Smith #endif
125367234432SDmitry Karpeev   /*
125467234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
125567234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
125667234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
125767234432SDmitry Karpeev   */
125867234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
125904102261SBarry Smith 
12602d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12612d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12622d53ad75SBarry Smith #endif
12632d53ad75SBarry Smith 
12642d53ad75SBarry Smith 
1265e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1266dff31646SBarry Smith   flg = PETSC_FALSE;
1267c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1268d5649816SBarry Smith   if (flg) {
1269e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1270d5649816SBarry Smith   }
1271d5649816SBarry Smith #endif
1272d5649816SBarry Smith 
1273681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1274681455b2SBarry Smith   flg1 = PETSC_FALSE;
1275c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1276681455b2SBarry Smith   if (flg1) {
1277681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1278681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1279681455b2SBarry Smith   }
1280681455b2SBarry Smith #endif
1281681455b2SBarry Smith 
128267584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1283c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1284e5c89e4eSSatish Balay   if (!flg2) {
128590d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1286c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1287e5c89e4eSSatish Balay   }
1288e5c89e4eSSatish Balay   if (flg2) {
12890841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1290e5c89e4eSSatish Balay   }
129167584ceeSBarry Smith #endif
1292e5c89e4eSSatish Balay 
1293e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
129490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1295c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1296e5c89e4eSSatish Balay   if (flg1) {
1297e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1298205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1299e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1300e5c89e4eSSatish Balay   }
1301e5c89e4eSSatish Balay #endif
1302e5c89e4eSSatish Balay 
1303e5c89e4eSSatish Balay 
1304e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1305e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1306e5c89e4eSSatish Balay   mname[0] = 0;
1307c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1308e5c89e4eSSatish Balay   if (flg1) {
1309e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1310e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1311e5c89e4eSSatish Balay   }
1312e5c89e4eSSatish Balay #endif
1313356e5837SBarry Smith #endif
1314a297a907SKarl Rupp 
1315dd710f27SBarry Smith   /*
1316dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1317dd710f27SBarry Smith   */
1318dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1319dd710f27SBarry Smith 
1320356e5837SBarry Smith #if defined(PETSC_USE_LOG)
132187c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1322f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
132387c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
132487c3beb6SLisandro Dalcin 
1325356e5837SBarry Smith   mname[0] = 0;
1326c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1327e5c89e4eSSatish Balay   if (flg1) {
132891eabc43SBarry Smith     PetscViewer viewer;
132920a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
133091eabc43SBarry Smith     if (mname[0]) {
133191eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
133291eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13336bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
133433f85c2fSBarry Smith     } else {
133533f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13369a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
133733f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13389a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
133933f85c2fSBarry Smith     }
1340e5c89e4eSSatish Balay   }
1341a297a907SKarl Rupp 
1342dd710f27SBarry Smith   /*
1343dd710f27SBarry Smith      Free any objects created by the last block of code.
1344dd710f27SBarry Smith   */
1345dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1346dd710f27SBarry Smith 
1347dd710f27SBarry Smith   mname[0] = 0;
1348c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1349c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
13507ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1351e5c89e4eSSatish Balay #endif
135210463e74SBarry Smith 
135333f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
135410463e74SBarry Smith 
135590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1356c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1357e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
135890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1359c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1360e5c89e4eSSatish Balay   if (flg1) {
1361e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1362e5c89e4eSSatish Balay   }
136390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
136490d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13658bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1366c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1367c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1368c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1369e4c476e2SSatish Balay 
1370e5c89e4eSSatish Balay   if (flg2) {
1371be56827dSJed Brown     PetscViewer viewer;
137202ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
137302ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1374c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1375be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1376e5c89e4eSSatish Balay   }
1377e5c89e4eSSatish Balay 
1378e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1379c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1380c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1381e5c89e4eSSatish Balay 
138233fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1383c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
1384c5929fdfSBarry Smith   ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
1385e5c89e4eSSatish Balay   if (flg3) {
1386e5c89e4eSSatish Balay     if (!flg2) { /* have not yet printed the options */
1387be56827dSJed Brown       PetscViewer viewer;
138802ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
138902ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1390c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1391be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1392e5c89e4eSSatish Balay     }
1393e5c89e4eSSatish Balay     if (!nopt) {
1394e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1395e5c89e4eSSatish Balay     } else if (nopt == 1) {
1396e5c89e4eSSatish Balay       ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1397e5c89e4eSSatish Balay     } else {
13987582186dSLisandro Dalcin       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1399e5c89e4eSSatish Balay     }
1400df12ba86SBarry Smith   }
1401e5c89e4eSSatish Balay #if defined(PETSC_USE_DEBUG)
1402da8b8a77SBarry Smith   if (nopt && !flg3 && !flg1) {
1403e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
1404e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
1405c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1406e5c89e4eSSatish Balay   } else if (nopt && flg3) {
1407e5c89e4eSSatish Balay #else
1408e5c89e4eSSatish Balay   if (nopt && flg3) {
1409e5c89e4eSSatish Balay #endif
1410c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1411e5c89e4eSSatish Balay   }
1412e5c89e4eSSatish Balay 
1413e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1414d45a07a7SBarry Smith   if (!PetscGlobalRank) {
141587f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
141616ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1417d45a07a7SBarry Smith   }
1418ec957eceSBarry Smith #endif
1419ec957eceSBarry Smith 
14204097062eSBarry Smith #if defined(PETSC_USE_LOG)
142110463e74SBarry Smith   /*
1422dbc8283eSBarry Smith        List all objects the user may have forgot to free
14232eff7a51SBarry Smith   */
142405df10baSBarry Smith   if (PetscObjectsLog) {
1425c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1426a64a8e02SBarry Smith     if (flg1) {
1427a64a8e02SBarry Smith       MPI_Comm local_comm;
14287eb1d149SBarry Smith       char     string[64];
1429a64a8e02SBarry Smith 
1430c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1431a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1432a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14337eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1434a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1435a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14360a1571b3SBarry Smith     }
143705df10baSBarry Smith   }
14384097062eSBarry Smith #endif
14394097062eSBarry Smith 
14404097062eSBarry Smith #if defined(PETSC_USE_LOG)
1441dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1442dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1443a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14444097062eSBarry Smith #endif
14452eff7a51SBarry Smith 
144633f85c2fSBarry Smith   /*
144733f85c2fSBarry Smith      Destroy any packages that registered a finalize
144833f85c2fSBarry Smith   */
144933f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
145033f85c2fSBarry Smith 
1451101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1452fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1453101409b8SToby Isaac #endif
1454101409b8SToby Isaac 
145533f85c2fSBarry Smith   /*
145648dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
145748dd1dffSBarry Smith 
145837e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
145948dd1dffSBarry Smith   */
146037e93019SBarry Smith 
14614028d114SSatish Balay   if (petsc_history) {
1462f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
1463e5c89e4eSSatish Balay     petsc_history = 0;
1464e5c89e4eSSatish Balay   }
14659de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1466e5c89e4eSSatish Balay 
14670298fd71SBarry Smith   ierr = PetscInfoAllow(PETSC_FALSE,NULL);CHKERRQ(ierr);
1468e5c89e4eSSatish Balay 
146967584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
14708bb29257SSatish Balay   {
1471e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
1472e5c89e4eSSatish Balay     FILE *fd;
1473ed9cf6e9SBarry Smith     int  err;
1474e5c89e4eSSatish Balay 
1475e5c89e4eSSatish Balay     fname[0] = 0;
1476a297a907SKarl Rupp 
1477c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1478dc92acbaSJed Brown     flg2 = PETSC_FALSE;
1479c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);
14808bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
1481dc92acbaSJed Brown     if (PETSC_RUNNING_ON_VALGRIND) flg2 = PETSC_FALSE;
1482dc92acbaSJed Brown #else
1483dc92acbaSJed Brown     flg2 = PETSC_FALSE;         /* Skip reporting for optimized builds regardless of -malloc_test */
1484dc92acbaSJed Brown #endif
1485e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1486e5c89e4eSSatish Balay       char sname[PETSC_MAX_PATH_LEN];
1487e5c89e4eSSatish Balay 
14882e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1489e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1490e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1491ed9cf6e9SBarry Smith       err  = fclose(fd);
1492e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1493dc92acbaSJed Brown     } else if (flg1 || flg2) {
1494e5c89e4eSSatish Balay       MPI_Comm local_comm;
1495e5c89e4eSSatish Balay 
1496e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1497e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1498e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1499e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1500e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1501e5c89e4eSSatish Balay     }
1502e5c89e4eSSatish Balay   }
1503a64a8e02SBarry Smith 
15048bb29257SSatish Balay   {
1505e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
15060298fd71SBarry Smith     FILE *fd = NULL;
1507e5c89e4eSSatish Balay 
1508e5c89e4eSSatish Balay     fname[0] = 0;
1509a297a907SKarl Rupp 
1510c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_log",fname,250,&flg1);CHKERRQ(ierr);
1511c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-malloc_log_threshold",&flg2);CHKERRQ(ierr);
1512e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1513ed9cf6e9SBarry Smith       int err;
1514e5c89e4eSSatish Balay 
1515574034a9SJed Brown       if (!rank) {
1516574034a9SJed Brown         fd = fopen(fname,"w");
1517574034a9SJed Brown         if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",fname);
1518574034a9SJed Brown       }
1519e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(fd);CHKERRQ(ierr);
1520574034a9SJed Brown       if (fd) {
1521ed9cf6e9SBarry Smith         err = fclose(fd);
1522e32f2f54SBarry Smith         if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1523574034a9SJed Brown       }
1524574034a9SJed Brown     } else if (flg1 || flg2) {
1525e5c89e4eSSatish Balay       ierr = PetscMallocDumpLog(stdout);CHKERRQ(ierr);
1526e5c89e4eSSatish Balay     }
1527e5c89e4eSSatish Balay   }
152867584ceeSBarry Smith #endif
152920e2c332SMatthew G. Knepley 
15305486ca60SMatthew G. Knepley   /*
15315486ca60SMatthew G. Knepley      Close any open dynamic libraries
15325486ca60SMatthew G. Knepley   */
15335486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15345486ca60SMatthew G. Knepley 
1535e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15364416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1537e5c89e4eSSatish Balay 
1538e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
1539e5c89e4eSSatish Balay   PetscGlobalArgs = 0;
1540e5c89e4eSSatish Balay 
1541008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1542e5c89e4eSSatish Balay 
1543dbc8283eSBarry Smith   /*
1544efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1545efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1546efb80d3cSBarry Smith 
1547efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1548efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1549dbc8283eSBarry Smith  */
1550b770b1f6SSatish Balay   {
1551dbc8283eSBarry Smith     PetscCommCounter *counter;
1552dbc8283eSBarry Smith     PetscMPIInt      flg;
1553dbc8283eSBarry Smith     MPI_Comm         icomm;
1554265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
155547435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1556dbc8283eSBarry Smith     if (flg) {
1557265f3f35SJed Brown       icomm = ucomm.comm;
155847435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1559dbc8283eSBarry 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");
1560dbc8283eSBarry Smith 
156147435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
156247435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1563efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1564dbc8283eSBarry Smith     }
156547435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1566dbc8283eSBarry Smith     if (flg) {
1567265f3f35SJed Brown       icomm = ucomm.comm;
156847435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1569dbc8283eSBarry 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");
1570dbc8283eSBarry Smith 
157147435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
157247435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1573efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1574dbc8283eSBarry Smith     }
1575b770b1f6SSatish Balay   }
1576dbc8283eSBarry Smith 
157747435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
157847435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
157947435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
15805f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1581480cf27aSJed Brown 
15825ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15835ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15845ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1585ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1586ef19f930SBarry Smith 
1587e5c89e4eSSatish Balay   if (PetscBeganMPI) {
158899608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
158999b1327fSBarry Smith     PetscMPIInt flag;
159099b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1591e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
159299608316SBarry Smith #endif
1593e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1594e5c89e4eSSatish Balay   }
1595e5c89e4eSSatish Balay /*
1596e5c89e4eSSatish Balay 
1597e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1598e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1599e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1600e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1601e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16020e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1603e5c89e4eSSatish Balay    memory was not freed.
1604e5c89e4eSSatish Balay 
1605e5c89e4eSSatish Balay */
16061d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1607a297a907SKarl Rupp 
1608e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1609e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
16103db9a53dSBarry Smith   PetscFunctionReturn(0);
1611e5c89e4eSSatish Balay }
1612e5c89e4eSSatish Balay 
161343db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16148cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
161543db4dbbSBarry Smith {
161643db4dbbSBarry Smith   if (*a == *b) return 1;
161743db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
161843db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
161943db4dbbSBarry Smith   return 0;
162043db4dbbSBarry Smith }
1621a70650f6SBarry Smith #endif
162243db4dbbSBarry Smith 
162343db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16248cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
162543db4dbbSBarry Smith {
162643db4dbbSBarry Smith   if (*a == *b) return 1;
162743db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
162843db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
162943db4dbbSBarry Smith   return 0;
163043db4dbbSBarry Smith }
163143db4dbbSBarry Smith #endif
1632