xref: /petsc/src/sys/objects/pinit.c (revision 05342407404a810c6afec3612a6850394c3ef881)
1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file defines the initialization of PETSc, including PetscInitialize()
4e5c89e4eSSatish Balay */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>        /*I  "petscsys.h"   I*/
6022afb99SBarry Smith #include <petscvalgrind.h>
7665c2dedSJed Brown #include <petscviewer.h>
88101f56cSMatthew Knepley 
9a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
11a9f03627SSatish Balay #endif
12f2d66bcaSShri Abhyankar 
132d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1495c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData;
152d53ad75SBarry Smith PetscFPT PetscFPTData = 0;
162d53ad75SBarry Smith #endif
172d53ad75SBarry Smith 
18a6790183SBarry Smith #if defined(PETSC_HAVE_SAWS)
1916ad0300SBarry Smith #include <petscviewersaws.h>
20a6790183SBarry Smith #endif
21eb02082bSJunchao Zhang 
22e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/
23e5c89e4eSSatish Balay 
2495c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history;
25e5c89e4eSSatish Balay 
2695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
2795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
2895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
2995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
3095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
3195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);
320069ddf5SShri Abhyankar 
33e5c89e4eSSatish Balay /* user may set this BEFORE calling PetscInitialize() */
34e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
35e5c89e4eSSatish Balay 
36480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
37480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
38480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
395f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
40480cf27aSJed Brown 
41e5c89e4eSSatish Balay /*
42e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
43e5c89e4eSSatish Balay */
4402c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
4502c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
46e5c89e4eSSatish Balay 
47ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
48ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
490f8e0872SSatish Balay 
50a2f94806SJed Brown PetscInt PetscHotRegionDepth;
51a2f94806SJed Brown 
52b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
53b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
54b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
55b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
56b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
57b22622e2STadeu Manoel #endif
58b22622e2STadeu Manoel 
59e5c89e4eSSatish Balay /*
60e5c89e4eSSatish Balay        Checks the options database for initializations related to the
61e5c89e4eSSatish Balay     PETSc components
62e5c89e4eSSatish Balay */
6395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode  PetscOptionsCheckInitial_Components(void)
64e5c89e4eSSatish Balay {
652d747510SLisandro Dalcin   PetscBool      flg;
66e5c89e4eSSatish Balay   PetscErrorCode ierr;
67e5c89e4eSSatish Balay 
68e5c89e4eSSatish Balay   PetscFunctionBegin;
692d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
702d747510SLisandro Dalcin   if (flg) {
71e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
72e8e7597cSSatish Balay     MPI_Comm comm = PETSC_COMM_WORLD;
73e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"------Additional PETSc component options--------\n");CHKERRQ(ierr);
748e81d068SLisandro Dalcin     ierr = (*PetscHelpPrintf)(comm," -log_exclude: <vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr);
75e5c89e4eSSatish Balay     ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr);
76e5c89e4eSSatish Balay #endif
77e5c89e4eSSatish Balay   }
78e5c89e4eSSatish Balay   PetscFunctionReturn(0);
79e5c89e4eSSatish Balay }
80e5c89e4eSSatish Balay 
810f11a792SBarry Smith /*
82945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
8372a42c3cSBarry Smith 
8472a42c3cSBarry Smith    Collective
8572a42c3cSBarry Smith 
8672a42c3cSBarry Smith    Level: advanced
8772a42c3cSBarry Smith 
8895452b02SPatrick Sanan     Notes:
89a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
900f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
91a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
920f11a792SBarry Smith 
93a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
941ea5a559SBarry Smith 
9572a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
960f11a792SBarry Smith */
97945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
9872a42c3cSBarry Smith {
9972a42c3cSBarry Smith   PetscErrorCode ierr;
10072a42c3cSBarry Smith   int            myargc   = argc;
10172a42c3cSBarry Smith   char           **myargs = args;
10272a42c3cSBarry Smith 
10372a42c3cSBarry Smith   PetscFunctionBegin;
104c52ac268SRichard Tran Mills   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr;
1051ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
106df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
10772a42c3cSBarry Smith   PetscFunctionReturn(ierr);
10872a42c3cSBarry Smith }
10972a42c3cSBarry Smith 
110f0865b08SBarry Smith /*
111a56f64adSBarry Smith       Used by Julia interface to get communicator
112f0865b08SBarry Smith */
113945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
114f0865b08SBarry Smith {
115f0865b08SBarry Smith   PetscFunctionBegin;
116f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
117f0865b08SBarry Smith   PetscFunctionReturn(0);
118f0865b08SBarry Smith }
119f0865b08SBarry Smith 
120e5c89e4eSSatish Balay /*@C
121e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
122e5c89e4eSSatish Balay         the command line arguments.
123e5c89e4eSSatish Balay 
124e5c89e4eSSatish Balay    Collective
125e5c89e4eSSatish Balay 
126e5c89e4eSSatish Balay    Level: advanced
127e5c89e4eSSatish Balay 
128e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
129e5c89e4eSSatish Balay @*/
1307087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
131e5c89e4eSSatish Balay {
132e5c89e4eSSatish Balay   PetscErrorCode ierr;
133e5c89e4eSSatish Balay   int            argc   = 0;
13402c9f0b5SLisandro Dalcin   char           **args = NULL;
135e5c89e4eSSatish Balay 
136e5c89e4eSSatish Balay   PetscFunctionBegin;
1370298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
138e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
139e5c89e4eSSatish Balay }
140e5c89e4eSSatish Balay 
141e5c89e4eSSatish Balay /*@
142e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
143e5c89e4eSSatish Balay 
14493b6d2d1SJed Brown    Level: beginner
145e5c89e4eSSatish Balay 
146e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
147e5c89e4eSSatish Balay @*/
1487087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool  *isInitialized)
149e5c89e4eSSatish Balay {
150e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
15193b6d2d1SJed Brown   return 0;
152e5c89e4eSSatish Balay }
153e5c89e4eSSatish Balay 
154e5c89e4eSSatish Balay /*@
155e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
156e5c89e4eSSatish Balay 
157e5c89e4eSSatish Balay    Level: developer
158e5c89e4eSSatish Balay 
159e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
160e5c89e4eSSatish Balay @*/
1617087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
162e5c89e4eSSatish Balay {
163e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
16493b6d2d1SJed Brown   return 0;
165e5c89e4eSSatish Balay }
166e5c89e4eSSatish Balay 
16795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void);
168e5c89e4eSSatish Balay 
169e5c89e4eSSatish Balay /*
170e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
171e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
172e5c89e4eSSatish Balay */
173367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
174e5c89e4eSSatish Balay 
175367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
176e5c89e4eSSatish Balay {
177e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
178e5c89e4eSSatish Balay 
179e5c89e4eSSatish Balay   PetscFunctionBegin;
180e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
181e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
18241e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
183e5c89e4eSSatish Balay   }
184e5c89e4eSSatish Balay 
185e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
186e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
187e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
188e5c89e4eSSatish Balay   }
189812af9f3SBarry Smith   PetscFunctionReturnVoid();
190e5c89e4eSSatish Balay }
191e5c89e4eSSatish Balay 
192e5c89e4eSSatish Balay /*
193e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
194e5c89e4eSSatish Balay sum of the second entry.
195b693b147SBarry Smith 
19676f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
197367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
198b693b147SBarry Smith there would be no place to store the both needed results.
199e5c89e4eSSatish Balay */
20076ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
201e5c89e4eSSatish Balay {
202e5c89e4eSSatish Balay   PetscErrorCode ierr;
203e5c89e4eSSatish Balay 
204e5c89e4eSSatish Balay   PetscFunctionBegin;
205d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
206d6e4c47cSJed Brown   {
207d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
208367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
209d6e4c47cSJed Brown     *max = work.max;
210d6e4c47cSJed Brown     *sum = work.sum;
211d6e4c47cSJed Brown   }
212d6e4c47cSJed Brown #else
213d6e4c47cSJed Brown   {
214d6e4c47cSJed Brown     PetscMPIInt    size,rank;
215d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
216e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
217e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
218785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
219367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2206ac3741eSJed Brown     *max = work[rank].max;
2216ac3741eSJed Brown     *sum = work[rank].sum;
222e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
223d6e4c47cSJed Brown   }
224d6e4c47cSJed Brown #endif
225e5c89e4eSSatish Balay   PetscFunctionReturn(0);
226e5c89e4eSSatish Balay }
227e5c89e4eSSatish Balay 
228e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
229e5c89e4eSSatish Balay 
230570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
23106a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
232e5c89e4eSSatish Balay 
2338cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
234e5c89e4eSSatish Balay {
235e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
236e5c89e4eSSatish Balay 
237e5c89e4eSSatish Balay   PetscFunctionBegin;
2387c2de775SJed Brown   if (*datatype == MPIU_REAL) {
239e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
240a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2417c2de775SJed Brown   }
2427c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2437c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2447c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
245a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2467c2de775SJed Brown   }
2477c2de775SJed Brown #endif
2487c2de775SJed Brown   else {
2497c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
25041e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
251e2e03761SBarry Smith   }
252812af9f3SBarry Smith   PetscFunctionReturnVoid();
253e5c89e4eSSatish Balay }
254e5c89e4eSSatish Balay #endif
255e5c89e4eSSatish Balay 
256570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
257d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
258d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
259d9822059SBarry Smith 
2608cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
261d9822059SBarry Smith {
262d9822059SBarry Smith   PetscInt i,count = *cnt;
263d9822059SBarry Smith 
264d9822059SBarry Smith   PetscFunctionBegin;
2657c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2668c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
267a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2687c2de775SJed Brown   }
2697c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2707c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2717c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2727c2de775SJed Brown     for (i=0; i<count; i++) {
2737c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2747c2de775SJed Brown     }
2757c2de775SJed Brown   }
2767c2de775SJed Brown #endif
2777c2de775SJed Brown   else {
2787c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
27941e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2808c764dc5SJose Roman   }
281d9822059SBarry Smith   PetscFunctionReturnVoid();
282d9822059SBarry Smith }
283d9822059SBarry Smith 
2848cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
285d9822059SBarry Smith {
286d9822059SBarry Smith   PetscInt    i,count = *cnt;
287d9822059SBarry Smith 
288d9822059SBarry Smith   PetscFunctionBegin;
2897c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2908c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
291a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2927c2de775SJed Brown   }
2937c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2947c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2957c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2967c2de775SJed Brown     for (i=0; i<count; i++) {
2977c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2987c2de775SJed Brown     }
2997c2de775SJed Brown   }
3007c2de775SJed Brown #endif
3017c2de775SJed Brown   else {
3028c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
30341e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
3048c764dc5SJose Roman   }
305d9822059SBarry Smith   PetscFunctionReturnVoid();
306d9822059SBarry Smith }
307d9822059SBarry Smith #endif
308d9822059SBarry Smith 
309480cf27aSJed Brown /*
310480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
311480cf27aSJed Brown 
312ff0e51ddSBarry 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.
313480cf27aSJed Brown 
31412801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
315480cf27aSJed Brown 
316480cf27aSJed Brown */
3178cc058d9SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
318480cf27aSJed Brown {
319480cf27aSJed Brown   PetscErrorCode   ierr;
320*05342407SJunchao Zhang   PetscCommCounter *counter=(PetscCommCounter*)count_val;
321480cf27aSJed Brown 
322480cf27aSJed Brown   PetscFunctionBegin;
32302c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
324*05342407SJunchao Zhang   ierr = PetscFree(counter->iflags);CHKERRMPI(ierr);
325*05342407SJunchao Zhang   ierr = PetscFree(counter);CHKERRMPI(ierr);
326480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
327480cf27aSJed Brown }
328480cf27aSJed Brown 
329480cf27aSJed Brown /*
33047435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
331da3039f7SJed Brown   calls MPI_Comm_free().
332da3039f7SJed Brown 
333da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
334480cf27aSJed Brown 
335ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
336480cf27aSJed Brown 
33712801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
338480cf27aSJed Brown 
339480cf27aSJed Brown */
340da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
341480cf27aSJed Brown {
342480cf27aSJed Brown   PetscErrorCode                    ierr;
343b89831e5SBarry Smith   PetscMPIInt                       flg;
344265f3f35SJed Brown   union {MPI_Comm comm; void *ptr;} icomm,ocomm;
345480cf27aSJed Brown 
346480cf27aSJed Brown   PetscFunctionBegin;
34712801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
348ec4fadc2SJed Brown   icomm.ptr = attr_val;
349da3039f7SJed Brown 
35047435625SJed Brown   ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
35112801b39SBarry Smith   if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm");
35212801b39SBarry Smith   if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm");
35347435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
35402c9f0b5SLisandro 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);
355da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
356b89831e5SBarry Smith }
357da3039f7SJed Brown 
358da3039f7SJed Brown /*
35947435625SJed 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.
360da3039f7SJed Brown  */
361da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
362da3039f7SJed Brown {
363da3039f7SJed Brown   PetscErrorCode ierr;
364da3039f7SJed Brown 
365da3039f7SJed Brown   PetscFunctionBegin;
36602c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
367480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
368480cf27aSJed Brown }
369480cf27aSJed Brown 
3705f7487a0SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Shm(MPI_Comm,PetscMPIInt,void *,void *);
37142218b76SBarry Smith 
372951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3738cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3748cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3758cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
376e39fd77fSBarry Smith #endif
377e39fd77fSBarry Smith 
3780f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
37912801b39SBarry Smith 
380eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
381eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3826ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
38302c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
384dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
385e5c89e4eSSatish Balay 
386dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
387051e4cf2SJed Brown {
388051e4cf2SJed Brown   PetscErrorCode ierr;
389051e4cf2SJed Brown 
390051e4cf2SJed Brown   PetscFunctionBegin;
391051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
392a1601671SBarry 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);
393051e4cf2SJed 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);
394051e4cf2SJed Brown   PetscFunctionReturn(0);
395051e4cf2SJed Brown }
396e5c89e4eSSatish Balay 
3972d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3982d747510SLisandro Dalcin 
3992d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
4002d747510SLisandro Dalcin {
4012d747510SLisandro Dalcin   PetscErrorCode ierr;
4022d747510SLisandro Dalcin 
4032d747510SLisandro Dalcin   PetscFunctionBegin;
4042d747510SLisandro Dalcin   ierr  = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
4052d747510SLisandro Dalcin   PetscFunctionReturn(0);
4062d747510SLisandro Dalcin }
4072d747510SLisandro Dalcin 
4082d747510SLisandro Dalcin /*@C
4092d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
4102d747510SLisandro Dalcin 
4112d747510SLisandro Dalcin     Not Collective
4122d747510SLisandro Dalcin 
4132d747510SLisandro Dalcin     Input Parameter:
4142d747510SLisandro Dalcin .   len - length of the string name
4152d747510SLisandro Dalcin 
4162d747510SLisandro Dalcin     Output Parameter:
4172d747510SLisandro Dalcin .   name - the name of the running program
4182d747510SLisandro Dalcin 
4192d747510SLisandro Dalcin    Level: advanced
4202d747510SLisandro Dalcin 
4212d747510SLisandro Dalcin     Notes:
4222d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4232d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4242d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4252d747510SLisandro Dalcin @*/
4262d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4272d747510SLisandro Dalcin {
4282d747510SLisandro Dalcin   PetscErrorCode ierr;
4292d747510SLisandro Dalcin 
4302d747510SLisandro Dalcin   PetscFunctionBegin;
4312d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4322d747510SLisandro Dalcin   PetscFunctionReturn(0);
4332d747510SLisandro Dalcin }
4342d747510SLisandro Dalcin 
435e5c89e4eSSatish Balay /*@C
436e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
437e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
438e5c89e4eSSatish Balay 
439e5c89e4eSSatish Balay    Not Collective
440e5c89e4eSSatish Balay 
441e5c89e4eSSatish Balay    Output Parameters:
442e5c89e4eSSatish Balay +  argc - count of number of command line arguments
443e5c89e4eSSatish Balay -  args - the command line arguments
444e5c89e4eSSatish Balay 
445e5c89e4eSSatish Balay    Level: intermediate
446e5c89e4eSSatish Balay 
447e5c89e4eSSatish Balay    Notes:
448e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
449e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
450e5c89e4eSSatish Balay 
451f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
452f177e3b1SBarry Smith 
453793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
454e5c89e4eSSatish Balay 
455e5c89e4eSSatish Balay @*/
4567087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
457e5c89e4eSSatish Balay {
458e5c89e4eSSatish Balay   PetscFunctionBegin;
45917186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
460e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
461e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
462e5c89e4eSSatish Balay   PetscFunctionReturn(0);
463e5c89e4eSSatish Balay }
464e5c89e4eSSatish Balay 
465793721a6SBarry Smith /*@C
466793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
467793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
468793721a6SBarry Smith 
469793721a6SBarry Smith    Not Collective
470793721a6SBarry Smith 
471793721a6SBarry Smith    Output Parameters:
472793721a6SBarry Smith .  args - the command line arguments
473793721a6SBarry Smith 
474793721a6SBarry Smith    Level: intermediate
475793721a6SBarry Smith 
476793721a6SBarry Smith    Notes:
477793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
478793721a6SBarry Smith 
479793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
480793721a6SBarry Smith 
481793721a6SBarry Smith @*/
4827087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
483793721a6SBarry Smith {
484793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
485793721a6SBarry Smith   PetscErrorCode ierr;
486793721a6SBarry Smith 
487793721a6SBarry Smith   PetscFunctionBegin;
48817186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4892d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
490785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
491793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
492793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
493793721a6SBarry Smith   }
4942d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
495793721a6SBarry Smith   PetscFunctionReturn(0);
496793721a6SBarry Smith }
497793721a6SBarry Smith 
498793721a6SBarry Smith /*@C
499793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
500793721a6SBarry Smith 
501793721a6SBarry Smith    Not Collective
502793721a6SBarry Smith 
503793721a6SBarry Smith    Output Parameters:
504793721a6SBarry Smith .  args - the command line arguments
505793721a6SBarry Smith 
506793721a6SBarry Smith    Level: intermediate
507793721a6SBarry Smith 
508793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
509793721a6SBarry Smith 
510793721a6SBarry Smith @*/
5117087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
512793721a6SBarry Smith {
513793721a6SBarry Smith   PetscInt       i = 0;
514793721a6SBarry Smith   PetscErrorCode ierr;
515793721a6SBarry Smith 
516793721a6SBarry Smith   PetscFunctionBegin;
517a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
518793721a6SBarry Smith   while (args[i]) {
519793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
520793721a6SBarry Smith     i++;
521793721a6SBarry Smith   }
522793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
523793721a6SBarry Smith   PetscFunctionReturn(0);
524793721a6SBarry Smith }
525793721a6SBarry Smith 
52611525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
52730befbd2SBarry Smith #include <petscconfiginfo.h>
52830befbd2SBarry Smith 
52995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
53011525c0dSBarry Smith {
53111525c0dSBarry Smith   if (!PetscGlobalRank) {
53230befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
53311525c0dSBarry Smith     int            port;
534ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
53511525c0dSBarry Smith     size_t         applinelen,introlen;
53611525c0dSBarry Smith     PetscErrorCode ierr;
537ffbd1cfbSBarry Smith     char           sawsurl[256];
53811525c0dSBarry Smith 
539c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
54011525c0dSBarry Smith     if (flg) {
54111525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
54211525c0dSBarry Smith 
543c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
54411525c0dSBarry Smith       if (sawslog[0]) {
54511525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
54611525c0dSBarry Smith       } else {
54711525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
54811525c0dSBarry Smith       }
54911525c0dSBarry Smith     }
550c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
55111525c0dSBarry Smith     if (flg) {
55211525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
55311525c0dSBarry Smith     }
554c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
555ffbd1cfbSBarry Smith     if (selectport) {
556ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
557ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
558ffbd1cfbSBarry Smith     } else {
559c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
56011525c0dSBarry Smith       if (flg) {
56111525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
56211525c0dSBarry Smith       }
563ffbd1cfbSBarry Smith     }
564c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
56511525c0dSBarry Smith     if (flg) {
56611525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
56711525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5689c1e0ce8SBarry Smith     } else {
569c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5709c1e0ce8SBarry Smith       if (flg) {
5713c01dfcfSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
5729c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5739c1e0ce8SBarry Smith       }
57411525c0dSBarry Smith     }
575c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
57611525c0dSBarry Smith     if (flg2) {
57711525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
57811525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
57911525c0dSBarry Smith       ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr);
58011525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
58111525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
58243da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
58311525c0dSBarry Smith     }
58411525c0dSBarry Smith     ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr);
58511525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
58611525c0dSBarry Smith     introlen   = 4096 + applinelen;
58730a8c9c0SSurtai Han     applinelen += 1024;
58811525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
58911525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
59011525c0dSBarry Smith 
59111525c0dSBarry Smith     if (rootlocal) {
59211525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
59311525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
59411525c0dSBarry Smith     }
59576a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
59611525c0dSBarry Smith     if (rootlocal && help) {
597928bb9adSStefano 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);
59811525c0dSBarry Smith     } else if (help) {
599928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
60011525c0dSBarry Smith     } else {
601928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
60211525c0dSBarry Smith     }
603b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
60430befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
60511525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
606a8d69d7bSBarry 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"
607df62c222SBarry 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"
608928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
60943da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
61011525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
61111525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
612ffbd1cfbSBarry Smith     if (selectport) {
613aa573868SBarry Smith       PetscBool silent;
6147d812c46SBarry Smith 
6157d812c46SBarry Smith       ierr = SAWs_Initialize();
6167d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6177d812c46SBarry Smith       while (ierr) {
6187d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6197d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6207d812c46SBarry Smith         ierr = SAWs_Initialize();
6217d812c46SBarry Smith       }
6227d812c46SBarry Smith 
623aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
624aa573868SBarry Smith       if (!silent) {
625ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
626ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
627ffbd1cfbSBarry Smith       }
6287d812c46SBarry Smith     } else {
6297d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
630aa573868SBarry Smith     }
6310af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6320af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6330af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6340af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6350af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
63611525c0dSBarry Smith   }
63711525c0dSBarry Smith   PetscFunctionReturn(0);
63811525c0dSBarry Smith }
63911525c0dSBarry Smith #endif
64011525c0dSBarry Smith 
6414dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6424dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6434dfee713SSatish Balay {
6444dfee713SSatish Balay   PetscFunctionBegin;
6454dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6464dfee713SSatish Balay     /* see MPI.py for details on this bug */
6474dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6484dfee713SSatish Balay #endif
6494dfee713SSatish Balay   PetscFunctionReturn(0);
6504dfee713SSatish Balay }
6514dfee713SSatish Balay 
652a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
653a56f64adSBarry Smith #include <adios.h>
65422580e64SBarry Smith #include <adios_read.h>
6557b56e58cSSatish Balay int64_t Petsc_adios_group;
656a56f64adSBarry Smith #endif
65755d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
65855d657eeSBarry Smith #include <adios2_c.h>
65955d657eeSBarry Smith #endif
660cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
661cd1b99a6SBarry Smith #include <omp.h>
662f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
663cd1b99a6SBarry Smith #endif
664a56f64adSBarry Smith 
665bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H)
666bc8a24ecSBarry Smith #include <dlfcn.h>
667bc8a24ecSBarry Smith #endif
668bc8a24ecSBarry Smith 
669e5c89e4eSSatish Balay /*@C
670e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
671e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
672e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
673e5c89e4eSSatish Balay    your program -- usually the very first line!
674e5c89e4eSSatish Balay 
675e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
676e5c89e4eSSatish Balay 
677e5c89e4eSSatish Balay    Input Parameters:
678e5c89e4eSSatish Balay +  argc - count of number of command line arguments
679e5c89e4eSSatish Balay .  args - the command line arguments
6800298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
681fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6820298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
683e5c89e4eSSatish Balay 
68405827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
68505827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
68605827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
68705827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
68805827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
689e5c89e4eSSatish Balay 
690e5c89e4eSSatish Balay    Options Database Keys:
6917ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6927ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
693e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
6948a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
6958a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
6968a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
6978a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
6988a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
699e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
700e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
701e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
702e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
70379dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
70479dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
70579dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
706aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
70792f119d6SBarry 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
70892f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
70992f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
71092f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
711e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
712e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
713e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
714e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
715e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7160841954dSBarry Smith -  -memory_view - Print memory usage at end of run
717e5c89e4eSSatish Balay 
718e5c89e4eSSatish Balay    Options Database Keys for Profiling:
719a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
720e94e781bSJacob Faibussowitsch +  -info [optional filename][:[~]optional list,of,classnames][:[~]optional "self"] - Prints verbose information to the screen. See PetscInfo().
721217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
722217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
723495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
724e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7259a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
72679dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
7279a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
728495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
729f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
730495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
731495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
732c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
73387c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
734c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
735495fc317SBarry Smith 
736609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
737e5c89e4eSSatish Balay 
738ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
739ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
740ffbd1cfbSBarry 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
741ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
742ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
743ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
744ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
745ffbd1cfbSBarry Smith 
746e5c89e4eSSatish Balay    Environmental Variables:
747e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
748e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
749e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
750e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
751e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
752e5c89e4eSSatish Balay 
753e5c89e4eSSatish Balay 
754e5c89e4eSSatish Balay    Level: beginner
755e5c89e4eSSatish Balay 
756e5c89e4eSSatish Balay    Notes:
757e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
758e5c89e4eSSatish Balay    it before PetscInitialize().
759e5c89e4eSSatish Balay 
760e5c89e4eSSatish Balay    Fortran Version:
761e5c89e4eSSatish Balay    In Fortran this routine has the format
762e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
763e5c89e4eSSatish Balay 
764e5c89e4eSSatish Balay +   ierr - error return code
7650eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
766fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
767e5c89e4eSSatish Balay 
768e5c89e4eSSatish Balay    Important Fortran Note:
7690eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7700298fd71SBarry Smith    null character string; you CANNOT just use NULL as
771a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
772e5c89e4eSSatish Balay 
77301cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
77401cb0274SBarry Smith    calling PetscInitialize().
775e5c89e4eSSatish Balay 
77601cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
777e5c89e4eSSatish Balay 
778e5c89e4eSSatish Balay @*/
7797087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
780e5c89e4eSSatish Balay {
781e5c89e4eSSatish Balay   PetscErrorCode ierr;
7824bb5149bSJed Brown   PetscMPIInt    flag, size;
78319c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
784e5c89e4eSSatish Balay   char           hostname[256];
785e5c89e4eSSatish Balay 
786e5c89e4eSSatish Balay   PetscFunctionBegin;
787e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7883d96e996SBarry Smith   /*
7893d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7903d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7913d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7923d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7933d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7943d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7953d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7963d96e996SBarry Smith       listed above and since that time are compatible.
7973d96e996SBarry Smith 
7983d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
7993d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
8003d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
8013d96e996SBarry Smith       to perform the checking.
8023d96e996SBarry Smith 
8033d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
8043d96e996SBarry Smith 
8053d96e996SBarry Smith       Questions:
8063d96e996SBarry Smith 
8073d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
8083d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
8093d96e996SBarry Smith   */
8103d96e996SBarry Smith 
81119c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
81219c5658aSBarry Smith   {
81319c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
81419c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
81519c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8163d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
81719c5658aSBarry Smith #if defined(MPICH_VERSION)
8183d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
81919c5658aSBarry Smith     {
82019c5658aSBarry Smith       char *ver,*lf;
82119c5658aSBarry Smith       flg = PETSC_FALSE;
82219c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
82319c5658aSBarry Smith       if (ver) {
82419c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
82519c5658aSBarry Smith         if (lf) {
82619c5658aSBarry Smith           *lf = 0;
82719c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
82819c5658aSBarry Smith         }
82919c5658aSBarry Smith       }
83019c5658aSBarry Smith       if (!flg) {
83119c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8323d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
83319c5658aSBarry Smith       }
83419c5658aSBarry Smith     }
8353d96e996SBarry Smith #endif
8363d96e996SBarry 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?) */
83719c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
83819c5658aSBarry Smith     {
83919c5658aSBarry Smith       char *ver,bs[32],*bsf;
84019c5658aSBarry Smith       flg = PETSC_FALSE;
84119c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
84219c5658aSBarry Smith       if (ver) {
8432e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
84419c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
84519c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
84619c5658aSBarry Smith       }
84719c5658aSBarry Smith       if (!flg) {
84819c5658aSBarry 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);
8493d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
85019c5658aSBarry Smith       }
85119c5658aSBarry Smith     }
85219c5658aSBarry Smith #endif
85319c5658aSBarry Smith   }
85419c5658aSBarry Smith #endif
85519c5658aSBarry Smith 
856bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM)
857bc8a24ecSBarry Smith   {
858bc8a24ecSBarry Smith     PetscInt cnt = 0;
859bc8a24ecSBarry 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 */
860bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++;
861bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++;
862bc8a24ecSBarry Smith     if (cnt > 1) {
863bc8a24ecSBarry Smith       fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
864bc8a24ecSBarry Smith       return PETSC_ERR_MPI_LIB_INCOMP;
865bc8a24ecSBarry Smith     }
866bc8a24ecSBarry Smith   }
867bc8a24ecSBarry Smith #endif
868bc8a24ecSBarry Smith 
869ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
870d89683f4Sbcordonn   PETSC_STDOUT = stdout;
871ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
872e5c89e4eSSatish Balay 
8730c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8740c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8750c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8760c30907bSSatish Balay #endif
8770c30907bSSatish Balay 
8784416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
879e5c89e4eSSatish Balay 
880e5c89e4eSSatish Balay   /*
881e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
882e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
883e5c89e4eSSatish Balay   */
884e5c89e4eSSatish Balay   if (argc && *argc) {
885e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
886e5c89e4eSSatish Balay   } else {
887e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
888e5c89e4eSSatish Balay   }
889e5c89e4eSSatish Balay 
890e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
891e5c89e4eSSatish Balay   if (!flag) {
892e32f2f54SBarry 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");
8934dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
8945e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8955e765c61SJed Brown     {
8965e765c61SJed Brown       PetscMPIInt provided;
8975e765c61SJed Brown       ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr);
8985e765c61SJed Brown     }
8995e765c61SJed Brown #else
900e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
9015e765c61SJed Brown #endif
902e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
903e5c89e4eSSatish Balay   }
9044dfee713SSatish Balay 
905e5c89e4eSSatish Balay   if (argc && args) {
906e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
907e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
908e5c89e4eSSatish Balay   }
909e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
9105ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
9115ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
9125ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
913ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
914e5c89e4eSSatish Balay 
915a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
916d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
917e5c89e4eSSatish Balay 
9180f9be574SPeter Hill   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
91912801b39SBarry Smith     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
92012801b39SBarry Smith     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
9210f9be574SPeter Hill   }
92212801b39SBarry Smith 
923e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
924e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
925e5c89e4eSSatish Balay 
926e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
927e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
928e5c89e4eSSatish Balay 
9298ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9308ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9317cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
932e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
933e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
934e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
935e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
936e316c87fSJed Brown #endif
937e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9388ad47952SJed Brown 
939e5c89e4eSSatish Balay   /*
940e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
941e5c89e4eSSatish Balay      shared libraries the constructors for global variables
942e5c89e4eSSatish Balay      are not called; at least on IRIX.
943e5c89e4eSSatish Balay   */
944886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
945e5c89e4eSSatish Balay   {
946252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
94750f81f78SJed Brown     PetscComplex ic(0.0,1.0);
948e5c89e4eSSatish Balay     PETSC_i = ic;
949252ecd31SSatish Balay #else
95050f81f78SJed Brown     PETSC_i = _Complex_I;
951b7940d39SSatish Balay #endif
952762437b8SSatish Balay   }
953762437b8SSatish Balay 
9542c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
955e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
956500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
957500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
958500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9592c876bd9SBarry Smith #endif
960886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
961e5c89e4eSSatish Balay 
962e5c89e4eSSatish Balay   /*
963e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
964e5c89e4eSSatish Balay      half of the entries and maxes the second half.
965e5c89e4eSSatish Balay   */
966367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
967e5c89e4eSSatish Balay 
968ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
969c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
970c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9717c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9728c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9738c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9748c764dc5SJose Roman #endif
975d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
976d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
977570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
978570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
979570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
980570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
981570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
982c90a1750SBarry Smith #endif
983c90a1750SBarry Smith 
984570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
985cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
986cca4cb22SSatish Balay #endif
987cca4cb22SSatish Balay 
988e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
989e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
990e5c89e4eSSatish Balay 
99140df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
992e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
993e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
99444041f26SJed Brown #endif
995e5c89e4eSSatish Balay 
996ec957eceSBarry Smith 
997e5c89e4eSSatish Balay   /*
998480cf27aSJed Brown      Attributes to be set on PETSc communicators
999480cf27aSJed Brown   */
100012801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
100112801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
100212801b39SBarry Smith   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
10035f7487a0SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Shm,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
1004480cf27aSJed Brown 
1005480cf27aSJed Brown   /*
1006e8fb0fc0SBarry Smith      Build the options database
1007e5c89e4eSSatish Balay   */
1008c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
1009e5c89e4eSSatish Balay 
1010703f3eceSBarry Smith   /* call a second time so it can look in the options database */
1011703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
10126dc8fec2Sbcordonn 
1013e5c89e4eSSatish Balay   /*
1014e5c89e4eSSatish Balay      Print main application help message
1015e5c89e4eSSatish Balay   */
10162d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr);
1017e5c89e4eSSatish Balay   if (help && flg) {
1018e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr);
1019e5c89e4eSSatish Balay   }
1020e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr);
1021e5c89e4eSSatish Balay 
1022d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
1023d45a07a7SBarry Smith 
1024e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
102511525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1026f4202a44SBarry Smith #endif
1027f4202a44SBarry Smith 
1028e5c89e4eSSatish Balay   /*
1029e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1030e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1031e5c89e4eSSatish Balay   */
1032e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1033e5c89e4eSSatish Balay 
1034e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
103502c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1036e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
103702c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1038cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1039101491d0SBarry Smith   {
1040101491d0SBarry Smith     PetscBool omp_view_flag;
1041cd1b99a6SBarry Smith     char      *threads = getenv("OMP_NUM_THREADS");
1042e5c89e4eSSatish Balay 
1043cd1b99a6SBarry Smith    if (threads) {
104402c9f0b5SLisandro Dalcin      ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
1045f90b075cSBarry Smith      (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
1046cd1b99a6SBarry Smith    } else {
1047cd1b99a6SBarry Smith #define NMAX  10000
1048101491d0SBarry Smith      int          i;
1049cd1b99a6SBarry Smith       PetscScalar *x;
1050cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1051cd1b99a6SBarry Smith #pragma omp parallel for
1052cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1053cd1b99a6SBarry Smith         x[i] = 0.0;
1054f90b075cSBarry Smith         PetscNumOMPThreads  = (PetscInt) omp_get_num_threads();
1055cd1b99a6SBarry Smith       }
1056cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
105702c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1058101491d0SBarry Smith     }
1059101491d0SBarry Smith     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
1060f90b075cSBarry 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);
1061101491d0SBarry Smith     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
1062101491d0SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1063101491d0SBarry Smith     if (flg) {
106402c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1065f90b075cSBarry Smith       omp_set_num_threads((int)PetscNumOMPThreads);
1066101491d0SBarry Smith     }
1067101491d0SBarry Smith     if (omp_view_flag) {
1068f90b075cSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
1069cd1b99a6SBarry Smith     }
1070cd1b99a6SBarry Smith   }
1071cd1b99a6SBarry Smith #endif
1072e5c89e4eSSatish Balay   ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr);
1073ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1074c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1075ef6c6fedSBoyana Norris 
1076951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1077e39fd77fSBarry Smith   /*
1078e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1079e39fd77fSBarry Smith 
1080e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1081e39fd77fSBarry Smith   */
108230815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
10830298fd71SBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
108430815ce0SLisandro Dalcin   }
1085951e3c8eSBarry Smith #endif
1086e39fd77fSBarry Smith 
108741c0b4b3SShri Abhyankar   /*
108841c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
108941c0b4b3SShri Abhyankar   */
1090ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1091e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1092e1167bb9SShri Abhyankar #endif
1093e1167bb9SShri Abhyankar 
10942d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10952d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10962d53ad75SBarry Smith #endif
10972d53ad75SBarry Smith 
10985e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
109987c3beb6SLisandro Dalcin   {
110087c3beb6SLisandro Dalcin     PetscViewer viewer;
110122e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
11025e71baefSBarry Smith     if (flg) {
11035e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
110487c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
110587c3beb6SLisandro Dalcin     }
11065e71baefSBarry Smith   }
11075e71baefSBarry Smith #endif
1108dff31646SBarry Smith 
110987c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
111087c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
111187c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
111287c3beb6SLisandro Dalcin 
1113a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1114a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
11157b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1116a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
11179fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1118a56f64adSBarry Smith #endif
111955d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
112055d657eeSBarry Smith #endif
1121a56f64adSBarry Smith 
1122301d30feSBarry Smith   /*
1123301d30feSBarry Smith       Once we are completedly initialized then we can set this variables
1124301d30feSBarry Smith   */
1125301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
11262db0e300SLisandro Dalcin 
11272db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
11282db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1129301d30feSBarry Smith   PetscFunctionReturn(0);
1130e5c89e4eSSatish Balay }
1131e5c89e4eSSatish Balay 
11324097062eSBarry Smith #if defined(PETSC_USE_LOG)
113395c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1134ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1135ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
113695c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11374097062eSBarry Smith #endif
1138e5c89e4eSSatish Balay 
1139008a6e76SBarry Smith /*
1140008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1141008a6e76SBarry Smith */
1142008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1143008a6e76SBarry Smith {
1144008a6e76SBarry Smith   PetscErrorCode ierr;
1145008a6e76SBarry Smith 
1146008a6e76SBarry Smith   PetscFunctionBegin;
1147008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1148008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1149008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1150008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1151008a6e76SBarry Smith #endif
1152008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1153008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1154008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1155008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1156008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1157008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1158008a6e76SBarry Smith #endif
1159008a6e76SBarry Smith 
1160008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1161008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1162008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1163008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1164008a6e76SBarry Smith #endif
1165008a6e76SBarry Smith #endif
1166008a6e76SBarry Smith 
1167008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1168008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1169008a6e76SBarry Smith #endif
1170008a6e76SBarry Smith 
1171008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
117240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1173008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1174008a6e76SBarry Smith #endif
1175008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1176008a6e76SBarry Smith   PetscFunctionReturn(0);
1177008a6e76SBarry Smith }
1178008a6e76SBarry Smith 
1179e5c89e4eSSatish Balay /*@C
1180e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1181e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1182e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1183e5c89e4eSSatish Balay 
1184e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1185e5c89e4eSSatish Balay 
1186e5c89e4eSSatish Balay    Options Database Keys:
118726a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1188e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11897eb1d149SBarry 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
1190e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
119192f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1192e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
119392f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1194e5c89e4eSSatish Balay 
1195e5c89e4eSSatish Balay    Level: beginner
1196e5c89e4eSSatish Balay 
1197e5c89e4eSSatish Balay    Note:
1198e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1199e5c89e4eSSatish Balay 
120088c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1201e5c89e4eSSatish Balay @*/
12027087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1203e5c89e4eSSatish Balay {
1204e5c89e4eSSatish Balay   PetscErrorCode ierr;
12054bb5149bSJed Brown   PetscMPIInt    rank;
1206a8d2bbe5SBarry Smith   PetscInt       nopt;
12072bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1208dff31646SBarry Smith   PetscBool      flg;
120910463e74SBarry Smith #if defined(PETSC_USE_LOG)
121010463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
121110463e74SBarry Smith #endif
1212e5c89e4eSSatish Balay 
1213e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
12144b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
12153cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1216e5c89e4eSSatish Balay   }
12173cbbc5ffSBarry Smith   PetscFunctionBegin;
12180298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1219b022a5c1SBarry Smith 
12201f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1221a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
122222580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1223a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1224a56f64adSBarry Smith #endif
122555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
122655d657eeSBarry Smith #endif
1227c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1228dff31646SBarry Smith   if (flg) {
12291f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
12301f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
12311f817a21SBarry Smith 
1232c5929fdfSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
12331f817a21SBarry Smith     if (filename[0]) {
12341f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
12351f817a21SBarry Smith     }
1236dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1237dff31646SBarry Smith     cits[0] = 0;
1238dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
12391f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
12401f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12411f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
12421f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12431f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1244dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1245dff31646SBarry Smith   }
1246dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1247dff31646SBarry Smith 
1248c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
124904102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
125004102261SBarry Smith   {
125104102261SBarry Smith     PetscInt nmax = 2;
125204102261SBarry Smith     char     **buffs;
125304102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1254c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
125504102261SBarry Smith     if (flg1) {
125604102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
125704102261SBarry Smith       if (nmax == 1) {
125804102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
125904102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
126004102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
126104102261SBarry Smith       }
126204102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
126304102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
126404102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
126504102261SBarry Smith     }
126604102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
126704102261SBarry Smith   }
1268763ec7b1SBarry Smith   {
1269763ec7b1SBarry Smith     PetscInt nmax = 2;
1270763ec7b1SBarry Smith     char     **buffs;
1271763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1272763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1273763ec7b1SBarry Smith     if (flg1) {
1274763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1275763ec7b1SBarry Smith       if (nmax == 1) {
1276763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1277763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1278763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1279763ec7b1SBarry Smith       }
1280763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1281763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1282763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1283763ec7b1SBarry Smith     }
1284763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1285763ec7b1SBarry Smith   }
128604102261SBarry Smith #endif
128767234432SDmitry Karpeev   /*
128867234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
128967234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
129067234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
129167234432SDmitry Karpeev   */
129267234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
129304102261SBarry Smith 
12942d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12952d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12962d53ad75SBarry Smith #endif
12972d53ad75SBarry Smith 
12982d53ad75SBarry Smith 
1299e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1300dff31646SBarry Smith   flg = PETSC_FALSE;
1301c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1302d5649816SBarry Smith   if (flg) {
1303e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1304d5649816SBarry Smith   }
1305d5649816SBarry Smith #endif
1306d5649816SBarry Smith 
1307681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1308681455b2SBarry Smith   flg1 = PETSC_FALSE;
1309c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1310681455b2SBarry Smith   if (flg1) {
1311681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1312681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1313681455b2SBarry Smith   }
1314681455b2SBarry Smith #endif
1315681455b2SBarry Smith 
131667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1317c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1318e5c89e4eSSatish Balay   if (!flg2) {
131990d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1320c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1321e5c89e4eSSatish Balay   }
1322e5c89e4eSSatish Balay   if (flg2) {
13230841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1324e5c89e4eSSatish Balay   }
132567584ceeSBarry Smith #endif
1326e5c89e4eSSatish Balay 
1327e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
132890d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1329c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1330e5c89e4eSSatish Balay   if (flg1) {
1331e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1332205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1333e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1334e5c89e4eSSatish Balay   }
1335e5c89e4eSSatish Balay #endif
1336e5c89e4eSSatish Balay 
1337e5c89e4eSSatish Balay 
1338e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1339e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1340e5c89e4eSSatish Balay   mname[0] = 0;
1341c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1342e5c89e4eSSatish Balay   if (flg1) {
1343e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1344e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1345e5c89e4eSSatish Balay   }
1346e5c89e4eSSatish Balay #endif
1347356e5837SBarry Smith #endif
1348a297a907SKarl Rupp 
1349dd710f27SBarry Smith   /*
1350dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1351dd710f27SBarry Smith   */
1352dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1353dd710f27SBarry Smith 
1354356e5837SBarry Smith #if defined(PETSC_USE_LOG)
135587c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1356f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
135787c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
135887c3beb6SLisandro Dalcin 
1359356e5837SBarry Smith   mname[0] = 0;
1360c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1361e5c89e4eSSatish Balay   if (flg1) {
136291eabc43SBarry Smith     PetscViewer viewer;
136320a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
136491eabc43SBarry Smith     if (mname[0]) {
136591eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
136691eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13676bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
136833f85c2fSBarry Smith     } else {
136933f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13709a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
137133f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13729a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
137333f85c2fSBarry Smith     }
1374e5c89e4eSSatish Balay   }
1375a297a907SKarl Rupp 
1376dd710f27SBarry Smith   /*
1377dd710f27SBarry Smith      Free any objects created by the last block of code.
1378dd710f27SBarry Smith   */
1379dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1380dd710f27SBarry Smith 
1381dd710f27SBarry Smith   mname[0] = 0;
1382c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr);
1383c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr);
13847ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1385e5c89e4eSSatish Balay #endif
138610463e74SBarry Smith 
138733f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
138810463e74SBarry Smith 
138990d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1390c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1391e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
139290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1393c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1394e5c89e4eSSatish Balay   if (flg1) {
1395e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1396e5c89e4eSSatish Balay   }
139790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
139890d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13998bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1400c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1401c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1402c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1403e4c476e2SSatish Balay 
1404e5c89e4eSSatish Balay   if (flg2) {
1405be56827dSJed Brown     PetscViewer viewer;
140602ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
140702ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1408c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1409be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1410e5c89e4eSSatish Balay   }
1411e5c89e4eSSatish Balay 
1412e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1413c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1414c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1415e5c89e4eSSatish Balay 
141633fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1417c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
14183de2bfdfSBarry Smith #if defined(PETSC_USE_DEBUG)
14193de2bfdfSBarry Smith   if (!flg1) flg3 = PETSC_TRUE;
14203de2bfdfSBarry Smith #endif
1421e5c89e4eSSatish Balay   if (flg3) {
14223de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1423be56827dSJed Brown       PetscViewer viewer;
142402ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
142502ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1426c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1427be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1428e5c89e4eSSatish Balay     }
14293de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
14303de2bfdfSBarry Smith     if (nopt) {
14313de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
14323de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
14333de2bfdfSBarry Smith       if (nopt == 1) {
1434e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1435e5c89e4eSSatish Balay       } else {
14367582186dSLisandro Dalcin         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1437e5c89e4eSSatish Balay       }
14383de2bfdfSBarry Smith     } else if (flg3 && flg1) {
14393de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1440df12ba86SBarry Smith     }
1441c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1442e5c89e4eSSatish Balay   }
1443e5c89e4eSSatish Balay 
1444e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1445d45a07a7SBarry Smith   if (!PetscGlobalRank) {
144687f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
144716ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1448d45a07a7SBarry Smith   }
1449ec957eceSBarry Smith #endif
1450ec957eceSBarry Smith 
14514097062eSBarry Smith #if defined(PETSC_USE_LOG)
145210463e74SBarry Smith   /*
1453dbc8283eSBarry Smith        List all objects the user may have forgot to free
14542eff7a51SBarry Smith   */
145505df10baSBarry Smith   if (PetscObjectsLog) {
1456c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1457a64a8e02SBarry Smith     if (flg1) {
1458a64a8e02SBarry Smith       MPI_Comm local_comm;
14597eb1d149SBarry Smith       char     string[64];
1460a64a8e02SBarry Smith 
1461c5929fdfSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
1462a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1463a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14647eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1465a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1466a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14670a1571b3SBarry Smith     }
146805df10baSBarry Smith   }
14694097062eSBarry Smith #endif
14704097062eSBarry Smith 
14714097062eSBarry Smith #if defined(PETSC_USE_LOG)
1472dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1473dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1474a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14754097062eSBarry Smith #endif
14762eff7a51SBarry Smith 
147733f85c2fSBarry Smith   /*
147833f85c2fSBarry Smith      Destroy any packages that registered a finalize
147933f85c2fSBarry Smith   */
148033f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
148133f85c2fSBarry Smith 
1482101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1483fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1484101409b8SToby Isaac #endif
1485101409b8SToby Isaac 
148633f85c2fSBarry Smith   /*
148748dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
148848dd1dffSBarry Smith 
148937e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
149048dd1dffSBarry Smith   */
149137e93019SBarry Smith 
14924028d114SSatish Balay   if (petsc_history) {
1493f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
149402c9f0b5SLisandro Dalcin     petsc_history = NULL;
1495e5c89e4eSSatish Balay   }
14969de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1497e94e781bSJacob Faibussowitsch   ierr = PetscInfoDestroy();CHKERRQ(ierr);
1498e5c89e4eSSatish Balay 
149967584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
150092f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1501e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
150292f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1503e5c89e4eSSatish Balay     FILE *fd;
1504ed9cf6e9SBarry Smith     int  err;
1505e5c89e4eSSatish Balay 
1506dc92acbaSJed Brown     flg2 = PETSC_FALSE;
150792f119d6SBarry Smith     flg3 = PETSC_FALSE;
15088bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG)
150992f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);
1510dc92acbaSJed Brown #endif
151192f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
151292f119d6SBarry Smith     fname[0] = 0;
151392f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr);
1514e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1515e5c89e4eSSatish Balay 
15162e924ca5SSatish Balay       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
1517e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1518e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1519ed9cf6e9SBarry Smith       err  = fclose(fd);
1520e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
152192f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1522e5c89e4eSSatish Balay       MPI_Comm local_comm;
1523e5c89e4eSSatish Balay 
1524e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1525e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1526e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1527e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1528e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1529e5c89e4eSSatish Balay     }
1530e5c89e4eSSatish Balay     fname[0] = 0;
153192f119d6SBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,250,&flg1);CHKERRQ(ierr);
1532e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1533e5c89e4eSSatish Balay 
153492f119d6SBarry Smith       PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank);
153592f119d6SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
153692f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1537ed9cf6e9SBarry Smith       err  = fclose(fd);
1538e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
153992f119d6SBarry Smith     } else if (flg1) {
154092f119d6SBarry Smith       MPI_Comm local_comm;
154192f119d6SBarry Smith 
154292f119d6SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
154392f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
154492f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
154592f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
154692f119d6SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1547e5c89e4eSSatish Balay     }
1548e5c89e4eSSatish Balay   }
154967584ceeSBarry Smith #endif
155020e2c332SMatthew G. Knepley 
15515486ca60SMatthew G. Knepley   /*
15525486ca60SMatthew G. Knepley      Close any open dynamic libraries
15535486ca60SMatthew G. Knepley   */
15545486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15555486ca60SMatthew G. Knepley 
1556e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15574416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1558e5c89e4eSSatish Balay 
1559e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
156002c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1561e5c89e4eSSatish Balay 
1562008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1563e5c89e4eSSatish Balay 
1564dbc8283eSBarry Smith   /*
1565efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1566efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1567efb80d3cSBarry Smith 
1568efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1569efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1570dbc8283eSBarry Smith  */
1571b770b1f6SSatish Balay   {
1572dbc8283eSBarry Smith     PetscCommCounter *counter;
1573dbc8283eSBarry Smith     PetscMPIInt      flg;
1574dbc8283eSBarry Smith     MPI_Comm         icomm;
1575265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
157647435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1577dbc8283eSBarry Smith     if (flg) {
1578265f3f35SJed Brown       icomm = ucomm.comm;
157947435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1580dbc8283eSBarry 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");
1581dbc8283eSBarry Smith 
158247435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
158347435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1584efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1585dbc8283eSBarry Smith     }
158647435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1587dbc8283eSBarry Smith     if (flg) {
1588265f3f35SJed Brown       icomm = ucomm.comm;
158947435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1590dbc8283eSBarry 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");
1591dbc8283eSBarry Smith 
159247435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
159347435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1594efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1595dbc8283eSBarry Smith     }
1596b770b1f6SSatish Balay   }
1597dbc8283eSBarry Smith 
159847435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
159947435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
160047435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
16015f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1602480cf27aSJed Brown 
16035ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
16045ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
16055ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1606ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1607ef19f930SBarry Smith 
1608e5c89e4eSSatish Balay   if (PetscBeganMPI) {
160999608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
161099b1327fSBarry Smith     PetscMPIInt flag;
161199b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1612e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
161399608316SBarry Smith #endif
1614e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1615e5c89e4eSSatish Balay   }
1616e5c89e4eSSatish Balay /*
1617e5c89e4eSSatish Balay 
1618e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1619e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1620e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1621e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1622e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16230e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1624e5c89e4eSSatish Balay    memory was not freed.
1625e5c89e4eSSatish Balay 
1626e5c89e4eSSatish Balay */
16271d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1628a297a907SKarl Rupp 
1629e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1630e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
16313db9a53dSBarry Smith   PetscFunctionReturn(0);
1632e5c89e4eSSatish Balay }
1633e5c89e4eSSatish Balay 
163443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16358cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
163643db4dbbSBarry Smith {
163743db4dbbSBarry Smith   if (*a == *b) return 1;
163843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
163943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
164043db4dbbSBarry Smith   return 0;
164143db4dbbSBarry Smith }
1642a70650f6SBarry Smith #endif
164343db4dbbSBarry Smith 
164443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16458cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
164643db4dbbSBarry Smith {
164743db4dbbSBarry Smith   if (*a == *b) return 1;
164843db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
164943db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
165043db4dbbSBarry Smith   return 0;
165143db4dbbSBarry Smith }
165243db4dbbSBarry Smith #endif
1653