xref: /petsc/src/sys/objects/pinit.c (revision bc8a24ec4bca81e0b462cc7ae851949ff3b58895)
1*bc8a24ecSBarry 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);
848e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -info_exclude: <null,vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
85e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
86e5c89e4eSSatish Balay #endif
87e5c89e4eSSatish Balay   }
88e5c89e4eSSatish Balay   PetscFunctionReturn(0);
89e5c89e4eSSatish Balay }
90e5c89e4eSSatish Balay 
910f11a792SBarry Smith /*
92945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
9372a42c3cSBarry Smith 
9472a42c3cSBarry Smith    Collective
9572a42c3cSBarry Smith 
9672a42c3cSBarry Smith    Level: advanced
9772a42c3cSBarry Smith 
9895452b02SPatrick Sanan     Notes:
99a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
1000f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
101a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
1020f11a792SBarry Smith 
103a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
1041ea5a559SBarry Smith 
10572a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
1060f11a792SBarry Smith */
107945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
10872a42c3cSBarry Smith {
10972a42c3cSBarry Smith   PetscErrorCode ierr;
11072a42c3cSBarry Smith   int            myargc   = argc;
11172a42c3cSBarry Smith   char           **myargs = args;
11272a42c3cSBarry Smith 
11372a42c3cSBarry Smith   PetscFunctionBegin;
114c52ac268SRichard Tran Mills   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr;
1151ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
116df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
11772a42c3cSBarry Smith   PetscFunctionReturn(ierr);
11872a42c3cSBarry Smith }
11972a42c3cSBarry Smith 
120f0865b08SBarry Smith /*
121a56f64adSBarry Smith       Used by Julia interface to get communicator
122f0865b08SBarry Smith */
123945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
124f0865b08SBarry Smith {
125f0865b08SBarry Smith   PetscFunctionBegin;
126f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
127f0865b08SBarry Smith   PetscFunctionReturn(0);
128f0865b08SBarry Smith }
129f0865b08SBarry Smith 
130e5c89e4eSSatish Balay /*@C
131e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
132e5c89e4eSSatish Balay         the command line arguments.
133e5c89e4eSSatish Balay 
134e5c89e4eSSatish Balay    Collective
135e5c89e4eSSatish Balay 
136e5c89e4eSSatish Balay    Level: advanced
137e5c89e4eSSatish Balay 
138e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
139e5c89e4eSSatish Balay @*/
1407087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
141e5c89e4eSSatish Balay {
142e5c89e4eSSatish Balay   PetscErrorCode ierr;
143e5c89e4eSSatish Balay   int            argc   = 0;
14402c9f0b5SLisandro Dalcin   char           **args = NULL;
145e5c89e4eSSatish Balay 
146e5c89e4eSSatish Balay   PetscFunctionBegin;
1470298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
148e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
149e5c89e4eSSatish Balay }
150e5c89e4eSSatish Balay 
151e5c89e4eSSatish Balay /*@
152e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
153e5c89e4eSSatish Balay 
15493b6d2d1SJed Brown    Level: beginner
155e5c89e4eSSatish Balay 
156e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
157e5c89e4eSSatish Balay @*/
1587087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool  *isInitialized)
159e5c89e4eSSatish Balay {
160e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
16193b6d2d1SJed Brown   return 0;
162e5c89e4eSSatish Balay }
163e5c89e4eSSatish Balay 
164e5c89e4eSSatish Balay /*@
165e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
166e5c89e4eSSatish Balay 
167e5c89e4eSSatish Balay    Level: developer
168e5c89e4eSSatish Balay 
169e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
170e5c89e4eSSatish Balay @*/
1717087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
172e5c89e4eSSatish Balay {
173e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
17493b6d2d1SJed Brown   return 0;
175e5c89e4eSSatish Balay }
176e5c89e4eSSatish Balay 
17795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void);
178e5c89e4eSSatish Balay 
179e5c89e4eSSatish Balay /*
180e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
181e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
182e5c89e4eSSatish Balay */
183367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
184e5c89e4eSSatish Balay 
185367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
186e5c89e4eSSatish Balay {
187e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
188e5c89e4eSSatish Balay 
189e5c89e4eSSatish Balay   PetscFunctionBegin;
190e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
191e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
19241e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
193e5c89e4eSSatish Balay   }
194e5c89e4eSSatish Balay 
195e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
196e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
197e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
198e5c89e4eSSatish Balay   }
199812af9f3SBarry Smith   PetscFunctionReturnVoid();
200e5c89e4eSSatish Balay }
201e5c89e4eSSatish Balay 
202e5c89e4eSSatish Balay /*
203e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
204e5c89e4eSSatish Balay sum of the second entry.
205b693b147SBarry Smith 
20676f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
207367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
208b693b147SBarry Smith there would be no place to store the both needed results.
209e5c89e4eSSatish Balay */
21076ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
211e5c89e4eSSatish Balay {
212e5c89e4eSSatish Balay   PetscErrorCode ierr;
213e5c89e4eSSatish Balay 
214e5c89e4eSSatish Balay   PetscFunctionBegin;
215d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
216d6e4c47cSJed Brown   {
217d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
218367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
219d6e4c47cSJed Brown     *max = work.max;
220d6e4c47cSJed Brown     *sum = work.sum;
221d6e4c47cSJed Brown   }
222d6e4c47cSJed Brown #else
223d6e4c47cSJed Brown   {
224d6e4c47cSJed Brown     PetscMPIInt    size,rank;
225d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
226e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
227e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
228785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
229367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2306ac3741eSJed Brown     *max = work[rank].max;
2316ac3741eSJed Brown     *sum = work[rank].sum;
232e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
233d6e4c47cSJed Brown   }
234d6e4c47cSJed Brown #endif
235e5c89e4eSSatish Balay   PetscFunctionReturn(0);
236e5c89e4eSSatish Balay }
237e5c89e4eSSatish Balay 
238e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
239e5c89e4eSSatish Balay 
240570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
24106a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
242e5c89e4eSSatish Balay 
2438cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
244e5c89e4eSSatish Balay {
245e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
246e5c89e4eSSatish Balay 
247e5c89e4eSSatish Balay   PetscFunctionBegin;
2487c2de775SJed Brown   if (*datatype == MPIU_REAL) {
249e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
250a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2517c2de775SJed Brown   }
2527c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2537c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2547c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
255a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2567c2de775SJed Brown   }
2577c2de775SJed Brown #endif
2587c2de775SJed Brown   else {
2597c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
26041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
261e2e03761SBarry Smith   }
262812af9f3SBarry Smith   PetscFunctionReturnVoid();
263e5c89e4eSSatish Balay }
264e5c89e4eSSatish Balay #endif
265e5c89e4eSSatish Balay 
266570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
267d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
268d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
269d9822059SBarry Smith 
2708cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
271d9822059SBarry Smith {
272d9822059SBarry Smith   PetscInt i,count = *cnt;
273d9822059SBarry Smith 
274d9822059SBarry Smith   PetscFunctionBegin;
2757c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2768c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
277a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2787c2de775SJed Brown   }
2797c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2807c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2817c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2827c2de775SJed Brown     for (i=0; i<count; i++) {
2837c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2847c2de775SJed Brown     }
2857c2de775SJed Brown   }
2867c2de775SJed Brown #endif
2877c2de775SJed Brown   else {
2887c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
28941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2908c764dc5SJose Roman   }
291d9822059SBarry Smith   PetscFunctionReturnVoid();
292d9822059SBarry Smith }
293d9822059SBarry Smith 
2948cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
295d9822059SBarry Smith {
296d9822059SBarry Smith   PetscInt    i,count = *cnt;
297d9822059SBarry Smith 
298d9822059SBarry Smith   PetscFunctionBegin;
2997c2de775SJed Brown   if (*datatype == MPIU_REAL) {
3008c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
301a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
3027c2de775SJed Brown   }
3037c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
3047c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
3057c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
3067c2de775SJed Brown     for (i=0; i<count; i++) {
3077c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
3087c2de775SJed Brown     }
3097c2de775SJed Brown   }
3107c2de775SJed Brown #endif
3117c2de775SJed Brown   else {
3128c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
31341e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
3148c764dc5SJose Roman   }
315d9822059SBarry Smith   PetscFunctionReturnVoid();
316d9822059SBarry Smith }
317d9822059SBarry Smith #endif
318d9822059SBarry Smith 
319480cf27aSJed Brown /*
320480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
321480cf27aSJed Brown 
322ff0e51ddSBarry 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.
323480cf27aSJed Brown 
32412801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
325480cf27aSJed Brown 
326480cf27aSJed Brown */
3278cc058d9SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
328480cf27aSJed Brown {
329480cf27aSJed Brown   PetscErrorCode ierr;
330480cf27aSJed Brown 
331480cf27aSJed Brown   PetscFunctionBegin;
33202c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
33312801b39SBarry Smith   ierr = PetscFree(count_val);CHKERRMPI(ierr);
334480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
335480cf27aSJed Brown }
336480cf27aSJed Brown 
337480cf27aSJed Brown /*
33847435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
339da3039f7SJed Brown   calls MPI_Comm_free().
340da3039f7SJed Brown 
341da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
342480cf27aSJed Brown 
343ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
344480cf27aSJed Brown 
34512801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
346480cf27aSJed Brown 
347480cf27aSJed Brown */
348da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
349480cf27aSJed Brown {
350480cf27aSJed Brown   PetscErrorCode                    ierr;
351b89831e5SBarry Smith   PetscMPIInt                       flg;
352265f3f35SJed Brown   union {MPI_Comm comm; void *ptr;} icomm,ocomm;
353480cf27aSJed Brown 
354480cf27aSJed Brown   PetscFunctionBegin;
35512801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
356ec4fadc2SJed Brown   icomm.ptr = attr_val;
357da3039f7SJed Brown 
35847435625SJed Brown   ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
35912801b39SBarry Smith   if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm");
36012801b39SBarry Smith   if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm");
36147435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
36202c9f0b5SLisandro 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);
363da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
364b89831e5SBarry Smith }
365da3039f7SJed Brown 
366da3039f7SJed Brown /*
36747435625SJed 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.
368da3039f7SJed Brown  */
369da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
370da3039f7SJed Brown {
371da3039f7SJed Brown   PetscErrorCode ierr;
372da3039f7SJed Brown 
373da3039f7SJed Brown   PetscFunctionBegin;
37402c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
375480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
376480cf27aSJed Brown }
377480cf27aSJed Brown 
3785f7487a0SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Shm(MPI_Comm,PetscMPIInt,void *,void *);
37942218b76SBarry Smith 
380951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3818cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3828cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3838cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
384e39fd77fSBarry Smith #endif
385e39fd77fSBarry Smith 
3860f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
38712801b39SBarry Smith 
388eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
389eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3906ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
39102c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
392dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
393e5c89e4eSSatish Balay 
394dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
395051e4cf2SJed Brown {
396051e4cf2SJed Brown   PetscErrorCode ierr;
397051e4cf2SJed Brown 
398051e4cf2SJed Brown   PetscFunctionBegin;
399051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
400a1601671SBarry 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);
401051e4cf2SJed 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);
402051e4cf2SJed Brown   PetscFunctionReturn(0);
403051e4cf2SJed Brown }
404e5c89e4eSSatish Balay 
4052d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
4062d747510SLisandro Dalcin 
4072d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4082d747510SLisandro Dalcin {
4092d747510SLisandro Dalcin   PetscErrorCode ierr;
4102d747510SLisandro Dalcin 
4112d747510SLisandro Dalcin   PetscFunctionBegin;
4122d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
4132d747510SLisandro Dalcin   PetscFunctionReturn(0);
4142d747510SLisandro Dalcin }
4152d747510SLisandro Dalcin 
4162d747510SLisandro Dalcin /*@C
4172d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4182d747510SLisandro Dalcin 
4192d747510SLisandro Dalcin     Not Collective
4202d747510SLisandro Dalcin 
4212d747510SLisandro Dalcin     Input Parameter:
4222d747510SLisandro Dalcin .   len - length of the string name
4232d747510SLisandro Dalcin 
4242d747510SLisandro Dalcin     Output Parameter:
4252d747510SLisandro Dalcin .   name - the name of the running program
4262d747510SLisandro Dalcin 
4272d747510SLisandro Dalcin    Level: advanced
4282d747510SLisandro Dalcin 
4292d747510SLisandro Dalcin     Notes:
4302d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4312d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4322d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4332d747510SLisandro Dalcin @*/
4342d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4352d747510SLisandro Dalcin {
4362d747510SLisandro Dalcin   PetscErrorCode ierr;
4372d747510SLisandro Dalcin 
4382d747510SLisandro Dalcin   PetscFunctionBegin;
4392d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4402d747510SLisandro Dalcin   PetscFunctionReturn(0);
4412d747510SLisandro Dalcin }
4422d747510SLisandro Dalcin 
443e5c89e4eSSatish Balay /*@C
444e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
445e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
446e5c89e4eSSatish Balay 
447e5c89e4eSSatish Balay    Not Collective
448e5c89e4eSSatish Balay 
449e5c89e4eSSatish Balay    Output Parameters:
450e5c89e4eSSatish Balay +  argc - count of number of command line arguments
451e5c89e4eSSatish Balay -  args - the command line arguments
452e5c89e4eSSatish Balay 
453e5c89e4eSSatish Balay    Level: intermediate
454e5c89e4eSSatish Balay 
455e5c89e4eSSatish Balay    Notes:
456e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
457e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
458e5c89e4eSSatish Balay 
459f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
460f177e3b1SBarry Smith 
461793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
462e5c89e4eSSatish Balay 
463e5c89e4eSSatish Balay @*/
4647087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
465e5c89e4eSSatish Balay {
466e5c89e4eSSatish Balay   PetscFunctionBegin;
46717186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
468e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
469e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
470e5c89e4eSSatish Balay   PetscFunctionReturn(0);
471e5c89e4eSSatish Balay }
472e5c89e4eSSatish Balay 
473793721a6SBarry Smith /*@C
474793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
475793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
476793721a6SBarry Smith 
477793721a6SBarry Smith    Not Collective
478793721a6SBarry Smith 
479793721a6SBarry Smith    Output Parameters:
480793721a6SBarry Smith .  args - the command line arguments
481793721a6SBarry Smith 
482793721a6SBarry Smith    Level: intermediate
483793721a6SBarry Smith 
484793721a6SBarry Smith    Notes:
485793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
486793721a6SBarry Smith 
487793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
488793721a6SBarry Smith 
489793721a6SBarry Smith @*/
4907087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
491793721a6SBarry Smith {
492793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
493793721a6SBarry Smith   PetscErrorCode ierr;
494793721a6SBarry Smith 
495793721a6SBarry Smith   PetscFunctionBegin;
49617186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4972d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
498785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
499793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
500793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
501793721a6SBarry Smith   }
5022d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
503793721a6SBarry Smith   PetscFunctionReturn(0);
504793721a6SBarry Smith }
505793721a6SBarry Smith 
506793721a6SBarry Smith /*@C
507793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
508793721a6SBarry Smith 
509793721a6SBarry Smith    Not Collective
510793721a6SBarry Smith 
511793721a6SBarry Smith    Output Parameters:
512793721a6SBarry Smith .  args - the command line arguments
513793721a6SBarry Smith 
514793721a6SBarry Smith    Level: intermediate
515793721a6SBarry Smith 
516793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
517793721a6SBarry Smith 
518793721a6SBarry Smith @*/
5197087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
520793721a6SBarry Smith {
521793721a6SBarry Smith   PetscInt       i = 0;
522793721a6SBarry Smith   PetscErrorCode ierr;
523793721a6SBarry Smith 
524793721a6SBarry Smith   PetscFunctionBegin;
525a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
526793721a6SBarry Smith   while (args[i]) {
527793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
528793721a6SBarry Smith     i++;
529793721a6SBarry Smith   }
530793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
531793721a6SBarry Smith   PetscFunctionReturn(0);
532793721a6SBarry Smith }
533793721a6SBarry Smith 
53411525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
53530befbd2SBarry Smith #include <petscconfiginfo.h>
53630befbd2SBarry Smith 
53795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
53811525c0dSBarry Smith {
53911525c0dSBarry Smith   if (!PetscGlobalRank) {
54030befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
54111525c0dSBarry Smith     int            port;
542ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
54311525c0dSBarry Smith     size_t         applinelen,introlen;
54411525c0dSBarry Smith     PetscErrorCode ierr;
545ffbd1cfbSBarry Smith     char           sawsurl[256];
54611525c0dSBarry Smith 
547c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
54811525c0dSBarry Smith     if (flg) {
54911525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
55011525c0dSBarry Smith 
551c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
55211525c0dSBarry Smith       if (sawslog[0]) {
55311525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
55411525c0dSBarry Smith       } else {
55511525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
55611525c0dSBarry Smith       }
55711525c0dSBarry Smith     }
558c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
55911525c0dSBarry Smith     if (flg) {
56011525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
56111525c0dSBarry Smith     }
562c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
563ffbd1cfbSBarry Smith     if (selectport) {
564ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
565ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
566ffbd1cfbSBarry Smith     } else {
567c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
56811525c0dSBarry Smith       if (flg) {
56911525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
57011525c0dSBarry Smith       }
571ffbd1cfbSBarry Smith     }
572c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
57311525c0dSBarry Smith     if (flg) {
57411525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
57511525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5769c1e0ce8SBarry Smith     } else {
577c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5789c1e0ce8SBarry Smith       if (flg) {
5793c01dfcfSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
5809c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5819c1e0ce8SBarry Smith       }
58211525c0dSBarry Smith     }
583c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
58411525c0dSBarry Smith     if (flg2) {
58511525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
58611525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
58711525c0dSBarry Smith       ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr);
58811525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
58911525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
59043da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
59111525c0dSBarry Smith     }
59211525c0dSBarry Smith     ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr);
59311525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
59411525c0dSBarry Smith     introlen   = 4096 + applinelen;
59530a8c9c0SSurtai Han     applinelen += 1024;
59611525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
59711525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
59811525c0dSBarry Smith 
59911525c0dSBarry Smith     if (rootlocal) {
60011525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
60111525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
60211525c0dSBarry Smith     }
60376a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
60411525c0dSBarry Smith     if (rootlocal && help) {
605928bb9adSStefano 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);
60611525c0dSBarry Smith     } else if (help) {
607928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
60811525c0dSBarry Smith     } else {
609928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
61011525c0dSBarry Smith     }
611b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
61230befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
61311525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
614a8d69d7bSBarry 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"
615df62c222SBarry 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"
616928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
61743da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
61811525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
61911525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
620ffbd1cfbSBarry Smith     if (selectport) {
621aa573868SBarry Smith       PetscBool silent;
6227d812c46SBarry Smith 
6237d812c46SBarry Smith       ierr = SAWs_Initialize();
6247d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6257d812c46SBarry Smith       while (ierr) {
6267d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6277d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6287d812c46SBarry Smith         ierr = SAWs_Initialize();
6297d812c46SBarry Smith       }
6307d812c46SBarry Smith 
631aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
632aa573868SBarry Smith       if (!silent) {
633ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
634ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
635ffbd1cfbSBarry Smith       }
6367d812c46SBarry Smith     } else {
6377d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
638aa573868SBarry Smith     }
6390af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6400af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6410af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6420af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6430af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
64411525c0dSBarry Smith   }
64511525c0dSBarry Smith   PetscFunctionReturn(0);
64611525c0dSBarry Smith }
64711525c0dSBarry Smith #endif
64811525c0dSBarry Smith 
6494dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6504dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6514dfee713SSatish Balay {
6524dfee713SSatish Balay   PetscFunctionBegin;
6534dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6544dfee713SSatish Balay     /* see MPI.py for details on this bug */
6554dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6564dfee713SSatish Balay #endif
6574dfee713SSatish Balay   PetscFunctionReturn(0);
6584dfee713SSatish Balay }
6594dfee713SSatish Balay 
660a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
661a56f64adSBarry Smith #include <adios.h>
66222580e64SBarry Smith #include <adios_read.h>
6637b56e58cSSatish Balay int64_t Petsc_adios_group;
664a56f64adSBarry Smith #endif
66555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
66655d657eeSBarry Smith #include <adios2_c.h>
66755d657eeSBarry Smith #endif
668cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
669cd1b99a6SBarry Smith #include <omp.h>
670f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
671cd1b99a6SBarry Smith #endif
672a56f64adSBarry Smith 
673*bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H)
674*bc8a24ecSBarry Smith #include <dlfcn.h>
675*bc8a24ecSBarry Smith #endif
676*bc8a24ecSBarry Smith 
677e5c89e4eSSatish Balay /*@C
678e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
679e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
680e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
681e5c89e4eSSatish Balay    your program -- usually the very first line!
682e5c89e4eSSatish Balay 
683e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
684e5c89e4eSSatish Balay 
685e5c89e4eSSatish Balay    Input Parameters:
686e5c89e4eSSatish Balay +  argc - count of number of command line arguments
687e5c89e4eSSatish Balay .  args - the command line arguments
6880298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
689fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6900298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
691e5c89e4eSSatish Balay 
69205827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
69305827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
69405827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
69505827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
69605827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
697e5c89e4eSSatish Balay 
698e5c89e4eSSatish Balay    Options Database Keys:
6997ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
7007ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
701e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
7028a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
7038a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
7048a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
7058a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
7068a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
707e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
708e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
709e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
710e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
71179dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
71279dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
71379dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
714aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
71592f119d6SBarry 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
71692f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
71792f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
71892f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
719e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
720e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
721e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
722e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
723e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7240841954dSBarry Smith -  -memory_view - Print memory usage at end of run
725e5c89e4eSSatish Balay 
726e5c89e4eSSatish Balay    Options Database Keys for Profiling:
727a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
728495fc317SBarry Smith +  -info <optional filename> - Prints verbose information to the screen
729495fc317SBarry Smith .  -info_exclude <null,vec,mat,pc,ksp,snes,ts> - Excludes some of the verbose messages
730217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
731217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
732495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
733e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7349a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
73579dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
7369a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
737495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
738f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
739495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
740495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
741c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
74287c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
743c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
744495fc317SBarry Smith 
745609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
746e5c89e4eSSatish Balay 
747ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
748ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
749ffbd1cfbSBarry 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
750ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
751ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
752ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
753ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
754ffbd1cfbSBarry Smith 
755e5c89e4eSSatish Balay    Environmental Variables:
756e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
757e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
758e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
759e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
760e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
761e5c89e4eSSatish Balay 
762e5c89e4eSSatish Balay 
763e5c89e4eSSatish Balay    Level: beginner
764e5c89e4eSSatish Balay 
765e5c89e4eSSatish Balay    Notes:
766e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
767e5c89e4eSSatish Balay    it before PetscInitialize().
768e5c89e4eSSatish Balay 
769e5c89e4eSSatish Balay    Fortran Version:
770e5c89e4eSSatish Balay    In Fortran this routine has the format
771e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
772e5c89e4eSSatish Balay 
773e5c89e4eSSatish Balay +   ierr - error return code
7740eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
775fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
776e5c89e4eSSatish Balay 
777e5c89e4eSSatish Balay    Important Fortran Note:
7780eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7790298fd71SBarry Smith    null character string; you CANNOT just use NULL as
780a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
781e5c89e4eSSatish Balay 
78201cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
78301cb0274SBarry Smith    calling PetscInitialize().
784e5c89e4eSSatish Balay 
78501cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
786e5c89e4eSSatish Balay 
787e5c89e4eSSatish Balay @*/
7887087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
789e5c89e4eSSatish Balay {
790e5c89e4eSSatish Balay   PetscErrorCode ierr;
7914bb5149bSJed Brown   PetscMPIInt    flag, size;
79219c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
793e5c89e4eSSatish Balay   char           hostname[256];
794e5c89e4eSSatish Balay 
795e5c89e4eSSatish Balay   PetscFunctionBegin;
796e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7973d96e996SBarry Smith   /*
7983d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7993d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
8003d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
8013d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
8023d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
8033d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
8043d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
8053d96e996SBarry Smith       listed above and since that time are compatible.
8063d96e996SBarry Smith 
8073d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
8083d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
8093d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
8103d96e996SBarry Smith       to perform the checking.
8113d96e996SBarry Smith 
8123d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
8133d96e996SBarry Smith 
8143d96e996SBarry Smith       Questions:
8153d96e996SBarry Smith 
8163d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
8173d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
8183d96e996SBarry Smith   */
8193d96e996SBarry Smith 
82019c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
82119c5658aSBarry Smith   {
82219c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
82319c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
82419c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8253d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
82619c5658aSBarry Smith #if defined(MPICH_VERSION)
8273d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
82819c5658aSBarry Smith     {
82919c5658aSBarry Smith       char *ver,*lf;
83019c5658aSBarry Smith       flg = PETSC_FALSE;
83119c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
83219c5658aSBarry Smith       if (ver) {
83319c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
83419c5658aSBarry Smith         if (lf) {
83519c5658aSBarry Smith           *lf = 0;
83619c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
83719c5658aSBarry Smith         }
83819c5658aSBarry Smith       }
83919c5658aSBarry Smith       if (!flg) {
84019c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8413d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
84219c5658aSBarry Smith       }
84319c5658aSBarry Smith     }
8443d96e996SBarry Smith #endif
8453d96e996SBarry 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?) */
84619c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
84719c5658aSBarry Smith     {
84819c5658aSBarry Smith       char *ver,bs[32],*bsf;
84919c5658aSBarry Smith       flg = PETSC_FALSE;
85019c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
85119c5658aSBarry Smith       if (ver) {
8522e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
85319c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
85419c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
85519c5658aSBarry Smith       }
85619c5658aSBarry Smith       if (!flg) {
85719c5658aSBarry 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);
8583d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
85919c5658aSBarry Smith       }
86019c5658aSBarry Smith     }
86119c5658aSBarry Smith #endif
86219c5658aSBarry Smith   }
86319c5658aSBarry Smith #endif
86419c5658aSBarry Smith 
865*bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM)
866*bc8a24ecSBarry Smith   {
867*bc8a24ecSBarry Smith     PetscInt cnt = 0;
868*bc8a24ecSBarry 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 */
869*bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++;
870*bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++;
871*bc8a24ecSBarry Smith     if (cnt > 1) {
872*bc8a24ecSBarry Smith       fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
873*bc8a24ecSBarry Smith       return PETSC_ERR_MPI_LIB_INCOMP;
874*bc8a24ecSBarry Smith     }
875*bc8a24ecSBarry Smith   }
876*bc8a24ecSBarry Smith #endif
877*bc8a24ecSBarry Smith 
878ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
879d89683f4Sbcordonn   PETSC_STDOUT = stdout;
880ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
881e5c89e4eSSatish Balay 
8820c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8830c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8840c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8850c30907bSSatish Balay #endif
8860c30907bSSatish Balay 
8874416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
888e5c89e4eSSatish Balay 
889e5c89e4eSSatish Balay   /*
890e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
891e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
892e5c89e4eSSatish Balay   */
893e5c89e4eSSatish Balay   if (argc && *argc) {
894e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
895e5c89e4eSSatish Balay   } else {
896e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
897e5c89e4eSSatish Balay   }
898e5c89e4eSSatish Balay 
899e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
900e5c89e4eSSatish Balay   if (!flag) {
901e32f2f54SBarry 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");
9024dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
9035e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
9045e765c61SJed Brown     {
9055e765c61SJed Brown       PetscMPIInt provided;
9065e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
9075e765c61SJed Brown     }
9085e765c61SJed Brown #else
909e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
9105e765c61SJed Brown #endif
911e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
912e5c89e4eSSatish Balay   }
9134dfee713SSatish Balay 
914e5c89e4eSSatish Balay   if (argc && args) {
915e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
916e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
917e5c89e4eSSatish Balay   }
918e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
9195ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
9205ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
9215ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
922ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
923e5c89e4eSSatish Balay 
924a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
925d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
926e5c89e4eSSatish Balay 
9270f9be574SPeter Hill   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
92812801b39SBarry Smith     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
92912801b39SBarry Smith     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
9300f9be574SPeter Hill   }
93112801b39SBarry Smith 
932e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
933e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
934e5c89e4eSSatish Balay 
935e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
936e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
937e5c89e4eSSatish Balay 
9388ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9398ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9407cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
941e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
942e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
943e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
944e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
945e316c87fSJed Brown #endif
946e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9478ad47952SJed Brown 
948e5c89e4eSSatish Balay   /*
949e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
950e5c89e4eSSatish Balay      shared libraries the constructors for global variables
951e5c89e4eSSatish Balay      are not called; at least on IRIX.
952e5c89e4eSSatish Balay   */
953886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
954e5c89e4eSSatish Balay   {
955252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
95650f81f78SJed Brown     PetscComplex ic(0.0,1.0);
957e5c89e4eSSatish Balay     PETSC_i = ic;
958252ecd31SSatish Balay #else
95950f81f78SJed Brown     PETSC_i = _Complex_I;
960b7940d39SSatish Balay #endif
961762437b8SSatish Balay   }
962762437b8SSatish Balay 
9632c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
964e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
965500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
966500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
967500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9682c876bd9SBarry Smith #endif
969886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
970e5c89e4eSSatish Balay 
971e5c89e4eSSatish Balay   /*
972e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
973e5c89e4eSSatish Balay      half of the entries and maxes the second half.
974e5c89e4eSSatish Balay   */
975367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
976e5c89e4eSSatish Balay 
977ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
978c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
979c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9807c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9818c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9828c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9838c764dc5SJose Roman #endif
984d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
985d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
986570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
987570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
988570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
989570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
990570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
991c90a1750SBarry Smith #endif
992c90a1750SBarry Smith 
993570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
994cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
995cca4cb22SSatish Balay #endif
996cca4cb22SSatish Balay 
997e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
998e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
999e5c89e4eSSatish Balay 
100040df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1001e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
1002e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
100344041f26SJed Brown #endif
1004e5c89e4eSSatish Balay 
1005ec957eceSBarry Smith 
1006e5c89e4eSSatish Balay   /*
1007480cf27aSJed Brown      Attributes to be set on PETSc communicators
1008480cf27aSJed Brown   */
100912801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
101012801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
101112801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
10125f7487a0SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Shm,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
1013480cf27aSJed Brown 
1014480cf27aSJed Brown   /*
1015e8fb0fc0SBarry Smith      Build the options database
1016e5c89e4eSSatish Balay   */
1017c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
1018e5c89e4eSSatish Balay 
1019703f3eceSBarry Smith   /* call a second time so it can look in the options database */
1020703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
10216dc8fec2Sbcordonn 
1022e5c89e4eSSatish Balay   /*
1023e5c89e4eSSatish Balay      Print main application help message
1024e5c89e4eSSatish Balay   */
10252d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
1026e5c89e4eSSatish Balay   if (help && flg) {
1027e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
1028e5c89e4eSSatish Balay   }
1029e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
1030e5c89e4eSSatish Balay 
1031d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
1032d45a07a7SBarry Smith 
1033e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
103411525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1035f4202a44SBarry Smith #endif
1036f4202a44SBarry Smith 
1037e5c89e4eSSatish Balay   /*
1038e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1039e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1040e5c89e4eSSatish Balay   */
1041e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1042e5c89e4eSSatish Balay 
1043e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
104402c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1045e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
104602c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1047cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1048101491d0SBarry Smith   {
1049101491d0SBarry Smith     PetscBool omp_view_flag;
1050cd1b99a6SBarry Smith     char      *threads = getenv("OMP_NUM_THREADS");
1051e5c89e4eSSatish Balay 
1052cd1b99a6SBarry Smith    if (threads) {
105302c9f0b5SLisandro Dalcin      ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
1054f90b075cSBarry Smith      (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
1055cd1b99a6SBarry Smith    } else {
1056cd1b99a6SBarry Smith #define NMAX  10000
1057101491d0SBarry Smith      int          i;
1058cd1b99a6SBarry Smith       PetscScalar *x;
1059cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1060cd1b99a6SBarry Smith #pragma omp parallel for
1061cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1062cd1b99a6SBarry Smith         x[i] = 0.0;
1063f90b075cSBarry Smith         PetscNumOMPThreads  = (PetscInt) omp_get_num_threads();
1064cd1b99a6SBarry Smith       }
1065cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
106602c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1067101491d0SBarry Smith     }
1068101491d0SBarry Smith     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
1069f90b075cSBarry 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);
1070101491d0SBarry Smith     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
1071101491d0SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1072101491d0SBarry Smith     if (flg) {
107302c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1074f90b075cSBarry Smith       omp_set_num_threads((int)PetscNumOMPThreads);
1075101491d0SBarry Smith     }
1076101491d0SBarry Smith     if (omp_view_flag) {
1077f90b075cSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
1078cd1b99a6SBarry Smith     }
1079cd1b99a6SBarry Smith   }
1080cd1b99a6SBarry Smith #endif
1081e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr);
1082ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1083c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1084ef6c6fedSBoyana Norris 
1085951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1086e39fd77fSBarry Smith   /*
1087e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1088e39fd77fSBarry Smith 
1089e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1090e39fd77fSBarry Smith   */
109130815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
10920298fd71SBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
109330815ce0SLisandro Dalcin   }
1094951e3c8eSBarry Smith #endif
1095e39fd77fSBarry Smith 
109641c0b4b3SShri Abhyankar   /*
109741c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
109841c0b4b3SShri Abhyankar   */
1099ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1100e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1101e1167bb9SShri Abhyankar #endif
1102e1167bb9SShri Abhyankar 
11032d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
11042d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
11052d53ad75SBarry Smith #endif
11062d53ad75SBarry Smith 
11075e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
110887c3beb6SLisandro Dalcin   {
110987c3beb6SLisandro Dalcin     PetscViewer viewer;
111022e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
11115e71baefSBarry Smith     if (flg) {
11125e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
111387c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
111487c3beb6SLisandro Dalcin     }
11155e71baefSBarry Smith   }
11165e71baefSBarry Smith #endif
1117dff31646SBarry Smith 
111887c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
111987c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
112087c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
112187c3beb6SLisandro Dalcin 
1122a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1123a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
11247b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1125a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
11269fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1127a56f64adSBarry Smith #endif
112855d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
112955d657eeSBarry Smith #endif
1130a56f64adSBarry Smith 
1131301d30feSBarry Smith   /*
1132301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1133301d30feSBarry Smith   */
1134301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
11352db0e300SLisandro Dalcin 
11362db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
11372db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1138eb02082bSJunchao Zhang 
1139eb02082bSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1140eb02082bSJunchao Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-use_gpu_aware_mpi",&use_gpu_aware_mpi,NULL);CHKERRQ(ierr);
1141eb02082bSJunchao Zhang #if defined(PETSC_HAVE_OPENMPI_MAJOR) && defined(MPIX_CUDA_AWARE_SUPPORT)
1142eb02082bSJunchao Zhang   /* OpenMPI supports compile time and runtime cuda support checking */
1143eb02082bSJunchao 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");
1144eb02082bSJunchao Zhang #endif
1145eb02082bSJunchao Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-sf_use_default_cuda_stream",&sf_use_default_cuda_stream,NULL);CHKERRQ(ierr);
1146eb02082bSJunchao Zhang #endif
1147eb02082bSJunchao Zhang 
1148301d30feSBarry Smith   PetscFunctionReturn(0);
1149e5c89e4eSSatish Balay }
1150e5c89e4eSSatish Balay 
11514097062eSBarry Smith #if defined(PETSC_USE_LOG)
115295c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1153ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1154ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
115595c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11564097062eSBarry Smith #endif
1157e5c89e4eSSatish Balay 
1158008a6e76SBarry Smith /*
1159008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1160008a6e76SBarry Smith */
1161008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1162008a6e76SBarry Smith {
1163008a6e76SBarry Smith   PetscErrorCode ierr;
1164008a6e76SBarry Smith 
1165008a6e76SBarry Smith   PetscFunctionBegin;
1166008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1167008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1168008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1169008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1170008a6e76SBarry Smith #endif
1171008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1172008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1173008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1174008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1175008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1176008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1177008a6e76SBarry Smith #endif
1178008a6e76SBarry Smith 
1179008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1180008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1181008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1182008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1183008a6e76SBarry Smith #endif
1184008a6e76SBarry Smith #endif
1185008a6e76SBarry Smith 
1186008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1187008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1188008a6e76SBarry Smith #endif
1189008a6e76SBarry Smith 
1190008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
119140df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1192008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1193008a6e76SBarry Smith #endif
1194008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1195008a6e76SBarry Smith   PetscFunctionReturn(0);
1196008a6e76SBarry Smith }
1197008a6e76SBarry Smith 
1198e5c89e4eSSatish Balay /*@C
1199e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1200e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1201e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1202e5c89e4eSSatish Balay 
1203e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1204e5c89e4eSSatish Balay 
1205e5c89e4eSSatish Balay    Options Database Keys:
120626a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1207e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
12087eb1d149SBarry 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
1209e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
121092f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1211e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
121292f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1213e5c89e4eSSatish Balay 
1214e5c89e4eSSatish Balay    Level: beginner
1215e5c89e4eSSatish Balay 
1216e5c89e4eSSatish Balay    Note:
1217e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1218e5c89e4eSSatish Balay 
121988c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1220e5c89e4eSSatish Balay @*/
12217087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1222e5c89e4eSSatish Balay {
1223e5c89e4eSSatish Balay   PetscErrorCode ierr;
12244bb5149bSJed Brown   PetscMPIInt    rank;
1225a8d2bbe5SBarry Smith   PetscInt       nopt;
12262bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1227dff31646SBarry Smith   PetscBool      flg;
122810463e74SBarry Smith #if defined(PETSC_USE_LOG)
122910463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
123010463e74SBarry Smith #endif
1231e5c89e4eSSatish Balay 
1232e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
12334b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
12343cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1235e5c89e4eSSatish Balay   }
12363cbbc5ffSBarry Smith   PetscFunctionBegin;
12370298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1238b022a5c1SBarry Smith 
12391f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1240a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
124122580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1242a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1243a56f64adSBarry Smith #endif
124455d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
124555d657eeSBarry Smith #endif
1246c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1247dff31646SBarry Smith   if (flg) {
12481f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
12491f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
12501f817a21SBarry Smith 
1251c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
12521f817a21SBarry Smith     if (filename[0]) {
12531f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
12541f817a21SBarry Smith     }
1255dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1256dff31646SBarry Smith     cits[0] = 0;
1257dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
12581f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
12591f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12601f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
12611f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12621f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1263dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1264dff31646SBarry Smith   }
1265dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1266dff31646SBarry Smith 
1267c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
126804102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
126904102261SBarry Smith   {
127004102261SBarry Smith     PetscInt nmax = 2;
127104102261SBarry Smith     char     **buffs;
127204102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1273c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
127404102261SBarry Smith     if (flg1) {
127504102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
127604102261SBarry Smith       if (nmax == 1) {
127704102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
127804102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
127904102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
128004102261SBarry Smith       }
128104102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
128204102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
128304102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
128404102261SBarry Smith     }
128504102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
128604102261SBarry Smith   }
1287763ec7b1SBarry Smith   {
1288763ec7b1SBarry Smith     PetscInt nmax = 2;
1289763ec7b1SBarry Smith     char     **buffs;
1290763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1291763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1292763ec7b1SBarry Smith     if (flg1) {
1293763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1294763ec7b1SBarry Smith       if (nmax == 1) {
1295763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1296763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1297763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1298763ec7b1SBarry Smith       }
1299763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1300763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1301763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1302763ec7b1SBarry Smith     }
1303763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1304763ec7b1SBarry Smith   }
130504102261SBarry Smith #endif
130667234432SDmitry Karpeev   /*
130767234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
130867234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
130967234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
131067234432SDmitry Karpeev   */
131167234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
131204102261SBarry Smith 
13132d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
13142d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
13152d53ad75SBarry Smith #endif
13162d53ad75SBarry Smith 
13172d53ad75SBarry Smith 
1318e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1319dff31646SBarry Smith   flg = PETSC_FALSE;
1320c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1321d5649816SBarry Smith   if (flg) {
1322e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1323d5649816SBarry Smith   }
1324d5649816SBarry Smith #endif
1325d5649816SBarry Smith 
1326681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1327681455b2SBarry Smith   flg1 = PETSC_FALSE;
1328c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1329681455b2SBarry Smith   if (flg1) {
1330681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1331681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1332681455b2SBarry Smith   }
1333681455b2SBarry Smith #endif
1334681455b2SBarry Smith 
133567584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1336c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1337e5c89e4eSSatish Balay   if (!flg2) {
133890d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1339c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1340e5c89e4eSSatish Balay   }
1341e5c89e4eSSatish Balay   if (flg2) {
13420841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1343e5c89e4eSSatish Balay   }
134467584ceeSBarry Smith #endif
1345e5c89e4eSSatish Balay 
1346e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
134790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1348c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1349e5c89e4eSSatish Balay   if (flg1) {
1350e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1351205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1352e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1353e5c89e4eSSatish Balay   }
1354e5c89e4eSSatish Balay #endif
1355e5c89e4eSSatish Balay 
1356e5c89e4eSSatish Balay 
1357e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1358e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1359e5c89e4eSSatish Balay   mname[0] = 0;
1360c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1361e5c89e4eSSatish Balay   if (flg1) {
1362e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1363e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1364e5c89e4eSSatish Balay   }
1365e5c89e4eSSatish Balay #endif
1366356e5837SBarry Smith #endif
1367a297a907SKarl Rupp 
1368dd710f27SBarry Smith   /*
1369dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1370dd710f27SBarry Smith   */
1371dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1372dd710f27SBarry Smith 
1373356e5837SBarry Smith #if defined(PETSC_USE_LOG)
137487c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1375f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
137687c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
137787c3beb6SLisandro Dalcin 
1378356e5837SBarry Smith   mname[0] = 0;
1379c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1380e5c89e4eSSatish Balay   if (flg1) {
138191eabc43SBarry Smith     PetscViewer viewer;
138220a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
138391eabc43SBarry Smith     if (mname[0]) {
138491eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
138591eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13866bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
138733f85c2fSBarry Smith     } else {
138833f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13899a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
139033f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13919a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
139233f85c2fSBarry Smith     }
1393e5c89e4eSSatish Balay   }
1394a297a907SKarl Rupp 
1395dd710f27SBarry Smith   /*
1396dd710f27SBarry Smith      Free any objects created by the last block of code.
1397dd710f27SBarry Smith   */
1398dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1399dd710f27SBarry Smith 
1400dd710f27SBarry Smith   mname[0] = 0;
1401c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1402c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
14037ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1404e5c89e4eSSatish Balay #endif
140510463e74SBarry Smith 
140633f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
140710463e74SBarry Smith 
140890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1409c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1410e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
141190d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1412c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1413e5c89e4eSSatish Balay   if (flg1) {
1414e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1415e5c89e4eSSatish Balay   }
141690d69ab7SBarry Smith   flg1 = PETSC_FALSE;
141790d69ab7SBarry Smith   flg2 = PETSC_FALSE;
14188bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1419c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1420c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1421c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1422e4c476e2SSatish Balay 
1423e5c89e4eSSatish Balay   if (flg2) {
1424be56827dSJed Brown     PetscViewer viewer;
142502ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
142602ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1427c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1428be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1429e5c89e4eSSatish Balay   }
1430e5c89e4eSSatish Balay 
1431e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1432c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1433c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1434e5c89e4eSSatish Balay 
143533fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1436c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
14373de2bfdfSBarry Smith #if defined(PETSC_USE_DEBUG)
14383de2bfdfSBarry Smith   if (!flg1) flg3 = PETSC_TRUE;
14393de2bfdfSBarry Smith #endif
1440e5c89e4eSSatish Balay   if (flg3) {
14413de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1442be56827dSJed Brown       PetscViewer viewer;
144302ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
144402ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1445c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1446be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1447e5c89e4eSSatish Balay     }
14483de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
14493de2bfdfSBarry Smith     if (nopt) {
14503de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
14513de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
14523de2bfdfSBarry Smith       if (nopt == 1) {
1453e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1454e5c89e4eSSatish Balay       } else {
14557582186dSLisandro Dalcin         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1456e5c89e4eSSatish Balay       }
14573de2bfdfSBarry Smith     } else if (flg3 && flg1) {
14583de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1459df12ba86SBarry Smith     }
1460c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1461e5c89e4eSSatish Balay   }
1462e5c89e4eSSatish Balay 
1463e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1464d45a07a7SBarry Smith   if (!PetscGlobalRank) {
146587f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
146616ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1467d45a07a7SBarry Smith   }
1468ec957eceSBarry Smith #endif
1469ec957eceSBarry Smith 
14704097062eSBarry Smith #if defined(PETSC_USE_LOG)
147110463e74SBarry Smith   /*
1472dbc8283eSBarry Smith        List all objects the user may have forgot to free
14732eff7a51SBarry Smith   */
147405df10baSBarry Smith   if (PetscObjectsLog) {
1475c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1476a64a8e02SBarry Smith     if (flg1) {
1477a64a8e02SBarry Smith       MPI_Comm local_comm;
14787eb1d149SBarry Smith       char     string[64];
1479a64a8e02SBarry Smith 
1480c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1481a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1482a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14837eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1484a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1485a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14860a1571b3SBarry Smith     }
148705df10baSBarry Smith   }
14884097062eSBarry Smith #endif
14894097062eSBarry Smith 
14904097062eSBarry Smith #if defined(PETSC_USE_LOG)
1491dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1492dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1493a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14944097062eSBarry Smith #endif
14952eff7a51SBarry Smith 
149633f85c2fSBarry Smith   /*
149733f85c2fSBarry Smith      Destroy any packages that registered a finalize
149833f85c2fSBarry Smith   */
149933f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
150033f85c2fSBarry Smith 
1501101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1502fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1503101409b8SToby Isaac #endif
1504101409b8SToby Isaac 
150533f85c2fSBarry Smith   /*
150648dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
150748dd1dffSBarry Smith 
150837e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
150948dd1dffSBarry Smith   */
151037e93019SBarry Smith 
15114028d114SSatish Balay   if (petsc_history) {
1512f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
151302c9f0b5SLisandro Dalcin     petsc_history = NULL;
1514e5c89e4eSSatish Balay   }
15159de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1516e5c89e4eSSatish Balay 
15170298fd71SBarry Smith   ierr = PetscInfoAllow(PETSC_FALSE,NULL);CHKERRQ(ierr);
1518e5c89e4eSSatish Balay 
151967584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
152092f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1521e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
152292f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1523e5c89e4eSSatish Balay     FILE *fd;
1524ed9cf6e9SBarry Smith     int  err;
1525e5c89e4eSSatish Balay 
1526dc92acbaSJed Brown     flg2 = PETSC_FALSE;
152792f119d6SBarry Smith     flg3 = PETSC_FALSE;
15288bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
152992f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);
1530dc92acbaSJed Brown #endif
153192f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
153292f119d6SBarry Smith     fname[0] = 0;
153392f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1534e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1535e5c89e4eSSatish Balay 
15362e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1537e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1538e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1539ed9cf6e9SBarry Smith       err  = fclose(fd);
1540e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
154192f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1542e5c89e4eSSatish Balay       MPI_Comm local_comm;
1543e5c89e4eSSatish Balay 
1544e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1545e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1546e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1547e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1548e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1549e5c89e4eSSatish Balay     }
1550e5c89e4eSSatish Balay     fname[0] = 0;
155192f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,250,&flg1);CHKERRQ(ierr);
1552e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1553e5c89e4eSSatish Balay 
155492f119d6SBarry Smith       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
155592f119d6SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
155692f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1557ed9cf6e9SBarry Smith       err  = fclose(fd);
1558e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
155992f119d6SBarry Smith     } else if (flg1) {
156092f119d6SBarry Smith       MPI_Comm local_comm;
156192f119d6SBarry Smith 
156292f119d6SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
156392f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
156492f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
156592f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
156692f119d6SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1567e5c89e4eSSatish Balay     }
1568e5c89e4eSSatish Balay   }
156967584ceeSBarry Smith #endif
157020e2c332SMatthew G. Knepley 
15715486ca60SMatthew G. Knepley   /*
15725486ca60SMatthew G. Knepley      Close any open dynamic libraries
15735486ca60SMatthew G. Knepley   */
15745486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15755486ca60SMatthew G. Knepley 
1576e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15774416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1578e5c89e4eSSatish Balay 
1579e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
158002c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1581e5c89e4eSSatish Balay 
1582008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1583e5c89e4eSSatish Balay 
1584dbc8283eSBarry Smith   /*
1585efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1586efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1587efb80d3cSBarry Smith 
1588efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1589efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1590dbc8283eSBarry Smith  */
1591b770b1f6SSatish Balay   {
1592dbc8283eSBarry Smith     PetscCommCounter *counter;
1593dbc8283eSBarry Smith     PetscMPIInt      flg;
1594dbc8283eSBarry Smith     MPI_Comm         icomm;
1595265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
159647435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1597dbc8283eSBarry Smith     if (flg) {
1598265f3f35SJed Brown       icomm = ucomm.comm;
159947435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1600dbc8283eSBarry 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");
1601dbc8283eSBarry Smith 
160247435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
160347435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1604efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1605dbc8283eSBarry Smith     }
160647435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1607dbc8283eSBarry Smith     if (flg) {
1608265f3f35SJed Brown       icomm = ucomm.comm;
160947435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1610dbc8283eSBarry 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");
1611dbc8283eSBarry Smith 
161247435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
161347435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1614efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1615dbc8283eSBarry Smith     }
1616b770b1f6SSatish Balay   }
1617dbc8283eSBarry Smith 
161847435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
161947435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
162047435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
16215f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1622480cf27aSJed Brown 
16235ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
16245ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
16255ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1626ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1627ef19f930SBarry Smith 
1628e5c89e4eSSatish Balay   if (PetscBeganMPI) {
162999608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
163099b1327fSBarry Smith     PetscMPIInt flag;
163199b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1632e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
163399608316SBarry Smith #endif
1634e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1635e5c89e4eSSatish Balay   }
1636e5c89e4eSSatish Balay /*
1637e5c89e4eSSatish Balay 
1638e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1639e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1640e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1641e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1642e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16430e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1644e5c89e4eSSatish Balay    memory was not freed.
1645e5c89e4eSSatish Balay 
1646e5c89e4eSSatish Balay */
16471d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1648a297a907SKarl Rupp 
1649e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1650e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
16513db9a53dSBarry Smith   PetscFunctionReturn(0);
1652e5c89e4eSSatish Balay }
1653e5c89e4eSSatish Balay 
165443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16558cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
165643db4dbbSBarry Smith {
165743db4dbbSBarry Smith   if (*a == *b) return 1;
165843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
165943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
166043db4dbbSBarry Smith   return 0;
166143db4dbbSBarry Smith }
1662a70650f6SBarry Smith #endif
166343db4dbbSBarry Smith 
166443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16658cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
166643db4dbbSBarry Smith {
166743db4dbbSBarry Smith   if (*a == *b) return 1;
166843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
166943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
167043db4dbbSBarry Smith   return 0;
167143db4dbbSBarry Smith }
167243db4dbbSBarry Smith #endif
1673