xref: /petsc/src/sys/objects/pinit.c (revision e94e781be4d0de67afa8d29cbcd676556dbc0369)
1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>        /*I  "petscsys.h"   I*/
6022afb99SBarry Smith #include <petscvalgrind.h>
7665c2dedSJed Brown #include <petscviewer.h>
88101f56cSMatthew Knepley 
9a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
11a9f03627SSatish Balay #endif
12f2d66bcaSShri Abhyankar 
132d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1495c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
152d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
162d53ad75SBarry Smith #endif
172d53ad75SBarry Smith 
18a6790183SBarry Smith #if defined(PETSC_HAVE_SAWS)
1916ad0300SBarry Smith #include <petscviewersaws.h>
20a6790183SBarry Smith #endif
21eb02082bSJunchao Zhang 
22eb02082bSJunchao Zhang #if defined(PETSC_HAVE_OPENMPI_MAJOR)
23eb02082bSJunchao Zhang #include "mpi-ext.h" /* Needed for CUDA-aware check */
24eb02082bSJunchao Zhang #endif
25e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
26e5c89e4eSSatish Balay 
2795c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
28e5c89e4eSSatish Balay 
2995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
3095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
3195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
3295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
3495c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
350069ddf5SShri Abhyankar 
36e5c89e4eSSatish Balay /* user may set this BEFORE calling PetscInitialize() */
37e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
38e5c89e4eSSatish Balay 
39480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
40480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
41480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
425f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
43480cf27aSJed Brown 
44eb02082bSJunchao Zhang /* Do not need put this in a guard like PETSC_HAVE_CUDA. Without configuring PETSc --with-cuda, users can still use option -use_gpu_aware_mpi */
45eb02082bSJunchao Zhang PetscBool use_gpu_aware_mpi = PETSC_FALSE;
46eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
47eb02082bSJunchao Zhang PetscBool sf_use_default_cuda_stream = PETSC_FALSE;
48eb02082bSJunchao Zhang #endif
49eb02082bSJunchao Zhang 
50e5c89e4eSSatish Balay /*
51e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
52e5c89e4eSSatish Balay */
5302c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
5402c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
55e5c89e4eSSatish Balay 
56ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
57ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
580f8e0872SSatish Balay 
59a2f94806SJed Brown PetscInt PetscHotRegionDepth;
60a2f94806SJed Brown 
61b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
62b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
63b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
64b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
65b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
66b22622e2STadeu Manoel #endif
67b22622e2STadeu Manoel 
68e5c89e4eSSatish Balay /*
69e5c89e4eSSatish Balay        Checks the options database for initializations related to the
70e5c89e4eSSatish Balay     PETSc components
71e5c89e4eSSatish Balay */
7295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode  PetscOptionsCheckInitial_Components(void)
73e5c89e4eSSatish Balay {
742d747510SLisandro Dalcin   PetscBool      flg;
75e5c89e4eSSatish Balay   PetscErrorCode ierr;
76e5c89e4eSSatish Balay 
77e5c89e4eSSatish Balay   PetscFunctionBegin;
782d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
792d747510SLisandro Dalcin   if (flg) {
80e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
81e8e7597cSSatish Balay     MPI_Comm comm = PETSC_COMM_WORLD;
82e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"------Additional PETSc component options--------\n");CHKERRQ(ierr);
838e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -log_exclude: <vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
84e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
85e5c89e4eSSatish Balay #endif
86e5c89e4eSSatish Balay   }
87e5c89e4eSSatish Balay   PetscFunctionReturn(0);
88e5c89e4eSSatish Balay }
89e5c89e4eSSatish Balay 
900f11a792SBarry Smith /*
91945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
9272a42c3cSBarry Smith 
9372a42c3cSBarry Smith    Collective
9472a42c3cSBarry Smith 
9572a42c3cSBarry Smith    Level: advanced
9672a42c3cSBarry Smith 
9795452b02SPatrick Sanan     Notes:
98a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
990f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
100a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
1010f11a792SBarry Smith 
102a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
1031ea5a559SBarry Smith 
10472a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
1050f11a792SBarry Smith */
106945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
10772a42c3cSBarry Smith {
10872a42c3cSBarry Smith   PetscErrorCode ierr;
10972a42c3cSBarry Smith   int            myargc   = argc;
11072a42c3cSBarry Smith   char           **myargs = args;
11172a42c3cSBarry Smith 
11272a42c3cSBarry Smith   PetscFunctionBegin;
113c52ac268SRichard Tran Mills   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr;
1141ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
115df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
11672a42c3cSBarry Smith   PetscFunctionReturn(ierr);
11772a42c3cSBarry Smith }
11872a42c3cSBarry Smith 
119f0865b08SBarry Smith /*
120a56f64adSBarry Smith       Used by Julia interface to get communicator
121f0865b08SBarry Smith */
122945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
123f0865b08SBarry Smith {
124f0865b08SBarry Smith   PetscFunctionBegin;
125f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
126f0865b08SBarry Smith   PetscFunctionReturn(0);
127f0865b08SBarry Smith }
128f0865b08SBarry Smith 
129e5c89e4eSSatish Balay /*@C
130e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
131e5c89e4eSSatish Balay         the command line arguments.
132e5c89e4eSSatish Balay 
133e5c89e4eSSatish Balay    Collective
134e5c89e4eSSatish Balay 
135e5c89e4eSSatish Balay    Level: advanced
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
138e5c89e4eSSatish Balay @*/
1397087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
140e5c89e4eSSatish Balay {
141e5c89e4eSSatish Balay   PetscErrorCode ierr;
142e5c89e4eSSatish Balay   int            argc   = 0;
14302c9f0b5SLisandro Dalcin   char           **args = NULL;
144e5c89e4eSSatish Balay 
145e5c89e4eSSatish Balay   PetscFunctionBegin;
1460298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
147e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
148e5c89e4eSSatish Balay }
149e5c89e4eSSatish Balay 
150e5c89e4eSSatish Balay /*@
151e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
152e5c89e4eSSatish Balay 
15393b6d2d1SJed Brown    Level: beginner
154e5c89e4eSSatish Balay 
155e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
156e5c89e4eSSatish Balay @*/
1577087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool  *isInitialized)
158e5c89e4eSSatish Balay {
159e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
16093b6d2d1SJed Brown   return 0;
161e5c89e4eSSatish Balay }
162e5c89e4eSSatish Balay 
163e5c89e4eSSatish Balay /*@
164e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
165e5c89e4eSSatish Balay 
166e5c89e4eSSatish Balay    Level: developer
167e5c89e4eSSatish Balay 
168e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
169e5c89e4eSSatish Balay @*/
1707087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
171e5c89e4eSSatish Balay {
172e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
17393b6d2d1SJed Brown   return 0;
174e5c89e4eSSatish Balay }
175e5c89e4eSSatish Balay 
17695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void);
177e5c89e4eSSatish Balay 
178e5c89e4eSSatish Balay /*
179e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
180e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
181e5c89e4eSSatish Balay */
182367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
183e5c89e4eSSatish Balay 
184367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
185e5c89e4eSSatish Balay {
186e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
187e5c89e4eSSatish Balay 
188e5c89e4eSSatish Balay   PetscFunctionBegin;
189e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
190e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
19141e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
192e5c89e4eSSatish Balay   }
193e5c89e4eSSatish Balay 
194e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
195e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
196e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
197e5c89e4eSSatish Balay   }
198812af9f3SBarry Smith   PetscFunctionReturnVoid();
199e5c89e4eSSatish Balay }
200e5c89e4eSSatish Balay 
201e5c89e4eSSatish Balay /*
202e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
203e5c89e4eSSatish Balay sum of the second entry.
204b693b147SBarry Smith 
20576f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
206367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
207b693b147SBarry Smith there would be no place to store the both needed results.
208e5c89e4eSSatish Balay */
20976ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
210e5c89e4eSSatish Balay {
211e5c89e4eSSatish Balay   PetscErrorCode ierr;
212e5c89e4eSSatish Balay 
213e5c89e4eSSatish Balay   PetscFunctionBegin;
214d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
215d6e4c47cSJed Brown   {
216d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
217367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
218d6e4c47cSJed Brown     *max = work.max;
219d6e4c47cSJed Brown     *sum = work.sum;
220d6e4c47cSJed Brown   }
221d6e4c47cSJed Brown #else
222d6e4c47cSJed Brown   {
223d6e4c47cSJed Brown     PetscMPIInt    size,rank;
224d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
225e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
226e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
227785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
228367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2296ac3741eSJed Brown     *max = work[rank].max;
2306ac3741eSJed Brown     *sum = work[rank].sum;
231e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
232d6e4c47cSJed Brown   }
233d6e4c47cSJed Brown #endif
234e5c89e4eSSatish Balay   PetscFunctionReturn(0);
235e5c89e4eSSatish Balay }
236e5c89e4eSSatish Balay 
237e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
238e5c89e4eSSatish Balay 
239570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
24006a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
241e5c89e4eSSatish Balay 
2428cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
243e5c89e4eSSatish Balay {
244e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
245e5c89e4eSSatish Balay 
246e5c89e4eSSatish Balay   PetscFunctionBegin;
2477c2de775SJed Brown   if (*datatype == MPIU_REAL) {
248e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
249a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2507c2de775SJed Brown   }
2517c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2527c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2537c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
254a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2557c2de775SJed Brown   }
2567c2de775SJed Brown #endif
2577c2de775SJed Brown   else {
2587c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
25941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
260e2e03761SBarry Smith   }
261812af9f3SBarry Smith   PetscFunctionReturnVoid();
262e5c89e4eSSatish Balay }
263e5c89e4eSSatish Balay #endif
264e5c89e4eSSatish Balay 
265570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
266d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
267d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
268d9822059SBarry Smith 
2698cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
270d9822059SBarry Smith {
271d9822059SBarry Smith   PetscInt i,count = *cnt;
272d9822059SBarry Smith 
273d9822059SBarry Smith   PetscFunctionBegin;
2747c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2758c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
276a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2777c2de775SJed Brown   }
2787c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2797c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2807c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2817c2de775SJed Brown     for (i=0; i<count; i++) {
2827c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2837c2de775SJed Brown     }
2847c2de775SJed Brown   }
2857c2de775SJed Brown #endif
2867c2de775SJed Brown   else {
2877c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
28841e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2898c764dc5SJose Roman   }
290d9822059SBarry Smith   PetscFunctionReturnVoid();
291d9822059SBarry Smith }
292d9822059SBarry Smith 
2938cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
294d9822059SBarry Smith {
295d9822059SBarry Smith   PetscInt    i,count = *cnt;
296d9822059SBarry Smith 
297d9822059SBarry Smith   PetscFunctionBegin;
2987c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2998c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
300a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
3017c2de775SJed Brown   }
3027c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
3037c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3047c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
3057c2de775SJed Brown     for (i=0; i<count; i++) {
3067c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3077c2de775SJed Brown     }
3087c2de775SJed Brown   }
3097c2de775SJed Brown #endif
3107c2de775SJed Brown   else {
3118c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
31241e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
3138c764dc5SJose Roman   }
314d9822059SBarry Smith   PetscFunctionReturnVoid();
315d9822059SBarry Smith }
316d9822059SBarry Smith #endif
317d9822059SBarry Smith 
318480cf27aSJed Brown /*
319480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
320480cf27aSJed Brown 
321ff0e51ddSBarry 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.
322480cf27aSJed Brown 
32312801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
324480cf27aSJed Brown 
325480cf27aSJed Brown */
3268cc058d9SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
327480cf27aSJed Brown {
328480cf27aSJed Brown   PetscErrorCode ierr;
329480cf27aSJed Brown 
330480cf27aSJed Brown   PetscFunctionBegin;
33102c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
33212801b39SBarry Smith   ierr = PetscFree(count_val);CHKERRMPI(ierr);
333480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
334480cf27aSJed Brown }
335480cf27aSJed Brown 
336480cf27aSJed Brown /*
33747435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
338da3039f7SJed Brown   calls MPI_Comm_free().
339da3039f7SJed Brown 
340da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
341480cf27aSJed Brown 
342ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
343480cf27aSJed Brown 
34412801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
345480cf27aSJed Brown 
346480cf27aSJed Brown */
347da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
348480cf27aSJed Brown {
349480cf27aSJed Brown   PetscErrorCode                    ierr;
350b89831e5SBarry Smith   PetscMPIInt                       flg;
351265f3f35SJed Brown   union {MPI_Comm comm; void *ptr;} icomm,ocomm;
352480cf27aSJed Brown 
353480cf27aSJed Brown   PetscFunctionBegin;
35412801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
355ec4fadc2SJed Brown   icomm.ptr = attr_val;
356da3039f7SJed Brown 
35747435625SJed Brown   ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
35812801b39SBarry Smith   if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm");
35912801b39SBarry Smith   if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm");
36047435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
36102c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"User MPI_Comm %ld is being freed after removing reference from inner PETSc comm to this outer comm\n",(long)comm);CHKERRMPI(ierr);
362da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
363b89831e5SBarry Smith }
364da3039f7SJed Brown 
365da3039f7SJed Brown /*
36647435625SJed 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.
367da3039f7SJed Brown  */
368da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
369da3039f7SJed Brown {
370da3039f7SJed Brown   PetscErrorCode ierr;
371da3039f7SJed Brown 
372da3039f7SJed Brown   PetscFunctionBegin;
37302c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
374480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
375480cf27aSJed Brown }
376480cf27aSJed Brown 
3775f7487a0SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Shm(MPI_Comm,PetscMPIInt,void *,void *);
37842218b76SBarry Smith 
379951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3808cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3818cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3828cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
383e39fd77fSBarry Smith #endif
384e39fd77fSBarry Smith 
3850f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
38612801b39SBarry Smith 
387eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
388eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3896ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
39002c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
391dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
392e5c89e4eSSatish Balay 
393dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
394051e4cf2SJed Brown {
395051e4cf2SJed Brown   PetscErrorCode ierr;
396051e4cf2SJed Brown 
397051e4cf2SJed Brown   PetscFunctionBegin;
398051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
399a1601671SBarry 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);
400051e4cf2SJed 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);
401051e4cf2SJed Brown   PetscFunctionReturn(0);
402051e4cf2SJed Brown }
403e5c89e4eSSatish Balay 
4042d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4052d747510SLisandro Dalcin 
4062d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4072d747510SLisandro Dalcin {
4082d747510SLisandro Dalcin   PetscErrorCode ierr;
4092d747510SLisandro Dalcin 
4102d747510SLisandro Dalcin   PetscFunctionBegin;
4112d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
4122d747510SLisandro Dalcin   PetscFunctionReturn(0);
4132d747510SLisandro Dalcin }
4142d747510SLisandro Dalcin 
4152d747510SLisandro Dalcin /*@C
4162d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4172d747510SLisandro Dalcin 
4182d747510SLisandro Dalcin     Not Collective
4192d747510SLisandro Dalcin 
4202d747510SLisandro Dalcin     Input Parameter:
4212d747510SLisandro Dalcin .   len - length of the string name
4222d747510SLisandro Dalcin 
4232d747510SLisandro Dalcin     Output Parameter:
4242d747510SLisandro Dalcin .   name - the name of the running program
4252d747510SLisandro Dalcin 
4262d747510SLisandro Dalcin    Level: advanced
4272d747510SLisandro Dalcin 
4282d747510SLisandro Dalcin     Notes:
4292d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4302d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4312d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4322d747510SLisandro Dalcin @*/
4332d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4342d747510SLisandro Dalcin {
4352d747510SLisandro Dalcin   PetscErrorCode ierr;
4362d747510SLisandro Dalcin 
4372d747510SLisandro Dalcin   PetscFunctionBegin;
4382d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4392d747510SLisandro Dalcin   PetscFunctionReturn(0);
4402d747510SLisandro Dalcin }
4412d747510SLisandro Dalcin 
442e5c89e4eSSatish Balay /*@C
443e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
444e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
445e5c89e4eSSatish Balay 
446e5c89e4eSSatish Balay    Not Collective
447e5c89e4eSSatish Balay 
448e5c89e4eSSatish Balay    Output Parameters:
449e5c89e4eSSatish Balay +  argc - count of number of command line arguments
450e5c89e4eSSatish Balay -  args - the command line arguments
451e5c89e4eSSatish Balay 
452e5c89e4eSSatish Balay    Level: intermediate
453e5c89e4eSSatish Balay 
454e5c89e4eSSatish Balay    Notes:
455e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
456e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
457e5c89e4eSSatish Balay 
458f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
459f177e3b1SBarry Smith 
460793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
461e5c89e4eSSatish Balay 
462e5c89e4eSSatish Balay @*/
4637087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
464e5c89e4eSSatish Balay {
465e5c89e4eSSatish Balay   PetscFunctionBegin;
46617186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
467e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
468e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
469e5c89e4eSSatish Balay   PetscFunctionReturn(0);
470e5c89e4eSSatish Balay }
471e5c89e4eSSatish Balay 
472793721a6SBarry Smith /*@C
473793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
474793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
475793721a6SBarry Smith 
476793721a6SBarry Smith    Not Collective
477793721a6SBarry Smith 
478793721a6SBarry Smith    Output Parameters:
479793721a6SBarry Smith .  args - the command line arguments
480793721a6SBarry Smith 
481793721a6SBarry Smith    Level: intermediate
482793721a6SBarry Smith 
483793721a6SBarry Smith    Notes:
484793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
485793721a6SBarry Smith 
486793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
487793721a6SBarry Smith 
488793721a6SBarry Smith @*/
4897087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
490793721a6SBarry Smith {
491793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
492793721a6SBarry Smith   PetscErrorCode ierr;
493793721a6SBarry Smith 
494793721a6SBarry Smith   PetscFunctionBegin;
49517186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4962d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
497785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
498793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
499793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
500793721a6SBarry Smith   }
5012d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
502793721a6SBarry Smith   PetscFunctionReturn(0);
503793721a6SBarry Smith }
504793721a6SBarry Smith 
505793721a6SBarry Smith /*@C
506793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
507793721a6SBarry Smith 
508793721a6SBarry Smith    Not Collective
509793721a6SBarry Smith 
510793721a6SBarry Smith    Output Parameters:
511793721a6SBarry Smith .  args - the command line arguments
512793721a6SBarry Smith 
513793721a6SBarry Smith    Level: intermediate
514793721a6SBarry Smith 
515793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
516793721a6SBarry Smith 
517793721a6SBarry Smith @*/
5187087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
519793721a6SBarry Smith {
520793721a6SBarry Smith   PetscInt       i = 0;
521793721a6SBarry Smith   PetscErrorCode ierr;
522793721a6SBarry Smith 
523793721a6SBarry Smith   PetscFunctionBegin;
524a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
525793721a6SBarry Smith   while (args[i]) {
526793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
527793721a6SBarry Smith     i++;
528793721a6SBarry Smith   }
529793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
530793721a6SBarry Smith   PetscFunctionReturn(0);
531793721a6SBarry Smith }
532793721a6SBarry Smith 
53311525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
53430befbd2SBarry Smith #include <petscconfiginfo.h>
53530befbd2SBarry Smith 
53695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
53711525c0dSBarry Smith {
53811525c0dSBarry Smith   if (!PetscGlobalRank) {
53930befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
54011525c0dSBarry Smith     int            port;
541ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
54211525c0dSBarry Smith     size_t         applinelen,introlen;
54311525c0dSBarry Smith     PetscErrorCode ierr;
544ffbd1cfbSBarry Smith     char           sawsurl[256];
54511525c0dSBarry Smith 
546c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
54711525c0dSBarry Smith     if (flg) {
54811525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
54911525c0dSBarry Smith 
550c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
55111525c0dSBarry Smith       if (sawslog[0]) {
55211525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
55311525c0dSBarry Smith       } else {
55411525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
55511525c0dSBarry Smith       }
55611525c0dSBarry Smith     }
557c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
55811525c0dSBarry Smith     if (flg) {
55911525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
56011525c0dSBarry Smith     }
561c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
562ffbd1cfbSBarry Smith     if (selectport) {
563ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
564ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
565ffbd1cfbSBarry Smith     } else {
566c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
56711525c0dSBarry Smith       if (flg) {
56811525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
56911525c0dSBarry Smith       }
570ffbd1cfbSBarry Smith     }
571c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
57211525c0dSBarry Smith     if (flg) {
57311525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
57411525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5759c1e0ce8SBarry Smith     } else {
576c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5779c1e0ce8SBarry Smith       if (flg) {
5783c01dfcfSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
5799c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5809c1e0ce8SBarry Smith       }
58111525c0dSBarry Smith     }
582c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
58311525c0dSBarry Smith     if (flg2) {
58411525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
58511525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
58611525c0dSBarry Smith       ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr);
58711525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
58811525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
58943da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
59011525c0dSBarry Smith     }
59111525c0dSBarry Smith     ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr);
59211525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
59311525c0dSBarry Smith     introlen   = 4096 + applinelen;
59430a8c9c0SSurtai Han     applinelen += 1024;
59511525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
59611525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
59711525c0dSBarry Smith 
59811525c0dSBarry Smith     if (rootlocal) {
59911525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
60011525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
60111525c0dSBarry Smith     }
60276a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
60311525c0dSBarry Smith     if (rootlocal && help) {
604928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running <a href=\"%s.c.html\">%s</a> %s</center><br><center><pre>%s</pre></center><br>\n",programname,programname,options,help);CHKERRQ(ierr);
60511525c0dSBarry Smith     } else if (help) {
606928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
60711525c0dSBarry Smith     } else {
608928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
60911525c0dSBarry Smith     }
610b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
61130befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
61211525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
613a8d69d7bSBarry 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"
614df62c222SBarry Smith                                     "<center>This is the default PETSc application dashboard, from it you can access any published PETSc objects or logging data</center><br><center>%s configured with %s</center><br>\n"
615928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
61643da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
61711525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
61811525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
619ffbd1cfbSBarry Smith     if (selectport) {
620aa573868SBarry Smith       PetscBool silent;
6217d812c46SBarry Smith 
6227d812c46SBarry Smith       ierr = SAWs_Initialize();
6237d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6247d812c46SBarry Smith       while (ierr) {
6257d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6267d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6277d812c46SBarry Smith         ierr = SAWs_Initialize();
6287d812c46SBarry Smith       }
6297d812c46SBarry Smith 
630aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
631aa573868SBarry Smith       if (!silent) {
632ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
633ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
634ffbd1cfbSBarry Smith       }
6357d812c46SBarry Smith     } else {
6367d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
637aa573868SBarry Smith     }
6380af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6390af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6400af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6410af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6420af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
64311525c0dSBarry Smith   }
64411525c0dSBarry Smith   PetscFunctionReturn(0);
64511525c0dSBarry Smith }
64611525c0dSBarry Smith #endif
64711525c0dSBarry Smith 
6484dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6494dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6504dfee713SSatish Balay {
6514dfee713SSatish Balay   PetscFunctionBegin;
6524dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6534dfee713SSatish Balay     /* see MPI.py for details on this bug */
6544dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6554dfee713SSatish Balay #endif
6564dfee713SSatish Balay   PetscFunctionReturn(0);
6574dfee713SSatish Balay }
6584dfee713SSatish Balay 
659a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
660a56f64adSBarry Smith #include <adios.h>
66122580e64SBarry Smith #include <adios_read.h>
6627b56e58cSSatish Balay int64_t Petsc_adios_group;
663a56f64adSBarry Smith #endif
66455d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
66555d657eeSBarry Smith #include <adios2_c.h>
66655d657eeSBarry Smith #endif
667cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
668cd1b99a6SBarry Smith #include <omp.h>
669f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
670cd1b99a6SBarry Smith #endif
671a56f64adSBarry Smith 
672bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H)
673bc8a24ecSBarry Smith #include <dlfcn.h>
674bc8a24ecSBarry Smith #endif
675bc8a24ecSBarry Smith 
676e5c89e4eSSatish Balay /*@C
677e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
678e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
679e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
680e5c89e4eSSatish Balay    your program -- usually the very first line!
681e5c89e4eSSatish Balay 
682e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
683e5c89e4eSSatish Balay 
684e5c89e4eSSatish Balay    Input Parameters:
685e5c89e4eSSatish Balay +  argc - count of number of command line arguments
686e5c89e4eSSatish Balay .  args - the command line arguments
6870298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
688fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6890298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
690e5c89e4eSSatish Balay 
69105827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
69205827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
69305827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
69405827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
69505827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
696e5c89e4eSSatish Balay 
697e5c89e4eSSatish Balay    Options Database Keys:
6987ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6997ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
700e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
7018a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
7028a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
7038a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
7048a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
7058a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
706e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
707e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
708e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
709e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
71079dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
71179dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
71279dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
713aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
71492f119d6SBarry Smith .  -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds, ignored in optimized build. May want to set in PETSC_OPTIONS environmental variable
71592f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
71692f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
71792f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
718e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
719e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
720e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
721e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
722e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7230841954dSBarry Smith -  -memory_view - Print memory usage at end of run
724e5c89e4eSSatish Balay 
725e5c89e4eSSatish Balay    Options Database Keys for Profiling:
726a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
727*e94e781bSJacob Faibussowitsch +  -info [optional filename][:[~]optional list,of,classnames][:[~]optional "self"] - Prints verbose information to the screen. See PetscInfo().
728217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
729217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
730495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
731e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7329a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
73379dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
7349a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
735495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
736f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
737495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
738495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
739c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
74087c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
741c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
742495fc317SBarry Smith 
743609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
744e5c89e4eSSatish Balay 
745ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
746ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
747ffbd1cfbSBarry 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
748ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
749ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
750ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
751ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
752ffbd1cfbSBarry Smith 
753e5c89e4eSSatish Balay    Environmental Variables:
754e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
755e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
756e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
757e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
758e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
759e5c89e4eSSatish Balay 
760e5c89e4eSSatish Balay 
761e5c89e4eSSatish Balay    Level: beginner
762e5c89e4eSSatish Balay 
763e5c89e4eSSatish Balay    Notes:
764e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
765e5c89e4eSSatish Balay    it before PetscInitialize().
766e5c89e4eSSatish Balay 
767e5c89e4eSSatish Balay    Fortran Version:
768e5c89e4eSSatish Balay    In Fortran this routine has the format
769e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
770e5c89e4eSSatish Balay 
771e5c89e4eSSatish Balay +   ierr - error return code
7720eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
773fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
774e5c89e4eSSatish Balay 
775e5c89e4eSSatish Balay    Important Fortran Note:
7760eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7770298fd71SBarry Smith    null character string; you CANNOT just use NULL as
778a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
779e5c89e4eSSatish Balay 
78001cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
78101cb0274SBarry Smith    calling PetscInitialize().
782e5c89e4eSSatish Balay 
78301cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
784e5c89e4eSSatish Balay 
785e5c89e4eSSatish Balay @*/
7867087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
787e5c89e4eSSatish Balay {
788e5c89e4eSSatish Balay   PetscErrorCode ierr;
7894bb5149bSJed Brown   PetscMPIInt    flag, size;
79019c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
791e5c89e4eSSatish Balay   char           hostname[256];
792e5c89e4eSSatish Balay 
793e5c89e4eSSatish Balay   PetscFunctionBegin;
794e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7953d96e996SBarry Smith   /*
7963d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7973d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7983d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7993d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
8003d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
8013d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
8023d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
8033d96e996SBarry Smith       listed above and since that time are compatible.
8043d96e996SBarry Smith 
8053d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
8063d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
8073d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
8083d96e996SBarry Smith       to perform the checking.
8093d96e996SBarry Smith 
8103d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
8113d96e996SBarry Smith 
8123d96e996SBarry Smith       Questions:
8133d96e996SBarry Smith 
8143d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
8153d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
8163d96e996SBarry Smith   */
8173d96e996SBarry Smith 
81819c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
81919c5658aSBarry Smith   {
82019c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
82119c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
82219c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8233d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
82419c5658aSBarry Smith #if defined(MPICH_VERSION)
8253d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
82619c5658aSBarry Smith     {
82719c5658aSBarry Smith       char *ver,*lf;
82819c5658aSBarry Smith       flg = PETSC_FALSE;
82919c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
83019c5658aSBarry Smith       if (ver) {
83119c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
83219c5658aSBarry Smith         if (lf) {
83319c5658aSBarry Smith           *lf = 0;
83419c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
83519c5658aSBarry Smith         }
83619c5658aSBarry Smith       }
83719c5658aSBarry Smith       if (!flg) {
83819c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8393d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
84019c5658aSBarry Smith       }
84119c5658aSBarry Smith     }
8423d96e996SBarry Smith #endif
8433d96e996SBarry 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?) */
84419c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
84519c5658aSBarry Smith     {
84619c5658aSBarry Smith       char *ver,bs[32],*bsf;
84719c5658aSBarry Smith       flg = PETSC_FALSE;
84819c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
84919c5658aSBarry Smith       if (ver) {
8502e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
85119c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
85219c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
85319c5658aSBarry Smith       }
85419c5658aSBarry Smith       if (!flg) {
85519c5658aSBarry 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);
8563d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
85719c5658aSBarry Smith       }
85819c5658aSBarry Smith     }
85919c5658aSBarry Smith #endif
86019c5658aSBarry Smith   }
86119c5658aSBarry Smith #endif
86219c5658aSBarry Smith 
863bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM)
864bc8a24ecSBarry Smith   {
865bc8a24ecSBarry Smith     PetscInt cnt = 0;
866bc8a24ecSBarry Smith     /* These symbols are currently in the OpenMPI and MPICH libraries; they may not always be, in that case the test will simply not detect the problem */
867bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++;
868bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++;
869bc8a24ecSBarry Smith     if (cnt > 1) {
870bc8a24ecSBarry Smith       fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
871bc8a24ecSBarry Smith       return PETSC_ERR_MPI_LIB_INCOMP;
872bc8a24ecSBarry Smith     }
873bc8a24ecSBarry Smith   }
874bc8a24ecSBarry Smith #endif
875bc8a24ecSBarry Smith 
876ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
877d89683f4Sbcordonn   PETSC_STDOUT = stdout;
878ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
879e5c89e4eSSatish Balay 
8800c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8810c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8820c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8830c30907bSSatish Balay #endif
8840c30907bSSatish Balay 
8854416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
886e5c89e4eSSatish Balay 
887e5c89e4eSSatish Balay   /*
888e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
889e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
890e5c89e4eSSatish Balay   */
891e5c89e4eSSatish Balay   if (argc && *argc) {
892e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
893e5c89e4eSSatish Balay   } else {
894e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
895e5c89e4eSSatish Balay   }
896e5c89e4eSSatish Balay 
897e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
898e5c89e4eSSatish Balay   if (!flag) {
899e32f2f54SBarry 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");
9004dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
9015e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
9025e765c61SJed Brown     {
9035e765c61SJed Brown       PetscMPIInt provided;
9045e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
9055e765c61SJed Brown     }
9065e765c61SJed Brown #else
907e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
9085e765c61SJed Brown #endif
909e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
910e5c89e4eSSatish Balay   }
9114dfee713SSatish Balay 
912e5c89e4eSSatish Balay   if (argc && args) {
913e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
914e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
915e5c89e4eSSatish Balay   }
916e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
9175ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
9185ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
9195ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
920ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
921e5c89e4eSSatish Balay 
922a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
923d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
924e5c89e4eSSatish Balay 
9250f9be574SPeter Hill   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
92612801b39SBarry Smith     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
92712801b39SBarry Smith     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
9280f9be574SPeter Hill   }
92912801b39SBarry Smith 
930e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
931e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
932e5c89e4eSSatish Balay 
933e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
934e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
935e5c89e4eSSatish Balay 
9368ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9378ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9387cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
939e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
940e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
941e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
942e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
943e316c87fSJed Brown #endif
944e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9458ad47952SJed Brown 
946e5c89e4eSSatish Balay   /*
947e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
948e5c89e4eSSatish Balay      shared libraries the constructors for global variables
949e5c89e4eSSatish Balay      are not called; at least on IRIX.
950e5c89e4eSSatish Balay   */
951886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
952e5c89e4eSSatish Balay   {
953252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
95450f81f78SJed Brown     PetscComplex ic(0.0,1.0);
955e5c89e4eSSatish Balay     PETSC_i = ic;
956252ecd31SSatish Balay #else
95750f81f78SJed Brown     PETSC_i = _Complex_I;
958b7940d39SSatish Balay #endif
959762437b8SSatish Balay   }
960762437b8SSatish Balay 
9612c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
962e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
963500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
964500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
965500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9662c876bd9SBarry Smith #endif
967886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
968e5c89e4eSSatish Balay 
969e5c89e4eSSatish Balay   /*
970e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
971e5c89e4eSSatish Balay      half of the entries and maxes the second half.
972e5c89e4eSSatish Balay   */
973367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
974e5c89e4eSSatish Balay 
975ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
976c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
977c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9787c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9798c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9808c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9818c764dc5SJose Roman #endif
982d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
983d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
984570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
985570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
986570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
987570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
988570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
989c90a1750SBarry Smith #endif
990c90a1750SBarry Smith 
991570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
992cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
993cca4cb22SSatish Balay #endif
994cca4cb22SSatish Balay 
995e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
996e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
997e5c89e4eSSatish Balay 
99840df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
999e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
1000e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
100144041f26SJed Brown #endif
1002e5c89e4eSSatish Balay 
1003ec957eceSBarry Smith 
1004e5c89e4eSSatish Balay   /*
1005480cf27aSJed Brown      Attributes to be set on PETSc communicators
1006480cf27aSJed Brown   */
100712801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
100812801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
100912801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
10105f7487a0SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Shm,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
1011480cf27aSJed Brown 
1012480cf27aSJed Brown   /*
1013e8fb0fc0SBarry Smith      Build the options database
1014e5c89e4eSSatish Balay   */
1015c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
1016e5c89e4eSSatish Balay 
1017703f3eceSBarry Smith   /* call a second time so it can look in the options database */
1018703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
10196dc8fec2Sbcordonn 
1020e5c89e4eSSatish Balay   /*
1021e5c89e4eSSatish Balay      Print main application help message
1022e5c89e4eSSatish Balay   */
10232d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
1024e5c89e4eSSatish Balay   if (help && flg) {
1025e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
1026e5c89e4eSSatish Balay   }
1027e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
1028e5c89e4eSSatish Balay 
1029d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
1030d45a07a7SBarry Smith 
1031e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
103211525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1033f4202a44SBarry Smith #endif
1034f4202a44SBarry Smith 
1035e5c89e4eSSatish Balay   /*
1036e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1037e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1038e5c89e4eSSatish Balay   */
1039e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1040e5c89e4eSSatish Balay 
1041e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
104202c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1043e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
104402c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1045cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1046101491d0SBarry Smith   {
1047101491d0SBarry Smith     PetscBool omp_view_flag;
1048cd1b99a6SBarry Smith     char      *threads = getenv("OMP_NUM_THREADS");
1049e5c89e4eSSatish Balay 
1050cd1b99a6SBarry Smith    if (threads) {
105102c9f0b5SLisandro Dalcin      ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
1052f90b075cSBarry Smith      (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
1053cd1b99a6SBarry Smith    } else {
1054cd1b99a6SBarry Smith #define NMAX  10000
1055101491d0SBarry Smith      int          i;
1056cd1b99a6SBarry Smith       PetscScalar *x;
1057cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1058cd1b99a6SBarry Smith #pragma omp parallel for
1059cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1060cd1b99a6SBarry Smith         x[i] = 0.0;
1061f90b075cSBarry Smith         PetscNumOMPThreads  = (PetscInt) omp_get_num_threads();
1062cd1b99a6SBarry Smith       }
1063cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
106402c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1065101491d0SBarry Smith     }
1066101491d0SBarry Smith     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
1067f90b075cSBarry Smith     ierr = PetscOptionsInt("-omp_num_threads","Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS","None",PetscNumOMPThreads,&PetscNumOMPThreads,&flg);CHKERRQ(ierr);
1068101491d0SBarry Smith     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
1069101491d0SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1070101491d0SBarry Smith     if (flg) {
107102c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1072f90b075cSBarry Smith       omp_set_num_threads((int)PetscNumOMPThreads);
1073101491d0SBarry Smith     }
1074101491d0SBarry Smith     if (omp_view_flag) {
1075f90b075cSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
1076cd1b99a6SBarry Smith     }
1077cd1b99a6SBarry Smith   }
1078cd1b99a6SBarry Smith #endif
1079e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr);
1080ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1081c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1082ef6c6fedSBoyana Norris 
1083951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1084e39fd77fSBarry Smith   /*
1085e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1086e39fd77fSBarry Smith 
1087e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1088e39fd77fSBarry Smith   */
108930815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
10900298fd71SBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
109130815ce0SLisandro Dalcin   }
1092951e3c8eSBarry Smith #endif
1093e39fd77fSBarry Smith 
109441c0b4b3SShri Abhyankar   /*
109541c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
109641c0b4b3SShri Abhyankar   */
1097ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1098e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1099e1167bb9SShri Abhyankar #endif
1100e1167bb9SShri Abhyankar 
11012d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
11022d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
11032d53ad75SBarry Smith #endif
11042d53ad75SBarry Smith 
11055e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
110687c3beb6SLisandro Dalcin   {
110787c3beb6SLisandro Dalcin     PetscViewer viewer;
110822e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
11095e71baefSBarry Smith     if (flg) {
11105e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
111187c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
111287c3beb6SLisandro Dalcin     }
11135e71baefSBarry Smith   }
11145e71baefSBarry Smith #endif
1115dff31646SBarry Smith 
111687c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
111787c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
111887c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
111987c3beb6SLisandro Dalcin 
1120a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1121a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
11227b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1123a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
11249fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1125a56f64adSBarry Smith #endif
112655d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
112755d657eeSBarry Smith #endif
1128a56f64adSBarry Smith 
1129301d30feSBarry Smith   /*
1130301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1131301d30feSBarry Smith   */
1132301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
11332db0e300SLisandro Dalcin 
11342db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
11352db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1136eb02082bSJunchao Zhang 
1137eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1138eb02082bSJunchao Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-use_gpu_aware_mpi",&use_gpu_aware_mpi,NULL);CHKERRQ(ierr);
1139eb02082bSJunchao Zhang #if defined(PETSC_HAVE_OPENMPI_MAJOR) && defined(MPIX_CUDA_AWARE_SUPPORT)
1140eb02082bSJunchao Zhang   /* OpenMPI supports compile time and runtime cuda support checking */
1141eb02082bSJunchao Zhang   if (use_gpu_aware_mpi && 1 != MPIX_Query_cuda_support()) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"OpenMPI you used does not support CUDA while you requested -use_gpu_aware_mpi");
1142eb02082bSJunchao Zhang #endif
1143eb02082bSJunchao Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-sf_use_default_cuda_stream",&sf_use_default_cuda_stream,NULL);CHKERRQ(ierr);
1144eb02082bSJunchao Zhang #endif
1145eb02082bSJunchao Zhang 
1146301d30feSBarry Smith   PetscFunctionReturn(0);
1147e5c89e4eSSatish Balay }
1148e5c89e4eSSatish Balay 
11494097062eSBarry Smith #if defined(PETSC_USE_LOG)
115095c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1151ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1152ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
115395c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11544097062eSBarry Smith #endif
1155e5c89e4eSSatish Balay 
1156008a6e76SBarry Smith /*
1157008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1158008a6e76SBarry Smith */
1159008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1160008a6e76SBarry Smith {
1161008a6e76SBarry Smith   PetscErrorCode ierr;
1162008a6e76SBarry Smith 
1163008a6e76SBarry Smith   PetscFunctionBegin;
1164008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1165008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1166008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1167008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1168008a6e76SBarry Smith #endif
1169008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1170008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1171008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1172008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1173008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1174008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1175008a6e76SBarry Smith #endif
1176008a6e76SBarry Smith 
1177008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1178008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1179008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1180008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1181008a6e76SBarry Smith #endif
1182008a6e76SBarry Smith #endif
1183008a6e76SBarry Smith 
1184008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1185008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1186008a6e76SBarry Smith #endif
1187008a6e76SBarry Smith 
1188008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
118940df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1190008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1191008a6e76SBarry Smith #endif
1192008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1193008a6e76SBarry Smith   PetscFunctionReturn(0);
1194008a6e76SBarry Smith }
1195008a6e76SBarry Smith 
1196e5c89e4eSSatish Balay /*@C
1197e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1198e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1199e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1200e5c89e4eSSatish Balay 
1201e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1202e5c89e4eSSatish Balay 
1203e5c89e4eSSatish Balay    Options Database Keys:
120426a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1205e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
12067eb1d149SBarry 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
1207e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
120892f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1209e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
121092f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1211e5c89e4eSSatish Balay 
1212e5c89e4eSSatish Balay    Level: beginner
1213e5c89e4eSSatish Balay 
1214e5c89e4eSSatish Balay    Note:
1215e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1216e5c89e4eSSatish Balay 
121788c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1218e5c89e4eSSatish Balay @*/
12197087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1220e5c89e4eSSatish Balay {
1221e5c89e4eSSatish Balay   PetscErrorCode ierr;
12224bb5149bSJed Brown   PetscMPIInt    rank;
1223a8d2bbe5SBarry Smith   PetscInt       nopt;
12242bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1225dff31646SBarry Smith   PetscBool      flg;
122610463e74SBarry Smith #if defined(PETSC_USE_LOG)
122710463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
122810463e74SBarry Smith #endif
1229e5c89e4eSSatish Balay 
1230e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
12314b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
12323cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1233e5c89e4eSSatish Balay   }
12343cbbc5ffSBarry Smith   PetscFunctionBegin;
12350298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1236b022a5c1SBarry Smith 
12371f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1238a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
123922580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1240a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1241a56f64adSBarry Smith #endif
124255d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
124355d657eeSBarry Smith #endif
1244c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1245dff31646SBarry Smith   if (flg) {
12461f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
12471f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
12481f817a21SBarry Smith 
1249c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
12501f817a21SBarry Smith     if (filename[0]) {
12511f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
12521f817a21SBarry Smith     }
1253dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1254dff31646SBarry Smith     cits[0] = 0;
1255dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
12561f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
12571f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12581f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
12591f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12601f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1261dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1262dff31646SBarry Smith   }
1263dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1264dff31646SBarry Smith 
1265c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
126604102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
126704102261SBarry Smith   {
126804102261SBarry Smith     PetscInt nmax = 2;
126904102261SBarry Smith     char     **buffs;
127004102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1271c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
127204102261SBarry Smith     if (flg1) {
127304102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
127404102261SBarry Smith       if (nmax == 1) {
127504102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
127604102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
127704102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
127804102261SBarry Smith       }
127904102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
128004102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
128104102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
128204102261SBarry Smith     }
128304102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
128404102261SBarry Smith   }
1285763ec7b1SBarry Smith   {
1286763ec7b1SBarry Smith     PetscInt nmax = 2;
1287763ec7b1SBarry Smith     char     **buffs;
1288763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1289763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1290763ec7b1SBarry Smith     if (flg1) {
1291763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1292763ec7b1SBarry Smith       if (nmax == 1) {
1293763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1294763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1295763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1296763ec7b1SBarry Smith       }
1297763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1298763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1299763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1300763ec7b1SBarry Smith     }
1301763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1302763ec7b1SBarry Smith   }
130304102261SBarry Smith #endif
130467234432SDmitry Karpeev   /*
130567234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
130667234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
130767234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
130867234432SDmitry Karpeev   */
130967234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
131004102261SBarry Smith 
13112d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
13122d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
13132d53ad75SBarry Smith #endif
13142d53ad75SBarry Smith 
13152d53ad75SBarry Smith 
1316e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1317dff31646SBarry Smith   flg = PETSC_FALSE;
1318c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1319d5649816SBarry Smith   if (flg) {
1320e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1321d5649816SBarry Smith   }
1322d5649816SBarry Smith #endif
1323d5649816SBarry Smith 
1324681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1325681455b2SBarry Smith   flg1 = PETSC_FALSE;
1326c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1327681455b2SBarry Smith   if (flg1) {
1328681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1329681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1330681455b2SBarry Smith   }
1331681455b2SBarry Smith #endif
1332681455b2SBarry Smith 
133367584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1334c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1335e5c89e4eSSatish Balay   if (!flg2) {
133690d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1337c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1338e5c89e4eSSatish Balay   }
1339e5c89e4eSSatish Balay   if (flg2) {
13400841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1341e5c89e4eSSatish Balay   }
134267584ceeSBarry Smith #endif
1343e5c89e4eSSatish Balay 
1344e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
134590d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1346c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1347e5c89e4eSSatish Balay   if (flg1) {
1348e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1349205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1350e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1351e5c89e4eSSatish Balay   }
1352e5c89e4eSSatish Balay #endif
1353e5c89e4eSSatish Balay 
1354e5c89e4eSSatish Balay 
1355e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1356e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1357e5c89e4eSSatish Balay   mname[0] = 0;
1358c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1359e5c89e4eSSatish Balay   if (flg1) {
1360e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1361e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1362e5c89e4eSSatish Balay   }
1363e5c89e4eSSatish Balay #endif
1364356e5837SBarry Smith #endif
1365a297a907SKarl Rupp 
1366dd710f27SBarry Smith   /*
1367dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1368dd710f27SBarry Smith   */
1369dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1370dd710f27SBarry Smith 
1371356e5837SBarry Smith #if defined(PETSC_USE_LOG)
137287c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1373f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
137487c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
137587c3beb6SLisandro Dalcin 
1376356e5837SBarry Smith   mname[0] = 0;
1377c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1378e5c89e4eSSatish Balay   if (flg1) {
137991eabc43SBarry Smith     PetscViewer viewer;
138020a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
138191eabc43SBarry Smith     if (mname[0]) {
138291eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
138391eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13846bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
138533f85c2fSBarry Smith     } else {
138633f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13879a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
138833f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13899a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
139033f85c2fSBarry Smith     }
1391e5c89e4eSSatish Balay   }
1392a297a907SKarl Rupp 
1393dd710f27SBarry Smith   /*
1394dd710f27SBarry Smith      Free any objects created by the last block of code.
1395dd710f27SBarry Smith   */
1396dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1397dd710f27SBarry Smith 
1398dd710f27SBarry Smith   mname[0] = 0;
1399c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1400c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
14017ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1402e5c89e4eSSatish Balay #endif
140310463e74SBarry Smith 
140433f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
140510463e74SBarry Smith 
140690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1407c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1408e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
140990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1410c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1411e5c89e4eSSatish Balay   if (flg1) {
1412e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1413e5c89e4eSSatish Balay   }
141490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
141590d69ab7SBarry Smith   flg2 = PETSC_FALSE;
14168bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1417c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1418c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1419c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1420e4c476e2SSatish Balay 
1421e5c89e4eSSatish Balay   if (flg2) {
1422be56827dSJed Brown     PetscViewer viewer;
142302ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
142402ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1425c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1426be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1427e5c89e4eSSatish Balay   }
1428e5c89e4eSSatish Balay 
1429e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1430c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1431c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1432e5c89e4eSSatish Balay 
143333fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1434c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
14353de2bfdfSBarry Smith #if defined(PETSC_USE_DEBUG)
14363de2bfdfSBarry Smith   if (!flg1) flg3 = PETSC_TRUE;
14373de2bfdfSBarry Smith #endif
1438e5c89e4eSSatish Balay   if (flg3) {
14393de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1440be56827dSJed Brown       PetscViewer viewer;
144102ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
144202ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1443c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1444be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1445e5c89e4eSSatish Balay     }
14463de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
14473de2bfdfSBarry Smith     if (nopt) {
14483de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
14493de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
14503de2bfdfSBarry Smith       if (nopt == 1) {
1451e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1452e5c89e4eSSatish Balay       } else {
14537582186dSLisandro Dalcin         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1454e5c89e4eSSatish Balay       }
14553de2bfdfSBarry Smith     } else if (flg3 && flg1) {
14563de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1457df12ba86SBarry Smith     }
1458c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1459e5c89e4eSSatish Balay   }
1460e5c89e4eSSatish Balay 
1461e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1462d45a07a7SBarry Smith   if (!PetscGlobalRank) {
146387f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
146416ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1465d45a07a7SBarry Smith   }
1466ec957eceSBarry Smith #endif
1467ec957eceSBarry Smith 
14684097062eSBarry Smith #if defined(PETSC_USE_LOG)
146910463e74SBarry Smith   /*
1470dbc8283eSBarry Smith        List all objects the user may have forgot to free
14712eff7a51SBarry Smith   */
147205df10baSBarry Smith   if (PetscObjectsLog) {
1473c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1474a64a8e02SBarry Smith     if (flg1) {
1475a64a8e02SBarry Smith       MPI_Comm local_comm;
14767eb1d149SBarry Smith       char     string[64];
1477a64a8e02SBarry Smith 
1478c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1479a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1480a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14817eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1482a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1483a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14840a1571b3SBarry Smith     }
148505df10baSBarry Smith   }
14864097062eSBarry Smith #endif
14874097062eSBarry Smith 
14884097062eSBarry Smith #if defined(PETSC_USE_LOG)
1489dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1490dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1491a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14924097062eSBarry Smith #endif
14932eff7a51SBarry Smith 
149433f85c2fSBarry Smith   /*
149533f85c2fSBarry Smith      Destroy any packages that registered a finalize
149633f85c2fSBarry Smith   */
149733f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
149833f85c2fSBarry Smith 
1499101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1500fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1501101409b8SToby Isaac #endif
1502101409b8SToby Isaac 
150333f85c2fSBarry Smith   /*
150448dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
150548dd1dffSBarry Smith 
150637e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
150748dd1dffSBarry Smith   */
150837e93019SBarry Smith 
15094028d114SSatish Balay   if (petsc_history) {
1510f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
151102c9f0b5SLisandro Dalcin     petsc_history = NULL;
1512e5c89e4eSSatish Balay   }
15139de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1514*e94e781bSJacob Faibussowitsch   ierr = PetscInfoDestroy();CHKERRQ(ierr);
1515e5c89e4eSSatish Balay 
151667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
151792f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1518e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
151992f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1520e5c89e4eSSatish Balay     FILE *fd;
1521ed9cf6e9SBarry Smith     int  err;
1522e5c89e4eSSatish Balay 
1523dc92acbaSJed Brown     flg2 = PETSC_FALSE;
152492f119d6SBarry Smith     flg3 = PETSC_FALSE;
15258bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
152692f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);
1527dc92acbaSJed Brown #endif
152892f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
152992f119d6SBarry Smith     fname[0] = 0;
153092f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1531e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1532e5c89e4eSSatish Balay 
15332e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1534e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1535e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1536ed9cf6e9SBarry Smith       err  = fclose(fd);
1537e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
153892f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1539e5c89e4eSSatish Balay       MPI_Comm local_comm;
1540e5c89e4eSSatish Balay 
1541e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1542e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1543e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1544e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1545e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1546e5c89e4eSSatish Balay     }
1547e5c89e4eSSatish Balay     fname[0] = 0;
154892f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,250,&flg1);CHKERRQ(ierr);
1549e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1550e5c89e4eSSatish Balay 
155192f119d6SBarry Smith       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
155292f119d6SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
155392f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1554ed9cf6e9SBarry Smith       err  = fclose(fd);
1555e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
155692f119d6SBarry Smith     } else if (flg1) {
155792f119d6SBarry Smith       MPI_Comm local_comm;
155892f119d6SBarry Smith 
155992f119d6SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
156092f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
156192f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
156292f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
156392f119d6SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1564e5c89e4eSSatish Balay     }
1565e5c89e4eSSatish Balay   }
156667584ceeSBarry Smith #endif
156720e2c332SMatthew G. Knepley 
15685486ca60SMatthew G. Knepley   /*
15695486ca60SMatthew G. Knepley      Close any open dynamic libraries
15705486ca60SMatthew G. Knepley   */
15715486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15725486ca60SMatthew G. Knepley 
1573e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15744416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1575e5c89e4eSSatish Balay 
1576e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
157702c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1578e5c89e4eSSatish Balay 
1579008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1580e5c89e4eSSatish Balay 
1581dbc8283eSBarry Smith   /*
1582efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1583efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1584efb80d3cSBarry Smith 
1585efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1586efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1587dbc8283eSBarry Smith  */
1588b770b1f6SSatish Balay   {
1589dbc8283eSBarry Smith     PetscCommCounter *counter;
1590dbc8283eSBarry Smith     PetscMPIInt      flg;
1591dbc8283eSBarry Smith     MPI_Comm         icomm;
1592265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
159347435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1594dbc8283eSBarry Smith     if (flg) {
1595265f3f35SJed Brown       icomm = ucomm.comm;
159647435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1597dbc8283eSBarry 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");
1598dbc8283eSBarry Smith 
159947435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
160047435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1601efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1602dbc8283eSBarry Smith     }
160347435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1604dbc8283eSBarry Smith     if (flg) {
1605265f3f35SJed Brown       icomm = ucomm.comm;
160647435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1607dbc8283eSBarry 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");
1608dbc8283eSBarry Smith 
160947435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
161047435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1611efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1612dbc8283eSBarry Smith     }
1613b770b1f6SSatish Balay   }
1614dbc8283eSBarry Smith 
161547435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
161647435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
161747435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
16185f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1619480cf27aSJed Brown 
16205ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
16215ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
16225ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1623ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1624ef19f930SBarry Smith 
1625e5c89e4eSSatish Balay   if (PetscBeganMPI) {
162699608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
162799b1327fSBarry Smith     PetscMPIInt flag;
162899b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1629e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
163099608316SBarry Smith #endif
1631e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1632e5c89e4eSSatish Balay   }
1633e5c89e4eSSatish Balay /*
1634e5c89e4eSSatish Balay 
1635e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1636e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1637e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1638e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1639e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16400e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1641e5c89e4eSSatish Balay    memory was not freed.
1642e5c89e4eSSatish Balay 
1643e5c89e4eSSatish Balay */
16441d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1645a297a907SKarl Rupp 
1646e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1647e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
16483db9a53dSBarry Smith   PetscFunctionReturn(0);
1649e5c89e4eSSatish Balay }
1650e5c89e4eSSatish Balay 
165143db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16528cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
165343db4dbbSBarry Smith {
165443db4dbbSBarry Smith   if (*a == *b) return 1;
165543db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
165643db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
165743db4dbbSBarry Smith   return 0;
165843db4dbbSBarry Smith }
1659a70650f6SBarry Smith #endif
166043db4dbbSBarry Smith 
166143db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16628cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
166343db4dbbSBarry Smith {
166443db4dbbSBarry Smith   if (*a == *b) return 1;
166543db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
166643db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
166743db4dbbSBarry Smith   return 0;
166843db4dbbSBarry Smith }
166943db4dbbSBarry Smith #endif
1670