xref: /petsc/src/sys/objects/pinit.c (revision 8ad20175d88643dafb91dc893b19fc7d09cd0ff4)
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 
336de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */
34e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
356de5d289SStefano Zampini #if defined(PETSC_HAVE_MPI_INIT_THREAD)
366de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
376de5d289SStefano Zampini #else
386de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
396de5d289SStefano Zampini #endif
40e5c89e4eSSatish Balay 
41480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
42480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
43480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
445f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;
45480cf27aSJed Brown 
46e5c89e4eSSatish Balay /*
47e5c89e4eSSatish Balay      Declare and set all the string names of the PETSc enums
48e5c89e4eSSatish Balay */
4902c9f0b5SLisandro Dalcin const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
5002c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};
51e5c89e4eSSatish Balay 
52ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE;
53ace3abfcSBarry Smith PetscBool PetscPreLoadingOn   = PETSC_FALSE;
540f8e0872SSatish Balay 
55a2f94806SJed Brown PetscInt PetscHotRegionDepth;
56a2f94806SJed Brown 
57b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY)
58b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen;
59b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout;
60b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr;
61b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock;
62b22622e2STadeu Manoel #endif
63b22622e2STadeu Manoel 
64e5c89e4eSSatish Balay /*
65945d1669SBarry Smith       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args
6672a42c3cSBarry Smith 
6772a42c3cSBarry Smith    Collective
6872a42c3cSBarry Smith 
6972a42c3cSBarry Smith    Level: advanced
7072a42c3cSBarry Smith 
7195452b02SPatrick Sanan     Notes:
72a56f64adSBarry Smith     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
730f11a792SBarry Smith      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
74a56f64adSBarry Smith      be called multiple times from Julia without the problem of trying to initialize MPI more than once.
750f11a792SBarry Smith 
76a56f64adSBarry Smith      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals
771ea5a559SBarry Smith 
7872a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
790f11a792SBarry Smith */
80945d1669SBarry Smith PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
8172a42c3cSBarry Smith {
8272a42c3cSBarry Smith   PetscErrorCode ierr;
8372a42c3cSBarry Smith   int            myargc   = argc;
8472a42c3cSBarry Smith   char           **myargs = args;
8572a42c3cSBarry Smith 
8672a42c3cSBarry Smith   PetscFunctionBegin;
87c52ac268SRichard Tran Mills   ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr;
881ea5a559SBarry Smith   ierr = PetscPopSignalHandler();CHKERRQ(ierr);
89df413903SBarry Smith   PetscBeganMPI = PETSC_FALSE;
9072a42c3cSBarry Smith   PetscFunctionReturn(ierr);
9172a42c3cSBarry Smith }
9272a42c3cSBarry Smith 
93f0865b08SBarry Smith /*
94a56f64adSBarry Smith       Used by Julia interface to get communicator
95f0865b08SBarry Smith */
96945d1669SBarry Smith PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
97f0865b08SBarry Smith {
98f0865b08SBarry Smith   PetscFunctionBegin;
99f0865b08SBarry Smith   *comm = PETSC_COMM_SELF;
100f0865b08SBarry Smith   PetscFunctionReturn(0);
101f0865b08SBarry Smith }
102f0865b08SBarry Smith 
103e5c89e4eSSatish Balay /*@C
104e5c89e4eSSatish Balay       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
105e5c89e4eSSatish Balay         the command line arguments.
106e5c89e4eSSatish Balay 
107e5c89e4eSSatish Balay    Collective
108e5c89e4eSSatish Balay 
109e5c89e4eSSatish Balay    Level: advanced
110e5c89e4eSSatish Balay 
111e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran()
112e5c89e4eSSatish Balay @*/
1137087cfbeSBarry Smith PetscErrorCode  PetscInitializeNoArguments(void)
114e5c89e4eSSatish Balay {
115e5c89e4eSSatish Balay   PetscErrorCode ierr;
116e5c89e4eSSatish Balay   int            argc   = 0;
11702c9f0b5SLisandro Dalcin   char           **args = NULL;
118e5c89e4eSSatish Balay 
119e5c89e4eSSatish Balay   PetscFunctionBegin;
1200298fd71SBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);
121e5c89e4eSSatish Balay   PetscFunctionReturn(ierr);
122e5c89e4eSSatish Balay }
123e5c89e4eSSatish Balay 
124e5c89e4eSSatish Balay /*@
125e5c89e4eSSatish Balay       PetscInitialized - Determine whether PETSc is initialized.
126e5c89e4eSSatish Balay 
12793b6d2d1SJed Brown    Level: beginner
128e5c89e4eSSatish Balay 
129e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
130e5c89e4eSSatish Balay @*/
1317087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized)
132e5c89e4eSSatish Balay {
133e5c89e4eSSatish Balay   *isInitialized = PetscInitializeCalled;
13493b6d2d1SJed Brown   return 0;
135e5c89e4eSSatish Balay }
136e5c89e4eSSatish Balay 
137e5c89e4eSSatish Balay /*@
138e5c89e4eSSatish Balay       PetscFinalized - Determine whether PetscFinalize() has been called yet
139e5c89e4eSSatish Balay 
140e5c89e4eSSatish Balay    Level: developer
141e5c89e4eSSatish Balay 
142e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
143e5c89e4eSSatish Balay @*/
1447087cfbeSBarry Smith PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
145e5c89e4eSSatish Balay {
146e5c89e4eSSatish Balay   *isFinalized = PetscFinalizeCalled;
14793b6d2d1SJed Brown   return 0;
148e5c89e4eSSatish Balay }
149e5c89e4eSSatish Balay 
15057171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []);
151e5c89e4eSSatish Balay 
152e5c89e4eSSatish Balay /*
153e5c89e4eSSatish Balay        This function is the MPI reduction operation used to compute the sum of the
154e5c89e4eSSatish Balay    first half of the datatype and the max of the second half.
155e5c89e4eSSatish Balay */
156367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0;
157e5c89e4eSSatish Balay 
158367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
159e5c89e4eSSatish Balay {
160e5c89e4eSSatish Balay   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;
161e5c89e4eSSatish Balay 
162e5c89e4eSSatish Balay   PetscFunctionBegin;
163e5c89e4eSSatish Balay   if (*datatype != MPIU_2INT) {
164e5c89e4eSSatish Balay     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
16541e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
166e5c89e4eSSatish Balay   }
167e5c89e4eSSatish Balay 
168e5c89e4eSSatish Balay   for (i=0; i<count; i++) {
169e5c89e4eSSatish Balay     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
170e5c89e4eSSatish Balay     xout[2*i+1] += xin[2*i+1];
171e5c89e4eSSatish Balay   }
172812af9f3SBarry Smith   PetscFunctionReturnVoid();
173e5c89e4eSSatish Balay }
174e5c89e4eSSatish Balay 
175e5c89e4eSSatish Balay /*
176e5c89e4eSSatish Balay     Returns the max of the first entry owned by this processor and the
177e5c89e4eSSatish Balay sum of the second entry.
178b693b147SBarry Smith 
17976f543a4SJed Brown     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
180367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
181b693b147SBarry Smith there would be no place to store the both needed results.
182e5c89e4eSSatish Balay */
18376ec1555SBarry Smith PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
184e5c89e4eSSatish Balay {
185e5c89e4eSSatish Balay   PetscErrorCode ierr;
186e5c89e4eSSatish Balay 
187e5c89e4eSSatish Balay   PetscFunctionBegin;
188d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
189d6e4c47cSJed Brown   {
190d6e4c47cSJed Brown     struct {PetscInt max,sum;} work;
191367daffbSBarry Smith     ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
192d6e4c47cSJed Brown     *max = work.max;
193d6e4c47cSJed Brown     *sum = work.sum;
194d6e4c47cSJed Brown   }
195d6e4c47cSJed Brown #else
196d6e4c47cSJed Brown   {
197d6e4c47cSJed Brown     PetscMPIInt    size,rank;
198d6e4c47cSJed Brown     struct {PetscInt max,sum;} *work;
199e5c89e4eSSatish Balay     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
200e5c89e4eSSatish Balay     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
201785e854fSJed Brown     ierr = PetscMalloc1(size,&work);CHKERRQ(ierr);
202367daffbSBarry Smith     ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr);
2036ac3741eSJed Brown     *max = work[rank].max;
2046ac3741eSJed Brown     *sum = work[rank].sum;
205e5c89e4eSSatish Balay     ierr = PetscFree(work);CHKERRQ(ierr);
206d6e4c47cSJed Brown   }
207d6e4c47cSJed Brown #endif
208e5c89e4eSSatish Balay   PetscFunctionReturn(0);
209e5c89e4eSSatish Balay }
210e5c89e4eSSatish Balay 
211e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/
212e5c89e4eSSatish Balay 
213570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
21406a205a8SBarry Smith MPI_Op MPIU_SUM = 0;
215e5c89e4eSSatish Balay 
2168cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
217e5c89e4eSSatish Balay {
218e5c89e4eSSatish Balay   PetscInt i,count = *cnt;
219e5c89e4eSSatish Balay 
220e5c89e4eSSatish Balay   PetscFunctionBegin;
2217c2de775SJed Brown   if (*datatype == MPIU_REAL) {
222e2e03761SBarry Smith     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
223a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2247c2de775SJed Brown   }
2257c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2267c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2277c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
228a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] += xin[i];
2297c2de775SJed Brown   }
2307c2de775SJed Brown #endif
2317c2de775SJed Brown   else {
2327c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
23341e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
234e2e03761SBarry Smith   }
235812af9f3SBarry Smith   PetscFunctionReturnVoid();
236e5c89e4eSSatish Balay }
237e5c89e4eSSatish Balay #endif
238e5c89e4eSSatish Balay 
239570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
240d9822059SBarry Smith MPI_Op MPIU_MAX = 0;
241d9822059SBarry Smith MPI_Op MPIU_MIN = 0;
242d9822059SBarry Smith 
2438cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
244d9822059SBarry Smith {
245d9822059SBarry Smith   PetscInt i,count = *cnt;
246d9822059SBarry Smith 
247d9822059SBarry Smith   PetscFunctionBegin;
2487c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2498c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
250a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
2517c2de775SJed Brown   }
2527c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2537c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2547c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2557c2de775SJed Brown     for (i=0; i<count; i++) {
2567c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2577c2de775SJed Brown     }
2587c2de775SJed Brown   }
2597c2de775SJed Brown #endif
2607c2de775SJed Brown   else {
2617c2de775SJed Brown     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
26241e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2638c764dc5SJose Roman   }
264d9822059SBarry Smith   PetscFunctionReturnVoid();
265d9822059SBarry Smith }
266d9822059SBarry Smith 
2678cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
268d9822059SBarry Smith {
269d9822059SBarry Smith   PetscInt    i,count = *cnt;
270d9822059SBarry Smith 
271d9822059SBarry Smith   PetscFunctionBegin;
2727c2de775SJed Brown   if (*datatype == MPIU_REAL) {
2738c764dc5SJose Roman     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
274a297a907SKarl Rupp     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
2757c2de775SJed Brown   }
2767c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
2777c2de775SJed Brown   else if (*datatype == MPIU_COMPLEX) {
2787c2de775SJed Brown     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
2797c2de775SJed Brown     for (i=0; i<count; i++) {
2807c2de775SJed Brown       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
2817c2de775SJed Brown     }
2827c2de775SJed Brown   }
2837c2de775SJed Brown #endif
2847c2de775SJed Brown   else {
2858c764dc5SJose Roman     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
28641e02c4dSJunchao Zhang     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
2878c764dc5SJose Roman   }
288d9822059SBarry Smith   PetscFunctionReturnVoid();
289d9822059SBarry Smith }
290d9822059SBarry Smith #endif
291d9822059SBarry Smith 
292480cf27aSJed Brown /*
293480cf27aSJed Brown    Private routine to delete internal tag/name counter storage when a communicator is freed.
294480cf27aSJed Brown 
295ff0e51ddSBarry 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.
296480cf27aSJed Brown 
29712801b39SBarry Smith    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
298480cf27aSJed Brown 
299480cf27aSJed Brown */
30033779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
301480cf27aSJed Brown {
302480cf27aSJed Brown   PetscErrorCode   ierr;
30305342407SJunchao Zhang   PetscCommCounter *counter=(PetscCommCounter*)count_val;
304480cf27aSJed Brown 
305480cf27aSJed Brown   PetscFunctionBegin;
30602c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
30705342407SJunchao Zhang   ierr = PetscFree(counter->iflags);CHKERRMPI(ierr);
30805342407SJunchao Zhang   ierr = PetscFree(counter);CHKERRMPI(ierr);
309480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
310480cf27aSJed Brown }
311480cf27aSJed Brown 
312480cf27aSJed Brown /*
31347435625SJed Brown   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
314da3039f7SJed Brown   calls MPI_Comm_free().
315da3039f7SJed Brown 
316da3039f7SJed Brown   This is the only entry point for breaking the links between inner and outer comms.
317480cf27aSJed Brown 
318ff0e51ddSBarry Smith   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.
319480cf27aSJed Brown 
32012801b39SBarry Smith   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()
321480cf27aSJed Brown 
322480cf27aSJed Brown */
32333779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
324480cf27aSJed Brown {
325480cf27aSJed Brown   PetscErrorCode                    ierr;
32633779a13SJunchao Zhang   union {MPI_Comm comm; void *ptr;} icomm;
327480cf27aSJed Brown 
328480cf27aSJed Brown   PetscFunctionBegin;
32912801b39SBarry Smith   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
330ec4fadc2SJed Brown   icomm.ptr = attr_val;
33176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
33233779a13SJunchao Zhang     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
33333779a13SJunchao Zhang     PetscMPIInt flg;
33433779a13SJunchao Zhang     union {MPI_Comm comm; void *ptr;} ocomm;
33547435625SJed Brown     ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr);
33633779a13SJunchao Zhang     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
33733779a13SJunchao Zhang     if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm's OuterComm attribute does not point to outer PETSc comm");
33833779a13SJunchao Zhang   }
33947435625SJed Brown   ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr);
34033779a13SJunchao Zhang   ierr = PetscInfo2(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);CHKERRMPI(ierr);
341da3039f7SJed Brown   PetscFunctionReturn(MPI_SUCCESS);
342b89831e5SBarry Smith }
343da3039f7SJed Brown 
344da3039f7SJed Brown /*
34533779a13SJunchao Zhang  * This is invoked on the inner comm when Petsc_InnerComm_Attr_Delete_Fn calls MPI_Comm_delete_attr().  It should not be reached any other way.
346da3039f7SJed Brown  */
34733779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
348da3039f7SJed Brown {
349da3039f7SJed Brown   PetscErrorCode ierr;
350da3039f7SJed Brown 
351da3039f7SJed Brown   PetscFunctionBegin;
35202c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr);
353480cf27aSJed Brown   PetscFunctionReturn(MPI_SUCCESS);
354480cf27aSJed Brown }
355480cf27aSJed Brown 
35633779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);
35742218b76SBarry Smith 
358951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
3598cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
3608cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
3618cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
362e39fd77fSBarry Smith #endif
363e39fd77fSBarry Smith 
3640f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;
36512801b39SBarry Smith 
366eb27c7c8SSatish Balay PETSC_INTERN int  PetscGlobalArgc;
367eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs;
3686ae9a8a6SBarry Smith int  PetscGlobalArgc   = 0;
36902c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL;
370dff31646SBarry Smith PetscSegBuffer PetscCitationsList;
371e5c89e4eSSatish Balay 
372dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void)
373051e4cf2SJed Brown {
374051e4cf2SJed Brown   PetscErrorCode ierr;
375051e4cf2SJed Brown 
376051e4cf2SJed Brown   PetscFunctionBegin;
377051e4cf2SJed Brown   ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr);
378a1601671SBarry 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);
379051e4cf2SJed 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);
380051e4cf2SJed Brown   PetscFunctionReturn(0);
381051e4cf2SJed Brown }
382e5c89e4eSSatish Balay 
3832d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */
3842d747510SLisandro Dalcin 
3852d747510SLisandro Dalcin PetscErrorCode  PetscSetProgramName(const char name[])
3862d747510SLisandro Dalcin {
3872d747510SLisandro Dalcin   PetscErrorCode ierr;
3882d747510SLisandro Dalcin 
3892d747510SLisandro Dalcin   PetscFunctionBegin;
390589a23caSBarry Smith   ierr  = PetscStrncpy(programname,name,sizeof(programname));CHKERRQ(ierr);
3912d747510SLisandro Dalcin   PetscFunctionReturn(0);
3922d747510SLisandro Dalcin }
3932d747510SLisandro Dalcin 
3942d747510SLisandro Dalcin /*@C
3952d747510SLisandro Dalcin     PetscGetProgramName - Gets the name of the running program.
3962d747510SLisandro Dalcin 
3972d747510SLisandro Dalcin     Not Collective
3982d747510SLisandro Dalcin 
3992d747510SLisandro Dalcin     Input Parameter:
4002d747510SLisandro Dalcin .   len - length of the string name
4012d747510SLisandro Dalcin 
4022d747510SLisandro Dalcin     Output Parameter:
4032d747510SLisandro Dalcin .   name - the name of the running program
4042d747510SLisandro Dalcin 
4052d747510SLisandro Dalcin    Level: advanced
4062d747510SLisandro Dalcin 
4072d747510SLisandro Dalcin     Notes:
4082d747510SLisandro Dalcin     The name of the program is copied into the user-provided character
4092d747510SLisandro Dalcin     array of length len.  On some machines the program name includes
4102d747510SLisandro Dalcin     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
4112d747510SLisandro Dalcin @*/
4122d747510SLisandro Dalcin PetscErrorCode  PetscGetProgramName(char name[],size_t len)
4132d747510SLisandro Dalcin {
4142d747510SLisandro Dalcin   PetscErrorCode ierr;
4152d747510SLisandro Dalcin 
4162d747510SLisandro Dalcin   PetscFunctionBegin;
4172d747510SLisandro Dalcin    ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr);
4182d747510SLisandro Dalcin   PetscFunctionReturn(0);
4192d747510SLisandro Dalcin }
4202d747510SLisandro Dalcin 
421e5c89e4eSSatish Balay /*@C
422e5c89e4eSSatish Balay    PetscGetArgs - Allows you to access the raw command line arguments anywhere
423e5c89e4eSSatish Balay      after PetscInitialize() is called but before PetscFinalize().
424e5c89e4eSSatish Balay 
425e5c89e4eSSatish Balay    Not Collective
426e5c89e4eSSatish Balay 
427e5c89e4eSSatish Balay    Output Parameters:
428e5c89e4eSSatish Balay +  argc - count of number of command line arguments
429e5c89e4eSSatish Balay -  args - the command line arguments
430e5c89e4eSSatish Balay 
431e5c89e4eSSatish Balay    Level: intermediate
432e5c89e4eSSatish Balay 
433e5c89e4eSSatish Balay    Notes:
434e5c89e4eSSatish Balay       This is usually used to pass the command line arguments into other libraries
435e5c89e4eSSatish Balay    that are called internally deep in PETSc or the application.
436e5c89e4eSSatish Balay 
437f177e3b1SBarry Smith       The first argument contains the program name as is normal for C arguments.
438f177e3b1SBarry Smith 
439793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()
440e5c89e4eSSatish Balay 
441e5c89e4eSSatish Balay @*/
4427087cfbeSBarry Smith PetscErrorCode  PetscGetArgs(int *argc,char ***args)
443e5c89e4eSSatish Balay {
444e5c89e4eSSatish Balay   PetscFunctionBegin;
44517186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
446e5c89e4eSSatish Balay   *argc = PetscGlobalArgc;
447e5c89e4eSSatish Balay   *args = PetscGlobalArgs;
448e5c89e4eSSatish Balay   PetscFunctionReturn(0);
449e5c89e4eSSatish Balay }
450e5c89e4eSSatish Balay 
451793721a6SBarry Smith /*@C
452793721a6SBarry Smith    PetscGetArguments - Allows you to access the  command line arguments anywhere
453793721a6SBarry Smith      after PetscInitialize() is called but before PetscFinalize().
454793721a6SBarry Smith 
455793721a6SBarry Smith    Not Collective
456793721a6SBarry Smith 
457793721a6SBarry Smith    Output Parameters:
458793721a6SBarry Smith .  args - the command line arguments
459793721a6SBarry Smith 
460793721a6SBarry Smith    Level: intermediate
461793721a6SBarry Smith 
462793721a6SBarry Smith    Notes:
463793721a6SBarry Smith       This does NOT start with the program name and IS null terminated (final arg is void)
464793721a6SBarry Smith 
465793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()
466793721a6SBarry Smith 
467793721a6SBarry Smith @*/
4687087cfbeSBarry Smith PetscErrorCode  PetscGetArguments(char ***args)
469793721a6SBarry Smith {
470793721a6SBarry Smith   PetscInt       i,argc = PetscGlobalArgc;
471793721a6SBarry Smith   PetscErrorCode ierr;
472793721a6SBarry Smith 
473793721a6SBarry Smith   PetscFunctionBegin;
47417186662SBarry Smith   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
4752d747510SLisandro Dalcin   if (!argc) {*args = NULL; PetscFunctionReturn(0);}
476785e854fSJed Brown   ierr = PetscMalloc1(argc,args);CHKERRQ(ierr);
477793721a6SBarry Smith   for (i=0; i<argc-1; i++) {
478793721a6SBarry Smith     ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr);
479793721a6SBarry Smith   }
4802d747510SLisandro Dalcin   (*args)[argc-1] = NULL;
481793721a6SBarry Smith   PetscFunctionReturn(0);
482793721a6SBarry Smith }
483793721a6SBarry Smith 
484793721a6SBarry Smith /*@C
485793721a6SBarry Smith    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()
486793721a6SBarry Smith 
487793721a6SBarry Smith    Not Collective
488793721a6SBarry Smith 
489793721a6SBarry Smith    Output Parameters:
490793721a6SBarry Smith .  args - the command line arguments
491793721a6SBarry Smith 
492793721a6SBarry Smith    Level: intermediate
493793721a6SBarry Smith 
494793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()
495793721a6SBarry Smith 
496793721a6SBarry Smith @*/
4977087cfbeSBarry Smith PetscErrorCode  PetscFreeArguments(char **args)
498793721a6SBarry Smith {
499793721a6SBarry Smith   PetscInt       i = 0;
500793721a6SBarry Smith   PetscErrorCode ierr;
501793721a6SBarry Smith 
502793721a6SBarry Smith   PetscFunctionBegin;
503a297a907SKarl Rupp   if (!args) PetscFunctionReturn(0);
504793721a6SBarry Smith   while (args[i]) {
505793721a6SBarry Smith     ierr = PetscFree(args[i]);CHKERRQ(ierr);
506793721a6SBarry Smith     i++;
507793721a6SBarry Smith   }
508793721a6SBarry Smith   ierr = PetscFree(args);CHKERRQ(ierr);
509793721a6SBarry Smith   PetscFunctionReturn(0);
510793721a6SBarry Smith }
511793721a6SBarry Smith 
51211525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS)
51330befbd2SBarry Smith #include <petscconfiginfo.h>
51430befbd2SBarry Smith 
51595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
51611525c0dSBarry Smith {
51711525c0dSBarry Smith   if (!PetscGlobalRank) {
51830befbd2SBarry Smith     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
51911525c0dSBarry Smith     int            port;
520ffbd1cfbSBarry Smith     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
52111525c0dSBarry Smith     size_t         applinelen,introlen;
52211525c0dSBarry Smith     PetscErrorCode ierr;
523ffbd1cfbSBarry Smith     char           sawsurl[256];
52411525c0dSBarry Smith 
525c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr);
52611525c0dSBarry Smith     if (flg) {
52711525c0dSBarry Smith       char  sawslog[PETSC_MAX_PATH_LEN];
52811525c0dSBarry Smith 
529589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL);CHKERRQ(ierr);
53011525c0dSBarry Smith       if (sawslog[0]) {
53111525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
53211525c0dSBarry Smith       } else {
53311525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
53411525c0dSBarry Smith       }
53511525c0dSBarry Smith     }
536589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg);CHKERRQ(ierr);
53711525c0dSBarry Smith     if (flg) {
53811525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
53911525c0dSBarry Smith     }
540c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr);
541ffbd1cfbSBarry Smith     if (selectport) {
542ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
543ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
544ffbd1cfbSBarry Smith     } else {
545c5929fdfSBarry Smith       ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr);
54611525c0dSBarry Smith       if (flg) {
54711525c0dSBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
54811525c0dSBarry Smith       }
549ffbd1cfbSBarry Smith     }
550589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg);CHKERRQ(ierr);
55111525c0dSBarry Smith     if (flg) {
55211525c0dSBarry Smith       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
55311525c0dSBarry Smith       ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr);
5549c1e0ce8SBarry Smith     } else {
555c5929fdfSBarry Smith       ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr);
5569c1e0ce8SBarry Smith       if (flg) {
557589a23caSBarry Smith         ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root));CHKERRQ(ierr);
5589c1e0ce8SBarry Smith         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr);
5599c1e0ce8SBarry Smith       }
56011525c0dSBarry Smith     }
561c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr);
56211525c0dSBarry Smith     if (flg2) {
56311525c0dSBarry Smith       char jsdir[PETSC_MAX_PATH_LEN];
56411525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
565589a23caSBarry Smith       ierr = PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root);CHKERRQ(ierr);
56611525c0dSBarry Smith       ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr);
56711525c0dSBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
56843da4ab2SBarry Smith       PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr);
56911525c0dSBarry Smith     }
570589a23caSBarry Smith     ierr = PetscGetProgramName(programname,sizeof(programname));CHKERRQ(ierr);
57111525c0dSBarry Smith     ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr);
57211525c0dSBarry Smith     introlen   = 4096 + applinelen;
57330a8c9c0SSurtai Han     applinelen += 1024;
57411525c0dSBarry Smith     ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr);
57511525c0dSBarry Smith     ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr);
57611525c0dSBarry Smith 
57711525c0dSBarry Smith     if (rootlocal) {
57811525c0dSBarry Smith       ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr);
57911525c0dSBarry Smith       ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr);
58011525c0dSBarry Smith     }
58176a34f28SBarry Smith     ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr);
58211525c0dSBarry Smith     if (rootlocal && help) {
583928bb9adSStefano 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);
58411525c0dSBarry Smith     } else if (help) {
585928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr);
58611525c0dSBarry Smith     } else {
587928bb9adSStefano Zampini       ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr);
58811525c0dSBarry Smith     }
589b0bb5815SBarry Smith     ierr = PetscFree(options);CHKERRQ(ierr);
59030befbd2SBarry Smith     ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
59111525c0dSBarry Smith     ierr = PetscSNPrintf(intro,introlen,"<body>\n"
592a8d69d7bSBarry 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"
593df62c222SBarry 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"
594928bb9adSStefano Zampini                                     "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr);
59543da4ab2SBarry Smith     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
59611525c0dSBarry Smith     ierr = PetscFree(intro);CHKERRQ(ierr);
59711525c0dSBarry Smith     ierr = PetscFree(appline);CHKERRQ(ierr);
598ffbd1cfbSBarry Smith     if (selectport) {
599aa573868SBarry Smith       PetscBool silent;
6007d812c46SBarry Smith 
6017d812c46SBarry Smith       ierr = SAWs_Initialize();
6027d812c46SBarry Smith       /* another process may have grabbed the port so keep trying */
6037d812c46SBarry Smith       while (ierr) {
6047d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
6057d812c46SBarry Smith         PetscStackCallSAWs(SAWs_Set_Port,(port));
6067d812c46SBarry Smith         ierr = SAWs_Initialize();
6077d812c46SBarry Smith       }
6087d812c46SBarry Smith 
609aa573868SBarry Smith       ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr);
610aa573868SBarry Smith       if (!silent) {
611ffbd1cfbSBarry Smith         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
612ffbd1cfbSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr);
613ffbd1cfbSBarry Smith       }
6147d812c46SBarry Smith     } else {
6157d812c46SBarry Smith       PetscStackCallSAWs(SAWs_Initialize,());
616aa573868SBarry Smith     }
6170af79b04SBarry Smith     ierr = PetscCitationsRegister("@TechReport{ saws,\n"
6180af79b04SBarry Smith                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
6190af79b04SBarry Smith                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
6200af79b04SBarry Smith                                   "  Institution = {Argonne National Laboratory},\n"
6210af79b04SBarry Smith                                   "  Year   = 2013\n}\n",NULL);CHKERRQ(ierr);
62211525c0dSBarry Smith   }
62311525c0dSBarry Smith   PetscFunctionReturn(0);
62411525c0dSBarry Smith }
62511525c0dSBarry Smith #endif
62611525c0dSBarry Smith 
6274dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
6284dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
6294dfee713SSatish Balay {
6304dfee713SSatish Balay   PetscFunctionBegin;
6314dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
6324dfee713SSatish Balay     /* see MPI.py for details on this bug */
6334dfee713SSatish Balay     (void) setenv("HWLOC_COMPONENTS","-x86",1);
6344dfee713SSatish Balay #endif
6354dfee713SSatish Balay   PetscFunctionReturn(0);
6364dfee713SSatish Balay }
6374dfee713SSatish Balay 
638a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
639a56f64adSBarry Smith #include <adios.h>
64022580e64SBarry Smith #include <adios_read.h>
6417b56e58cSSatish Balay int64_t Petsc_adios_group;
642a56f64adSBarry Smith #endif
64355d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
64455d657eeSBarry Smith #include <adios2_c.h>
64555d657eeSBarry Smith #endif
646cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
647cd1b99a6SBarry Smith #include <omp.h>
648f90b075cSBarry Smith PetscInt PetscNumOMPThreads;
649cd1b99a6SBarry Smith #endif
650a56f64adSBarry Smith 
651bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H)
652bc8a24ecSBarry Smith #include <dlfcn.h>
653bc8a24ecSBarry Smith #endif
654bc8a24ecSBarry Smith 
655e5c89e4eSSatish Balay /*@C
656e5c89e4eSSatish Balay    PetscInitialize - Initializes the PETSc database and MPI.
657e5c89e4eSSatish Balay    PetscInitialize() calls MPI_Init() if that has yet to be called,
658e5c89e4eSSatish Balay    so this routine should always be called near the beginning of
659e5c89e4eSSatish Balay    your program -- usually the very first line!
660e5c89e4eSSatish Balay 
661e5c89e4eSSatish Balay    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
662e5c89e4eSSatish Balay 
663e5c89e4eSSatish Balay    Input Parameters:
664e5c89e4eSSatish Balay +  argc - count of number of command line arguments
665e5c89e4eSSatish Balay .  args - the command line arguments
6660298fd71SBarry Smith .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
667fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
6680298fd71SBarry Smith -  help - [optional] Help message to print, use NULL for no message
669e5c89e4eSSatish Balay 
67005827820SBarry Smith    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
67105827820SBarry Smith    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
67205827820SBarry Smith    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
67305827820SBarry Smith    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
67405827820SBarry Smith    if different subcommunicators of the job are doing different things with PETSc.
675e5c89e4eSSatish Balay 
676e5c89e4eSSatish Balay    Options Database Keys:
6777ca660e7SBarry Smith +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
6787ca660e7SBarry Smith .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
679e5c89e4eSSatish Balay .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
6808a690491SBarry Smith .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
6818a690491SBarry Smith .  -on_error_abort - calls abort() when error detected (no traceback)
6828a690491SBarry Smith .  -on_error_mpiabort - calls MPI_abort() when error detected
6838a690491SBarry Smith .  -error_output_stderr - prints error messages to stderr instead of the default stdout
6848a690491SBarry Smith .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
685e5c89e4eSSatish Balay .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
686e5c89e4eSSatish Balay .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
687e5c89e4eSSatish Balay .  -stop_for_debugger - Print message on how to attach debugger manually to
688e5c89e4eSSatish Balay                         process and wait (-debugger_pause) seconds for attachment
68979dccf82SBarry Smith .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
69079dccf82SBarry Smith .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
69179dccf82SBarry Smith .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
692aee23540SBarry Smith .  -malloc_dump - prints a list of all unfreed memory at the end of the run
69392f119d6SBarry 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
69492f119d6SBarry Smith .  -malloc_view - show a list of all allocated memory during PetscFinalize()
69592f119d6SBarry Smith .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
69692f119d6SBarry Smith .  -fp_trap - Stops on floating point exceptions
697e5c89e4eSSatish Balay .  -no_signal_handler - Indicates not to trap error signals
698e5c89e4eSSatish Balay .  -shared_tmp - indicates /tmp directory is shared by all processors
699e5c89e4eSSatish Balay .  -not_shared_tmp - each processor has own /tmp
700e5c89e4eSSatish Balay .  -tmp - alternative name of /tmp directory
701e5c89e4eSSatish Balay .  -get_total_flops - returns total flops done by all processors
7020841954dSBarry Smith -  -memory_view - Print memory usage at end of run
703e5c89e4eSSatish Balay 
704e5c89e4eSSatish Balay    Options Database Keys for Profiling:
705a7f22e61SSatish Balay    See Users-Manual: ch_profiling for details.
706fe9b927eSVaclav Hapla +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
707217044c2SLisandro Dalcin .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
708217044c2SLisandro Dalcin         however it slows things down and gives a distorted view of the overall runtime.
709495fc317SBarry Smith .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
710e5c89e4eSSatish Balay         hangs without running in the debugger).  See PetscLogTraceBegin().
7119a9a5d4cSBarry Smith .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
71279dccf82SBarry Smith .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
7139a9a5d4cSBarry Smith .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
714495fc317SBarry Smith         summary is written to the file.  See PetscLogView().
715f5d6ab90SLisandro Dalcin .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
716495fc317SBarry Smith .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
717495fc317SBarry Smith .  -log [filename] - Logs basic profiline information  See PetscLogDump().
718c2f74817SBarry Smith .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
71987c3beb6SLisandro Dalcin .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
720c2f74817SBarry Smith -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code
721495fc317SBarry Smith 
722609bdbeeSBarry Smith     Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time
723e5c89e4eSSatish Balay 
724ffbd1cfbSBarry Smith    Options Database Keys for SAWs:
725ffbd1cfbSBarry Smith +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
726ffbd1cfbSBarry 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
727ffbd1cfbSBarry Smith                             this is useful when you are running many jobs that utilize SAWs at the same time
728ffbd1cfbSBarry Smith .  -saws_log <filename> - save a log of all SAWs communication
729ffbd1cfbSBarry Smith .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
730ffbd1cfbSBarry Smith -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files
731ffbd1cfbSBarry Smith 
732e5c89e4eSSatish Balay    Environmental Variables:
733e5c89e4eSSatish Balay +   PETSC_TMP - alternative tmp directory
734e5c89e4eSSatish Balay .   PETSC_SHARED_TMP - tmp is shared by all processes
735e5c89e4eSSatish Balay .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
7364a971ea4SToby Isaac .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
7374a971ea4SToby Isaac .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
738e5c89e4eSSatish Balay .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
739e5c89e4eSSatish Balay -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to
740e5c89e4eSSatish Balay 
741e5c89e4eSSatish Balay 
742e5c89e4eSSatish Balay    Level: beginner
743e5c89e4eSSatish Balay 
744e5c89e4eSSatish Balay    Notes:
745e5c89e4eSSatish Balay    If for some reason you must call MPI_Init() separately, call
746e5c89e4eSSatish Balay    it before PetscInitialize().
747e5c89e4eSSatish Balay 
748e5c89e4eSSatish Balay    Fortran Version:
749e5c89e4eSSatish Balay    In Fortran this routine has the format
750e5c89e4eSSatish Balay $       call PetscInitialize(file,ierr)
751e5c89e4eSSatish Balay 
752e5c89e4eSSatish Balay +   ierr - error return code
7530eb4c9c0SBarry Smith -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
754fc2bca9aSBarry Smith           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
755e5c89e4eSSatish Balay 
756e5c89e4eSSatish Balay    Important Fortran Note:
7570eb4c9c0SBarry Smith    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
7580298fd71SBarry Smith    null character string; you CANNOT just use NULL as
759a7f22e61SSatish Balay    in the C version. See Users-Manual: ch_fortran for details.
760e5c89e4eSSatish Balay 
76101cb0274SBarry Smith    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
76201cb0274SBarry Smith    calling PetscInitialize().
763e5c89e4eSSatish Balay 
76401cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
765e5c89e4eSSatish Balay 
766e5c89e4eSSatish Balay @*/
7677087cfbeSBarry Smith PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
768e5c89e4eSSatish Balay {
769e5c89e4eSSatish Balay   PetscErrorCode ierr;
7704bb5149bSJed Brown   PetscMPIInt    flag, size;
77119c5658aSBarry Smith   PetscBool      flg = PETSC_TRUE;
772e5c89e4eSSatish Balay   char           hostname[256];
773e5c89e4eSSatish Balay 
774e5c89e4eSSatish Balay   PetscFunctionBegin;
775e5c89e4eSSatish Balay   if (PetscInitializeCalled) PetscFunctionReturn(0);
7763d96e996SBarry Smith   /*
7773d96e996SBarry Smith       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
7783d96e996SBarry Smith       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
7793d96e996SBarry Smith         MPICH v3.1 (Released Feburary 2014)
7803d96e996SBarry Smith         IBM MPI v2.1 (December 2014)
7813d96e996SBarry Smith         Intel® MPI Library v5.0 (2014)
7823d96e996SBarry Smith         Cray MPT v7.0.0 (June 2014)
7833d96e996SBarry Smith       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
7843d96e996SBarry Smith       listed above and since that time are compatible.
7853d96e996SBarry Smith 
7863d96e996SBarry Smith       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
7873d96e996SBarry Smith       at compile time or runtime. Thus we will need to systematically track the allowed versions
7883d96e996SBarry Smith       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
7893d96e996SBarry Smith       to perform the checking.
7903d96e996SBarry Smith 
7913d96e996SBarry Smith       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).
7923d96e996SBarry Smith 
7933d96e996SBarry Smith       Questions:
7943d96e996SBarry Smith 
7953d96e996SBarry Smith         Should the checks for ABI incompatibility be only on the major version number below?
7963d96e996SBarry Smith         Presumably the output to stderr will be removed before a release.
7973d96e996SBarry Smith   */
7983d96e996SBarry Smith 
79919c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
80019c5658aSBarry Smith   {
80119c5658aSBarry Smith     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
80219c5658aSBarry Smith     PetscMPIInt mpilibraryversionlength;
80319c5658aSBarry Smith     ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr;
8043d96e996SBarry Smith     /* check for MPICH versions before MPI ABI initiative */
80519c5658aSBarry Smith #if defined(MPICH_VERSION)
8063d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000
80719c5658aSBarry Smith     {
80819c5658aSBarry Smith       char *ver,*lf;
80919c5658aSBarry Smith       flg = PETSC_FALSE;
81019c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr;
81119c5658aSBarry Smith       if (ver) {
81219c5658aSBarry Smith         ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr;
81319c5658aSBarry Smith         if (lf) {
81419c5658aSBarry Smith           *lf = 0;
81519c5658aSBarry Smith           ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr;
81619c5658aSBarry Smith         }
81719c5658aSBarry Smith       }
81819c5658aSBarry Smith       if (!flg) {
81919c5658aSBarry Smith         fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION);
8203d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
82119c5658aSBarry Smith       }
82219c5658aSBarry Smith     }
8233d96e996SBarry Smith #endif
8243d96e996SBarry 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?) */
82519c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION)
82619c5658aSBarry Smith     {
82719c5658aSBarry Smith       char *ver,bs[32],*bsf;
82819c5658aSBarry Smith       flg = PETSC_FALSE;
82919c5658aSBarry Smith       ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr;
83019c5658aSBarry Smith       if (ver) {
8312e924ca5SSatish Balay         PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
83219c5658aSBarry Smith         ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr;
83319c5658aSBarry Smith         if (bsf) flg = PETSC_TRUE;
83419c5658aSBarry Smith       }
83519c5658aSBarry Smith       if (!flg) {
83619c5658aSBarry 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);
8373d96e996SBarry Smith         return PETSC_ERR_MPI_LIB_INCOMP;
83819c5658aSBarry Smith       }
83919c5658aSBarry Smith     }
84019c5658aSBarry Smith #endif
84119c5658aSBarry Smith   }
84219c5658aSBarry Smith #endif
84319c5658aSBarry Smith 
844bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM)
845bc8a24ecSBarry Smith   {
846bc8a24ecSBarry Smith     PetscInt cnt = 0;
847bc8a24ecSBarry 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 */
848bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++;
849bc8a24ecSBarry Smith     if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++;
850bc8a24ecSBarry Smith     if (cnt > 1) {
851bc8a24ecSBarry Smith       fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
852bc8a24ecSBarry Smith       return PETSC_ERR_MPI_LIB_INCOMP;
853bc8a24ecSBarry Smith     }
854bc8a24ecSBarry Smith   }
855bc8a24ecSBarry Smith #endif
856bc8a24ecSBarry Smith 
857ae9b4142SLisandro Dalcin   /* these must be initialized in a routine, not as a constant declaration*/
858d89683f4Sbcordonn   PETSC_STDOUT = stdout;
859ae9b4142SLisandro Dalcin   PETSC_STDERR = stderr;
860e5c89e4eSSatish Balay 
861*8ad20175SVaclav Hapla   /* CHKERRQ can be used from now */
862*8ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_TRUE;
863*8ad20175SVaclav Hapla 
8640c30907bSSatish Balay   /* on Windows - set printf to default to printing 2 digit exponents */
8650c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
8660c30907bSSatish Balay   _set_output_format(_TWO_DIGIT_EXPONENT);
8670c30907bSSatish Balay #endif
8680c30907bSSatish Balay 
8694416b707SBarry Smith   ierr = PetscOptionsCreateDefault();CHKERRQ(ierr);
870e5c89e4eSSatish Balay 
871e5c89e4eSSatish Balay   /*
872e5c89e4eSSatish Balay      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
873e5c89e4eSSatish Balay      it that it sets args[0] on all processors to be args[0] on the first processor.
874e5c89e4eSSatish Balay   */
875e5c89e4eSSatish Balay   if (argc && *argc) {
876e5c89e4eSSatish Balay     ierr = PetscSetProgramName(**args);CHKERRQ(ierr);
877e5c89e4eSSatish Balay   } else {
878e5c89e4eSSatish Balay     ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr);
879e5c89e4eSSatish Balay   }
880e5c89e4eSSatish Balay 
881e5c89e4eSSatish Balay   ierr = MPI_Initialized(&flag);CHKERRQ(ierr);
882e5c89e4eSSatish Balay   if (!flag) {
883e32f2f54SBarry 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");
8844dfee713SSatish Balay     ierr = PetscPreMPIInit_Private();CHKERRQ(ierr);
8855e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD)
8865e765c61SJed Brown     {
8875e765c61SJed Brown       PetscMPIInt provided;
8886de5d289SStefano Zampini       ierr = MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided);CHKERRQ(ierr);
8895e765c61SJed Brown     }
8905e765c61SJed Brown #else
891e5c89e4eSSatish Balay     ierr = MPI_Init(argc,args);CHKERRQ(ierr);
8925e765c61SJed Brown #endif
893e5c89e4eSSatish Balay     PetscBeganMPI = PETSC_TRUE;
894e5c89e4eSSatish Balay   }
8954dfee713SSatish Balay 
896e5c89e4eSSatish Balay   if (argc && args) {
897e5c89e4eSSatish Balay     PetscGlobalArgc = *argc;
898e5c89e4eSSatish Balay     PetscGlobalArgs = *args;
899e5c89e4eSSatish Balay   }
900e5c89e4eSSatish Balay   PetscFinalizeCalled = PETSC_FALSE;
9015ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
9025ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
9035ad9ad5bSBarry Smith   ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
904ef19f930SBarry Smith   ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr);
905e5c89e4eSSatish Balay 
906a297a907SKarl Rupp   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
907d54338ecSKarl Rupp   ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr);
908e5c89e4eSSatish Balay 
9090f9be574SPeter Hill   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
91012801b39SBarry Smith     ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr);
91112801b39SBarry Smith     ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr);
9120f9be574SPeter Hill   }
91312801b39SBarry Smith 
914e5c89e4eSSatish Balay   /* Done after init due to a bug in MPICH-GM? */
915e5c89e4eSSatish Balay   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
916e5c89e4eSSatish Balay 
917e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr);
918e5c89e4eSSatish Balay   ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr);
919e5c89e4eSSatish Balay 
9208ad47952SJed Brown   MPIU_BOOL = MPI_INT;
9218ad47952SJed Brown   MPIU_ENUM = MPI_INT;
9227cdaf61dSJed Brown   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
923e316c87fSJed Brown   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
924e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
925e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
926e316c87fSJed Brown   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
927e316c87fSJed Brown #endif
928e316c87fSJed Brown   else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;}
9298ad47952SJed Brown 
930e5c89e4eSSatish Balay   /*
931e5c89e4eSSatish Balay      Initialized the global complex variable; this is because with
932e5c89e4eSSatish Balay      shared libraries the constructors for global variables
933e5c89e4eSSatish Balay      are not called; at least on IRIX.
934e5c89e4eSSatish Balay   */
935886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX)
936e5c89e4eSSatish Balay   {
937252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
93850f81f78SJed Brown     PetscComplex ic(0.0,1.0);
939e5c89e4eSSatish Balay     PETSC_i = ic;
940252ecd31SSatish Balay #else
94150f81f78SJed Brown     PETSC_i = _Complex_I;
942b7940d39SSatish Balay #endif
943762437b8SSatish Balay   }
944762437b8SSatish Balay 
9452c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
946e69cd0e6SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
947500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
948500d8756SSatish Balay   ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr);
949500d8756SSatish Balay   ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr);
9502c876bd9SBarry Smith #endif
951886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */
952e5c89e4eSSatish Balay 
953e5c89e4eSSatish Balay   /*
954e5c89e4eSSatish Balay      Create the PETSc MPI reduction operator that sums of the first
955e5c89e4eSSatish Balay      half of the entries and maxes the second half.
956e5c89e4eSSatish Balay   */
957367daffbSBarry Smith   ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr);
958e5c89e4eSSatish Balay 
959ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
960c90a1750SBarry Smith   ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr);
961c90a1750SBarry Smith   ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr);
9627c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX)
9638c764dc5SJose Roman   ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr);
9648c764dc5SJose Roman   ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr);
9658c764dc5SJose Roman #endif
966d9822059SBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
967d9822059SBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
968570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
969570b7f6dSBarry Smith   ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr);
970570b7f6dSBarry Smith   ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr);
971570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr);
972570b7f6dSBarry Smith   ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr);
973c90a1750SBarry Smith #endif
974c90a1750SBarry Smith 
975570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
976cca4cb22SSatish Balay   ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr);
977cca4cb22SSatish Balay #endif
978cca4cb22SSatish Balay 
979e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr);
980e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr);
981e5c89e4eSSatish Balay 
98240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
983e5c89e4eSSatish Balay   ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr);
984e5c89e4eSSatish Balay   ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr);
98544041f26SJed Brown #endif
986e5c89e4eSSatish Balay 
987e5c89e4eSSatish Balay   /*
988480cf27aSJed Brown      Attributes to be set on PETSc communicators
989480cf27aSJed Brown   */
99033779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr);
99133779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr);
99233779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr);
99333779a13SJunchao Zhang   ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr);
994480cf27aSJed Brown 
995480cf27aSJed Brown   /*
996e8fb0fc0SBarry Smith      Build the options database
997e5c89e4eSSatish Balay   */
998c5929fdfSBarry Smith   ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr);
999e5c89e4eSSatish Balay 
1000703f3eceSBarry Smith   /* call a second time so it can look in the options database */
1001703f3eceSBarry Smith   ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr);
10026dc8fec2Sbcordonn 
1003e5c89e4eSSatish Balay   /*
100457171095SVaclav Hapla      Check system options and print help
1005e5c89e4eSSatish Balay   */
100657171095SVaclav Hapla   ierr = PetscOptionsCheckInitial_Private(help);CHKERRQ(ierr);
1007e5c89e4eSSatish Balay 
1008d45a07a7SBarry Smith   ierr = PetscCitationsInitialize();CHKERRQ(ierr);
1009d45a07a7SBarry Smith 
1010e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
101111525c0dSBarry Smith   ierr = PetscInitializeSAWs(help);CHKERRQ(ierr);
1012f4202a44SBarry Smith #endif
1013f4202a44SBarry Smith 
1014e5c89e4eSSatish Balay   /*
1015e5c89e4eSSatish Balay      Load the dynamic libraries (on machines that support them), this registers all
1016e5c89e4eSSatish Balay      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
1017e5c89e4eSSatish Balay   */
1018e5c89e4eSSatish Balay   ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr);
1019e5c89e4eSSatish Balay 
1020e5c89e4eSSatish Balay   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
102102c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr);
1022e5c89e4eSSatish Balay   ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr);
102302c9f0b5SLisandro Dalcin   ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr);
1024cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP)
1025101491d0SBarry Smith   {
1026101491d0SBarry Smith     PetscBool omp_view_flag;
1027cd1b99a6SBarry Smith     char      *threads = getenv("OMP_NUM_THREADS");
1028e5c89e4eSSatish Balay 
1029cd1b99a6SBarry Smith    if (threads) {
103002c9f0b5SLisandro Dalcin      ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr);
1031f90b075cSBarry Smith      (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
1032cd1b99a6SBarry Smith    } else {
1033cd1b99a6SBarry Smith #define NMAX  10000
1034101491d0SBarry Smith      int          i;
1035cd1b99a6SBarry Smith       PetscScalar *x;
1036cd1b99a6SBarry Smith       ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr);
1037cd1b99a6SBarry Smith #pragma omp parallel for
1038cd1b99a6SBarry Smith       for (i=0; i<NMAX; i++) {
1039cd1b99a6SBarry Smith         x[i] = 0.0;
1040f90b075cSBarry Smith         PetscNumOMPThreads  = (PetscInt) omp_get_num_threads();
1041cd1b99a6SBarry Smith       }
1042cd1b99a6SBarry Smith       ierr = PetscFree(x);CHKERRQ(ierr);
104302c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1044101491d0SBarry Smith     }
1045101491d0SBarry Smith     ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr);
1046f90b075cSBarry 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);
1047101491d0SBarry Smith     ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr);
1048101491d0SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1049101491d0SBarry Smith     if (flg) {
105002c9f0b5SLisandro Dalcin       ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr);
1051f90b075cSBarry Smith       omp_set_num_threads((int)PetscNumOMPThreads);
1052101491d0SBarry Smith     }
1053101491d0SBarry Smith     if (omp_view_flag) {
1054f90b075cSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr);
1055cd1b99a6SBarry Smith     }
1056cd1b99a6SBarry Smith   }
1057cd1b99a6SBarry Smith #endif
1058ef6c6fedSBoyana Norris   /* Check the options database for options related to the options database itself */
1059c5929fdfSBarry Smith   ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr);
1060ef6c6fedSBoyana Norris 
1061951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
1062e39fd77fSBarry Smith   /*
1063e39fd77fSBarry Smith       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI
1064e39fd77fSBarry Smith 
1065e39fd77fSBarry Smith       Currently not used because it is not supported by MPICH.
1066e39fd77fSBarry Smith   */
106730815ce0SLisandro Dalcin   if (!PetscBinaryBigEndian()) {
10680298fd71SBarry Smith     ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr);
106930815ce0SLisandro Dalcin   }
1070951e3c8eSBarry Smith #endif
1071e39fd77fSBarry Smith 
107241c0b4b3SShri Abhyankar   /*
107341c0b4b3SShri Abhyankar       Setup building of stack frames for all function calls
107441c0b4b3SShri Abhyankar   */
1075ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1076e1167bb9SShri Abhyankar   ierr = PetscStackCreate();CHKERRQ(ierr);
1077e1167bb9SShri Abhyankar #endif
1078e1167bb9SShri Abhyankar 
10792d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
10802d53ad75SBarry Smith   ierr = PetscFPTCreate(10000);CHKERRQ(ierr);
10812d53ad75SBarry Smith #endif
10822d53ad75SBarry Smith 
10835e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC)
108487c3beb6SLisandro Dalcin   {
108587c3beb6SLisandro Dalcin     PetscViewer viewer;
108622e7e69cSBarry Smith     ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr);
10875e71baefSBarry Smith     if (flg) {
10885e71baefSBarry Smith       ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr);
108987c3beb6SLisandro Dalcin       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
109087c3beb6SLisandro Dalcin     }
10915e71baefSBarry Smith   }
10925e71baefSBarry Smith #endif
1093dff31646SBarry Smith 
109487c3beb6SLisandro Dalcin   flg = PETSC_TRUE;
109587c3beb6SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr);
109687c3beb6SLisandro Dalcin   if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);}
109787c3beb6SLisandro Dalcin 
1098a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
1099a56f64adSBarry Smith   ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr);
11007b56e58cSSatish Balay   ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr);
1101a56f64adSBarry Smith   ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr);
11029fc7e16cSBarry Smith   ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr);
1103a56f64adSBarry Smith #endif
110455d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
110555d657eeSBarry Smith #endif
1106a56f64adSBarry Smith 
1107301d30feSBarry Smith   /*
110857171095SVaclav Hapla       Set flag that we are completely initialized
1109301d30feSBarry Smith   */
1110301d30feSBarry Smith   PetscInitializeCalled = PETSC_TRUE;
11112db0e300SLisandro Dalcin 
11122db0e300SLisandro Dalcin   ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr);
11132db0e300SLisandro Dalcin   if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);}
1114301d30feSBarry Smith   PetscFunctionReturn(0);
1115e5c89e4eSSatish Balay }
1116e5c89e4eSSatish Balay 
11174097062eSBarry Smith #if defined(PETSC_USE_LOG)
111895c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects;
1119ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsCounts;
1120ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
112195c0884eSLisandro Dalcin PETSC_INTERN PetscBool   PetscObjectsLog;
11224097062eSBarry Smith #endif
1123e5c89e4eSSatish Balay 
1124008a6e76SBarry Smith /*
1125008a6e76SBarry Smith     Frees all the MPI types and operations that PETSc may have created
1126008a6e76SBarry Smith */
1127008a6e76SBarry Smith PetscErrorCode  PetscFreeMPIResources(void)
1128008a6e76SBarry Smith {
1129008a6e76SBarry Smith   PetscErrorCode ierr;
1130008a6e76SBarry Smith 
1131008a6e76SBarry Smith   PetscFunctionBegin;
1132008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128)
1133008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr);
1134008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1135008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr);
1136008a6e76SBarry Smith #endif
1137008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1138008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1139008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16)
1140008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr);
1141008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr);
1142008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr);
1143008a6e76SBarry Smith #endif
1144008a6e76SBarry Smith 
1145008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1146008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1147008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr);
1148008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr);
1149008a6e76SBarry Smith #endif
1150008a6e76SBarry Smith #endif
1151008a6e76SBarry Smith 
1152008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1153008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr);
1154008a6e76SBarry Smith #endif
1155008a6e76SBarry Smith 
1156008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr);
115740df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES)
1158008a6e76SBarry Smith   ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr);
1159008a6e76SBarry Smith #endif
1160008a6e76SBarry Smith   ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr);
1161008a6e76SBarry Smith   PetscFunctionReturn(0);
1162008a6e76SBarry Smith }
1163008a6e76SBarry Smith 
1164e5c89e4eSSatish Balay /*@C
1165e5c89e4eSSatish Balay    PetscFinalize - Checks for options to be called at the conclusion
1166e5c89e4eSSatish Balay    of the program. MPI_Finalize() is called only if the user had not
1167e5c89e4eSSatish Balay    called MPI_Init() before calling PetscInitialize().
1168e5c89e4eSSatish Balay 
1169e5c89e4eSSatish Balay    Collective on PETSC_COMM_WORLD
1170e5c89e4eSSatish Balay 
1171e5c89e4eSSatish Balay    Options Database Keys:
117226a7e8d4SBarry Smith +  -options_view - Calls PetscOptionsView()
1173e5c89e4eSSatish Balay .  -options_left - Prints unused options that remain in the database
11747eb1d149SBarry 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
1175e5c89e4eSSatish Balay .  -mpidump - Calls PetscMPIDump()
117692f119d6SBarry Smith .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1177e5c89e4eSSatish Balay .  -malloc_info - Prints total memory usage
117892f119d6SBarry Smith -  -malloc_view <optional filename> - Prints list of all memory allocated and where
1179e5c89e4eSSatish Balay 
1180e5c89e4eSSatish Balay    Level: beginner
1181e5c89e4eSSatish Balay 
1182e5c89e4eSSatish Balay    Note:
1183e5c89e4eSSatish Balay    See PetscInitialize() for more general runtime options.
1184e5c89e4eSSatish Balay 
118588c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1186e5c89e4eSSatish Balay @*/
11877087cfbeSBarry Smith PetscErrorCode  PetscFinalize(void)
1188e5c89e4eSSatish Balay {
1189e5c89e4eSSatish Balay   PetscErrorCode ierr;
11904bb5149bSJed Brown   PetscMPIInt    rank;
1191a8d2bbe5SBarry Smith   PetscInt       nopt;
11922bf49c77SBarry Smith   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1193dff31646SBarry Smith   PetscBool      flg;
119410463e74SBarry Smith #if defined(PETSC_USE_LOG)
119510463e74SBarry Smith   char           mname[PETSC_MAX_PATH_LEN];
119610463e74SBarry Smith #endif
1197e5c89e4eSSatish Balay 
1198e5c89e4eSSatish Balay   if (!PetscInitializeCalled) {
11994b09e917SBarry Smith     printf("PetscInitialize() must be called before PetscFinalize()\n");
12003cbbc5ffSBarry Smith     return(PETSC_ERR_ARG_WRONGSTATE);
1201e5c89e4eSSatish Balay   }
12023cbbc5ffSBarry Smith   PetscFunctionBegin;
12030298fd71SBarry Smith   ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr);
1204b022a5c1SBarry Smith 
12051f817a21SBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1206a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS)
120722580e64SBarry Smith   ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr);
1208a56f64adSBarry Smith   ierr = adios_finalize(rank);CHKERRQ(ierr);
1209a56f64adSBarry Smith #endif
121055d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2)
121155d657eeSBarry Smith #endif
1212c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr);
1213dff31646SBarry Smith   if (flg) {
12141f817a21SBarry Smith     char  *cits, filename[PETSC_MAX_PATH_LEN];
12151f817a21SBarry Smith     FILE  *fd = PETSC_STDOUT;
12161f817a21SBarry Smith 
1217589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL);CHKERRQ(ierr);
12181f817a21SBarry Smith     if (filename[0]) {
12191f817a21SBarry Smith       ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr);
12201f817a21SBarry Smith     }
1221dff31646SBarry Smith     ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr);
1222dff31646SBarry Smith     cits[0] = 0;
1223dff31646SBarry Smith     ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr);
12241f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr);
12251f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12261f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr);
12271f817a21SBarry Smith     ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr);
12281f817a21SBarry Smith     ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr);
1229dff31646SBarry Smith     ierr = PetscFree(cits);CHKERRQ(ierr);
1230dff31646SBarry Smith   }
1231dff31646SBarry Smith   ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr);
1232dff31646SBarry Smith 
1233c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
123404102261SBarry Smith   /* TextBelt is run for testing purposes only, please do not use this feature often */
123504102261SBarry Smith   {
123604102261SBarry Smith     PetscInt nmax = 2;
123704102261SBarry Smith     char     **buffs;
123804102261SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1239c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr);
124004102261SBarry Smith     if (flg1) {
124104102261SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
124204102261SBarry Smith       if (nmax == 1) {
124304102261SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
124404102261SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
124504102261SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
124604102261SBarry Smith       }
124704102261SBarry Smith       ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
124804102261SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
124904102261SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
125004102261SBarry Smith     }
125104102261SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
125204102261SBarry Smith   }
1253763ec7b1SBarry Smith   {
1254763ec7b1SBarry Smith     PetscInt nmax = 2;
1255763ec7b1SBarry Smith     char     **buffs;
1256763ec7b1SBarry Smith     ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr);
1257763ec7b1SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr);
1258763ec7b1SBarry Smith     if (flg1) {
1259763ec7b1SBarry Smith       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1260763ec7b1SBarry Smith       if (nmax == 1) {
1261763ec7b1SBarry Smith         ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr);
1262763ec7b1SBarry Smith         ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr);
1263763ec7b1SBarry Smith         ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr);
1264763ec7b1SBarry Smith       }
1265763ec7b1SBarry Smith       ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr);
1266763ec7b1SBarry Smith       ierr = PetscFree(buffs[0]);CHKERRQ(ierr);
1267763ec7b1SBarry Smith       ierr = PetscFree(buffs[1]);CHKERRQ(ierr);
1268763ec7b1SBarry Smith     }
1269763ec7b1SBarry Smith     ierr = PetscFree(buffs);CHKERRQ(ierr);
1270763ec7b1SBarry Smith   }
127104102261SBarry Smith #endif
127267234432SDmitry Karpeev   /*
127367234432SDmitry Karpeev     It should be safe to cancel the options monitors, since we don't expect to be setting options
127467234432SDmitry Karpeev     here (at least that are worth monitoring).  Monitors ought to be released so that they release
127567234432SDmitry Karpeev     whatever memory was allocated there before -malloc_dump reports unfreed memory.
127667234432SDmitry Karpeev   */
127767234432SDmitry Karpeev   ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr);
127804102261SBarry Smith 
12792d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
12802d53ad75SBarry Smith   ierr = PetscFPTDestroy();CHKERRQ(ierr);
12812d53ad75SBarry Smith #endif
12822d53ad75SBarry Smith 
12832d53ad75SBarry Smith 
1284e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1285dff31646SBarry Smith   flg = PETSC_FALSE;
1286c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr);
1287d5649816SBarry Smith   if (flg) {
1288e04113cfSBarry Smith     ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr);
1289d5649816SBarry Smith   }
1290d5649816SBarry Smith #endif
1291d5649816SBarry Smith 
1292681455b2SBarry Smith #if defined(PETSC_HAVE_X)
1293681455b2SBarry Smith   flg1 = PETSC_FALSE;
1294c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr);
1295681455b2SBarry Smith   if (flg1) {
1296681455b2SBarry Smith     /*  this is a crude hack, but better than nothing */
1297681455b2SBarry Smith     ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr);
1298681455b2SBarry Smith   }
1299681455b2SBarry Smith #endif
1300681455b2SBarry Smith 
130167584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
1302c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr);
1303e5c89e4eSSatish Balay   if (!flg2) {
130490d69ab7SBarry Smith     flg2 = PETSC_FALSE;
1305c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr);
1306e5c89e4eSSatish Balay   }
1307e5c89e4eSSatish Balay   if (flg2) {
13080841954dSBarry Smith     ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr);
1309e5c89e4eSSatish Balay   }
131067584ceeSBarry Smith #endif
1311e5c89e4eSSatish Balay 
1312e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
131390d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1314c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr);
1315e5c89e4eSSatish Balay   if (flg1) {
1316e5c89e4eSSatish Balay     PetscLogDouble flops = 0;
1317205a32c2SJed Brown     ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
1318e5c89e4eSSatish Balay     ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr);
1319e5c89e4eSSatish Balay   }
1320e5c89e4eSSatish Balay #endif
1321e5c89e4eSSatish Balay 
1322e5c89e4eSSatish Balay 
1323e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG)
1324e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE)
1325e5c89e4eSSatish Balay   mname[0] = 0;
1326589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1327e5c89e4eSSatish Balay   if (flg1) {
1328e5c89e4eSSatish Balay     if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);}
1329e5c89e4eSSatish Balay     else          {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);}
1330e5c89e4eSSatish Balay   }
1331e5c89e4eSSatish Balay #endif
1332356e5837SBarry Smith #endif
1333a297a907SKarl Rupp 
1334dd710f27SBarry Smith   /*
1335dd710f27SBarry Smith      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1336dd710f27SBarry Smith   */
1337dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1338dd710f27SBarry Smith 
1339356e5837SBarry Smith #if defined(PETSC_USE_LOG)
134087c3beb6SLisandro Dalcin   ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr);
1341f14045dbSBarry Smith   ierr = PetscLogViewFromOptions();CHKERRQ(ierr);
134287c3beb6SLisandro Dalcin   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
134387c3beb6SLisandro Dalcin 
1344356e5837SBarry Smith   mname[0] = 0;
1345589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1346e5c89e4eSSatish Balay   if (flg1) {
134791eabc43SBarry Smith     PetscViewer viewer;
134820a8bfc3SBarry Smith     ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr);
134991eabc43SBarry Smith     if (mname[0]) {
135091eabc43SBarry Smith       ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr);
135191eabc43SBarry Smith       ierr = PetscLogView(viewer);CHKERRQ(ierr);
13526bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
135333f85c2fSBarry Smith     } else {
135433f85c2fSBarry Smith       viewer = PETSC_VIEWER_STDOUT_WORLD;
13559a9a5d4cSBarry Smith       ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
135633f85c2fSBarry Smith       ierr   = PetscLogView(viewer);CHKERRQ(ierr);
13579a9a5d4cSBarry Smith       ierr   = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
135833f85c2fSBarry Smith     }
1359e5c89e4eSSatish Balay   }
1360a297a907SKarl Rupp 
1361dd710f27SBarry Smith   /*
1362dd710f27SBarry Smith      Free any objects created by the last block of code.
1363dd710f27SBarry Smith   */
1364dd710f27SBarry Smith   ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr);
1365dd710f27SBarry Smith 
1366dd710f27SBarry Smith   mname[0] = 0;
1367589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1);CHKERRQ(ierr);
1368589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2);CHKERRQ(ierr);
13697ff663adSLisandro Dalcin   if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);}
1370e5c89e4eSSatish Balay #endif
137110463e74SBarry Smith 
137233f85c2fSBarry Smith   ierr = PetscStackDestroy();CHKERRQ(ierr);
137310463e74SBarry Smith 
137490d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1375c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr);
1376e5c89e4eSSatish Balay   if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);}
137790d69ab7SBarry Smith   flg1 = PETSC_FALSE;
1378c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr);
1379e5c89e4eSSatish Balay   if (flg1) {
1380e5c89e4eSSatish Balay     ierr = PetscMPIDump(stdout);CHKERRQ(ierr);
1381e5c89e4eSSatish Balay   }
138290d69ab7SBarry Smith   flg1 = PETSC_FALSE;
138390d69ab7SBarry Smith   flg2 = PETSC_FALSE;
13848bb29257SSatish Balay   /* preemptive call to avoid listing this option in options table as unused */
1385c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr);
1386c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1387c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr);
1388e4c476e2SSatish Balay 
1389e5c89e4eSSatish Balay   if (flg2) {
1390be56827dSJed Brown     PetscViewer viewer;
139102ba9f54SBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
139202ba9f54SBarry Smith     ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1393c5929fdfSBarry Smith     ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1394be56827dSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1395e5c89e4eSSatish Balay   }
1396e5c89e4eSSatish Balay 
1397e5c89e4eSSatish Balay   /* to prevent PETSc -options_left from warning */
1398c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr);
1399c5929fdfSBarry Smith   ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr);
1400e5c89e4eSSatish Balay 
140133fc4174SSatish Balay   flg3 = PETSC_FALSE; /* default value is required */
1402c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr);
1403cf9c20a2SJed Brown   if (PetscDefined(USE_DEBUG) && !flg1) flg3 = PETSC_TRUE;
1404e5c89e4eSSatish Balay   if (flg3) {
14053de2bfdfSBarry Smith     if (!flg2 && flg1) { /* have not yet printed the options */
1406be56827dSJed Brown       PetscViewer viewer;
140702ba9f54SBarry Smith       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
140802ba9f54SBarry Smith       ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1409c5929fdfSBarry Smith       ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1410be56827dSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1411e5c89e4eSSatish Balay     }
14123de2bfdfSBarry Smith     ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr);
14133de2bfdfSBarry Smith     if (nopt) {
14143de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr);
14153de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr);
14163de2bfdfSBarry Smith       if (nopt == 1) {
1417e5c89e4eSSatish Balay         ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr);
1418e5c89e4eSSatish Balay       } else {
14197582186dSLisandro Dalcin         ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr);
1420e5c89e4eSSatish Balay       }
14213de2bfdfSBarry Smith     } else if (flg3 && flg1) {
14223de2bfdfSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr);
1423df12ba86SBarry Smith     }
1424c5929fdfSBarry Smith     ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr);
1425e5c89e4eSSatish Balay   }
1426e5c89e4eSSatish Balay 
1427e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1428d45a07a7SBarry Smith   if (!PetscGlobalRank) {
142987f587eeSBarry Smith     ierr = PetscStackSAWsViewOff();CHKERRQ(ierr);
143016ad0300SBarry Smith     PetscStackCallSAWs(SAWs_Finalize,());
1431d45a07a7SBarry Smith   }
1432ec957eceSBarry Smith #endif
1433ec957eceSBarry Smith 
14344097062eSBarry Smith #if defined(PETSC_USE_LOG)
143510463e74SBarry Smith   /*
1436dbc8283eSBarry Smith        List all objects the user may have forgot to free
14372eff7a51SBarry Smith   */
143805df10baSBarry Smith   if (PetscObjectsLog) {
1439c5929fdfSBarry Smith     ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
1440a64a8e02SBarry Smith     if (flg1) {
1441a64a8e02SBarry Smith       MPI_Comm local_comm;
14427eb1d149SBarry Smith       char     string[64];
1443a64a8e02SBarry Smith 
1444589a23caSBarry Smith       ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL);CHKERRQ(ierr);
1445a64a8e02SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1446a64a8e02SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
14477eb1d149SBarry Smith       ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
1448a64a8e02SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1449a64a8e02SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
14500a1571b3SBarry Smith     }
145105df10baSBarry Smith   }
14524097062eSBarry Smith #endif
14534097062eSBarry Smith 
14544097062eSBarry Smith #if defined(PETSC_USE_LOG)
1455dbc8283eSBarry Smith   PetscObjectsCounts    = 0;
1456dbc8283eSBarry Smith   PetscObjectsMaxCounts = 0;
1457a297a907SKarl Rupp   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);
14584097062eSBarry Smith #endif
14592eff7a51SBarry Smith 
146033f85c2fSBarry Smith   /*
146133f85c2fSBarry Smith      Destroy any packages that registered a finalize
146233f85c2fSBarry Smith   */
146333f85c2fSBarry Smith   ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr);
146433f85c2fSBarry Smith 
1465101409b8SToby Isaac #if defined(PETSC_USE_LOG)
1466fa2bb9feSLisandro Dalcin   ierr = PetscLogFinalize();CHKERRQ(ierr);
1467101409b8SToby Isaac #endif
1468101409b8SToby Isaac 
146933f85c2fSBarry Smith   /*
147048dd1dffSBarry Smith      Print PetscFunctionLists that have not been properly freed
147148dd1dffSBarry Smith 
147237e93019SBarry Smith   ierr = PetscFunctionListPrintAll();CHKERRQ(ierr);
147348dd1dffSBarry Smith   */
147437e93019SBarry Smith 
14754028d114SSatish Balay   if (petsc_history) {
1476f3dea69dSBarry Smith     ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr);
147702c9f0b5SLisandro Dalcin     petsc_history = NULL;
1478e5c89e4eSSatish Balay   }
14799de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
1480e94e781bSJacob Faibussowitsch   ierr = PetscInfoDestroy();CHKERRQ(ierr);
1481e5c89e4eSSatish Balay 
148267584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY)
148392f119d6SBarry Smith   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1484e5c89e4eSSatish Balay     char fname[PETSC_MAX_PATH_LEN];
148592f119d6SBarry Smith     char sname[PETSC_MAX_PATH_LEN];
1486e5c89e4eSSatish Balay     FILE *fd;
1487ed9cf6e9SBarry Smith     int  err;
1488e5c89e4eSSatish Balay 
1489dc92acbaSJed Brown     flg2 = PETSC_FALSE;
149092f119d6SBarry Smith     flg3 = PETSC_FALSE;
1491cf9c20a2SJed Brown     if (PetscDefined(USE_DEBUG)) {ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr);}
149292f119d6SBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr);
149392f119d6SBarry Smith     fname[0] = 0;
1494589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1495e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1496e5c89e4eSSatish Balay 
1497589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1498e32f2f54SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1499e5c89e4eSSatish Balay       ierr = PetscMallocDump(fd);CHKERRQ(ierr);
1500ed9cf6e9SBarry Smith       err  = fclose(fd);
1501e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
150292f119d6SBarry Smith     } else if (flg1 || flg2 || flg3) {
1503e5c89e4eSSatish Balay       MPI_Comm local_comm;
1504e5c89e4eSSatish Balay 
1505e5c89e4eSSatish Balay       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
1506e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
1507e5c89e4eSSatish Balay       ierr = PetscMallocDump(stdout);CHKERRQ(ierr);
1508e5c89e4eSSatish Balay       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
1509e5c89e4eSSatish Balay       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1510e5c89e4eSSatish Balay     }
1511e5c89e4eSSatish Balay     fname[0] = 0;
1512589a23caSBarry Smith     ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1);CHKERRQ(ierr);
1513e5c89e4eSSatish Balay     if (flg1 && fname[0]) {
1514e5c89e4eSSatish Balay 
1515589a23caSBarry Smith       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
151692f119d6SBarry Smith       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
151792f119d6SBarry Smith       ierr = PetscMallocView(fd);CHKERRQ(ierr);
1518ed9cf6e9SBarry Smith       err  = fclose(fd);
1519e32f2f54SBarry Smith       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
152092f119d6SBarry Smith     } else if (flg1) {
152192f119d6SBarry Smith       MPI_Comm local_comm;
152292f119d6SBarry Smith 
152392f119d6SBarry Smith       ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
152492f119d6SBarry Smith       ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
152592f119d6SBarry Smith       ierr = PetscMallocView(stdout);CHKERRQ(ierr);
152692f119d6SBarry Smith       ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
152792f119d6SBarry Smith       ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
1528e5c89e4eSSatish Balay     }
1529e5c89e4eSSatish Balay   }
153067584ceeSBarry Smith #endif
153120e2c332SMatthew G. Knepley 
15325486ca60SMatthew G. Knepley   /*
15335486ca60SMatthew G. Knepley      Close any open dynamic libraries
15345486ca60SMatthew G. Knepley   */
15355486ca60SMatthew G. Knepley   ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr);
15365486ca60SMatthew G. Knepley 
1537e5c89e4eSSatish Balay   /* Can be destroyed only after all the options are used */
15384416b707SBarry Smith   ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr);
1539e5c89e4eSSatish Balay 
1540e5c89e4eSSatish Balay   PetscGlobalArgc = 0;
154102c9f0b5SLisandro Dalcin   PetscGlobalArgs = NULL;
1542e5c89e4eSSatish Balay 
1543008a6e76SBarry Smith   ierr = PetscFreeMPIResources();CHKERRQ(ierr);
1544e5c89e4eSSatish Balay 
1545dbc8283eSBarry Smith   /*
1546efb80d3cSBarry Smith      Destroy any known inner MPI_Comm's and attributes pointing to them
1547efb80d3cSBarry Smith      Note this will not destroy any new communicators the user has created.
1548efb80d3cSBarry Smith 
1549efb80d3cSBarry Smith      If all PETSc objects were not destroyed those left over objects will have hanging references to
1550efb80d3cSBarry Smith      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1551dbc8283eSBarry Smith  */
1552b770b1f6SSatish Balay   {
1553dbc8283eSBarry Smith     PetscCommCounter *counter;
1554dbc8283eSBarry Smith     PetscMPIInt      flg;
1555dbc8283eSBarry Smith     MPI_Comm         icomm;
1556265f3f35SJed Brown     union {MPI_Comm comm; void *ptr;} ucomm;
155747435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1558dbc8283eSBarry Smith     if (flg) {
1559265f3f35SJed Brown       icomm = ucomm.comm;
156047435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1561dbc8283eSBarry 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");
1562dbc8283eSBarry Smith 
156347435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr);
156447435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1565efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1566dbc8283eSBarry Smith     }
156747435625SJed Brown     ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr);
1568dbc8283eSBarry Smith     if (flg) {
1569265f3f35SJed Brown       icomm = ucomm.comm;
157047435625SJed Brown       ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
1571dbc8283eSBarry 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");
1572dbc8283eSBarry Smith 
157347435625SJed Brown       ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr);
157447435625SJed Brown       ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr);
1575efb80d3cSBarry Smith       ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr);
1576dbc8283eSBarry Smith     }
1577b770b1f6SSatish Balay   }
1578dbc8283eSBarry Smith 
157947435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr);
158047435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr);
158147435625SJed Brown   ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr);
15825f7487a0SJunchao Zhang   ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr);
1583480cf27aSJed Brown 
15845ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr);
15855ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr);
15865ad9ad5bSBarry Smith   ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr);
1587ef19f930SBarry Smith   ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr);
1588ef19f930SBarry Smith 
1589e5c89e4eSSatish Balay   if (PetscBeganMPI) {
159099608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED)
159199b1327fSBarry Smith     PetscMPIInt flag;
159299b1327fSBarry Smith     ierr = MPI_Finalized(&flag);CHKERRQ(ierr);
1593e32f2f54SBarry Smith     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
159499608316SBarry Smith #endif
1595e5c89e4eSSatish Balay     ierr = MPI_Finalize();CHKERRQ(ierr);
1596e5c89e4eSSatish Balay   }
1597e5c89e4eSSatish Balay /*
1598e5c89e4eSSatish Balay 
1599e5c89e4eSSatish Balay      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1600e5c89e4eSSatish Balay    the communicator has some outstanding requests on it. Specifically if the
1601e5c89e4eSSatish Balay    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1602e5c89e4eSSatish Balay    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1603e5c89e4eSSatish Balay    is never freed as it should be. Thus one may obtain messages of the form
16040e5e90baSSatish Balay    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1605e5c89e4eSSatish Balay    memory was not freed.
1606e5c89e4eSSatish Balay 
1607e5c89e4eSSatish Balay */
16081d1a0024SBarry Smith   ierr = PetscMallocClear();CHKERRQ(ierr);
1609a297a907SKarl Rupp 
1610*8ad20175SVaclav Hapla   PetscErrorHandlingInitialized = PETSC_FALSE;
1611e5c89e4eSSatish Balay   PetscInitializeCalled = PETSC_FALSE;
1612e5c89e4eSSatish Balay   PetscFinalizeCalled   = PETSC_TRUE;
16133db9a53dSBarry Smith   PetscFunctionReturn(0);
1614e5c89e4eSSatish Balay }
1615e5c89e4eSSatish Balay 
161643db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_)
16178cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b)
161843db4dbbSBarry Smith {
161943db4dbbSBarry Smith   if (*a == *b) return 1;
162043db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
162143db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
162243db4dbbSBarry Smith   return 0;
162343db4dbbSBarry Smith }
1624a70650f6SBarry Smith #endif
162543db4dbbSBarry Smith 
162643db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame)
16278cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b)
162843db4dbbSBarry Smith {
162943db4dbbSBarry Smith   if (*a == *b) return 1;
163043db4dbbSBarry Smith   if (*a + 32 == *b) return 1;
163143db4dbbSBarry Smith   if (*a - 32 == *b) return 1;
163243db4dbbSBarry Smith   return 0;
163343db4dbbSBarry Smith }
163443db4dbbSBarry Smith #endif
1635