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*/ 6665c2dedSJed Brown #include <petscviewer.h> 7a0e72f99SJunchao Zhang 86edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS) 96edef35eSSatish Balay #include <petsc/private/valgrind/valgrind.h> 106edef35eSSatish Balay #endif 116edef35eSSatish Balay 1285649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN) 1385649d77SJunchao Zhang #include <petsc/private/fortranimpl.h> 1485649d77SJunchao Zhang #endif 1585649d77SJunchao Zhang 1656883f60SBarry Smith #if defined(PETSC_USE_GCOV) 1756883f60SBarry Smith EXTERN_C_BEGIN 18aaf3846cSSatish Balay #if defined(PETSC_HAVE___GCOV_DUMP) 19aaf3846cSSatish Balay #define __gcov_flush(x) __gcov_dump(x) 20aaf3846cSSatish Balay #endif 2156883f60SBarry Smith void __gcov_flush(void); 2256883f60SBarry Smith EXTERN_C_END 2356883f60SBarry Smith #endif 248101f56cSMatthew Knepley 252d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 2695c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData; 272d53ad75SBarry Smith PetscFPT PetscFPTData = 0; 282d53ad75SBarry Smith #endif 292d53ad75SBarry Smith 3027104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 3116ad0300SBarry Smith #include <petscviewersaws.h> 32a6790183SBarry Smith #endif 33eb02082bSJunchao Zhang 34e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/ 35e5c89e4eSSatish Balay 3695c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history; 37e5c89e4eSSatish Balay 3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void); 3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void); 4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void); 4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int); 4295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int); 4395c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**); 440069ddf5SShri Abhyankar 456de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */ 46e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL; 4727104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD) 486de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED; 496de5d289SStefano Zampini #else 506de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0; 516de5d289SStefano Zampini #endif 52e5c89e4eSSatish Balay 53480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval = MPI_KEYVAL_INVALID; 54480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID; 55480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID; 565f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval = MPI_KEYVAL_INVALID; 57480cf27aSJed Brown 58e5c89e4eSSatish Balay /* 59e5c89e4eSSatish Balay Declare and set all the string names of the PETSc enums 60e5c89e4eSSatish Balay */ 6102c9f0b5SLisandro Dalcin const char *const PetscBools[] = {"FALSE","TRUE","PetscBool","PETSC_",NULL}; 6202c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL}; 63e5c89e4eSSatish Balay 64ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE; 65ace3abfcSBarry Smith PetscBool PetscPreLoadingOn = PETSC_FALSE; 660f8e0872SSatish Balay 67a2f94806SJed Brown PetscInt PetscHotRegionDepth; 68a2f94806SJed Brown 696edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE; 706edef35eSSatish Balay 71b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY) 72b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen; 73b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout; 74b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr; 75b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock; 76b22622e2STadeu Manoel #endif 77b22622e2STadeu Manoel 78e5c89e4eSSatish Balay /* 79945d1669SBarry Smith PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args 8072a42c3cSBarry Smith 8172a42c3cSBarry Smith Collective 8272a42c3cSBarry Smith 8372a42c3cSBarry Smith Level: advanced 8472a42c3cSBarry Smith 8595452b02SPatrick Sanan Notes: 86a56f64adSBarry Smith this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to 870f11a792SBarry Smith indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to 88a56f64adSBarry Smith be called multiple times from Julia without the problem of trying to initialize MPI more than once. 890f11a792SBarry Smith 90a56f64adSBarry Smith Developer Note: Turns off PETSc signal handling to allow Julia to manage signals 911ea5a559SBarry Smith 9272a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments() 930f11a792SBarry Smith */ 94945d1669SBarry Smith PetscErrorCode PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help) 9572a42c3cSBarry Smith { 9672a42c3cSBarry Smith int myargc = argc; 9772a42c3cSBarry Smith char **myargs = args; 9872a42c3cSBarry Smith 9972a42c3cSBarry Smith PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&myargc,&myargs,filename,help)); 1019566063dSJacob Faibussowitsch PetscCall(PetscPopSignalHandler()); 102df413903SBarry Smith PetscBeganMPI = PETSC_FALSE; 10327104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 10472a42c3cSBarry Smith } 10572a42c3cSBarry Smith 106f0865b08SBarry Smith /* 107a56f64adSBarry Smith Used by Julia interface to get communicator 108f0865b08SBarry Smith */ 109945d1669SBarry Smith PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm) 110f0865b08SBarry Smith { 111f0865b08SBarry Smith PetscFunctionBegin; 11227104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscValidPointer(comm,1); 113f0865b08SBarry Smith *comm = PETSC_COMM_SELF; 114f0865b08SBarry Smith PetscFunctionReturn(0); 115f0865b08SBarry Smith } 116f0865b08SBarry Smith 117e5c89e4eSSatish Balay /*@C 118e5c89e4eSSatish Balay PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without 119e5c89e4eSSatish Balay the command line arguments. 120e5c89e4eSSatish Balay 121e5c89e4eSSatish Balay Collective 122e5c89e4eSSatish Balay 123e5c89e4eSSatish Balay Level: advanced 124e5c89e4eSSatish Balay 125e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran() 126e5c89e4eSSatish Balay @*/ 1277087cfbeSBarry Smith PetscErrorCode PetscInitializeNoArguments(void) 128e5c89e4eSSatish Balay { 129e5c89e4eSSatish Balay int argc = 0; 13002c9f0b5SLisandro Dalcin char **args = NULL; 131e5c89e4eSSatish Balay 132e5c89e4eSSatish Balay PetscFunctionBegin; 1339566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,NULL,NULL)); 13439a651e2SJacob Faibussowitsch PetscFunctionReturn(0); 135e5c89e4eSSatish Balay } 136e5c89e4eSSatish Balay 137e5c89e4eSSatish Balay /*@ 138e5c89e4eSSatish Balay PetscInitialized - Determine whether PETSc is initialized. 139e5c89e4eSSatish Balay 14093b6d2d1SJed Brown Level: beginner 141e5c89e4eSSatish Balay 142e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran() 143e5c89e4eSSatish Balay @*/ 1447087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized) 145e5c89e4eSSatish Balay { 14627104ee2SJacob Faibussowitsch PetscFunctionBegin; 14727104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized,1); 148e5c89e4eSSatish Balay *isInitialized = PetscInitializeCalled; 14927104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 150e5c89e4eSSatish Balay } 151e5c89e4eSSatish Balay 152e5c89e4eSSatish Balay /*@ 153e5c89e4eSSatish Balay PetscFinalized - Determine whether PetscFinalize() has been called yet 154e5c89e4eSSatish Balay 155e5c89e4eSSatish Balay Level: developer 156e5c89e4eSSatish Balay 157e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran() 158e5c89e4eSSatish Balay @*/ 1597087cfbeSBarry Smith PetscErrorCode PetscFinalized(PetscBool *isFinalized) 160e5c89e4eSSatish Balay { 16127104ee2SJacob Faibussowitsch PetscFunctionBegin; 16227104ee2SJacob Faibussowitsch if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized,1); 163e5c89e4eSSatish Balay *isFinalized = PetscFinalizeCalled; 16427104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 165e5c89e4eSSatish Balay } 166e5c89e4eSSatish Balay 16757171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []); 168e5c89e4eSSatish Balay 169e5c89e4eSSatish Balay /* 170e5c89e4eSSatish Balay This function is the MPI reduction operation used to compute the sum of the 171e5c89e4eSSatish Balay first half of the datatype and the max of the second half. 172e5c89e4eSSatish Balay */ 173367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0; 174e5c89e4eSSatish Balay 175367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype) 176e5c89e4eSSatish Balay { 177e5c89e4eSSatish Balay PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt; 178e5c89e4eSSatish Balay 179e5c89e4eSSatish Balay PetscFunctionBegin; 180e5c89e4eSSatish Balay if (*datatype != MPIU_2INT) { 181e5c89e4eSSatish Balay (*PetscErrorPrintf)("Can only handle MPIU_2INT data types"); 18241e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 183e5c89e4eSSatish Balay } 184e5c89e4eSSatish Balay 185e5c89e4eSSatish Balay for (i=0; i<count; i++) { 186e5c89e4eSSatish Balay xout[2*i] = PetscMax(xout[2*i],xin[2*i]); 187e5c89e4eSSatish Balay xout[2*i+1] += xin[2*i+1]; 188e5c89e4eSSatish Balay } 189812af9f3SBarry Smith PetscFunctionReturnVoid(); 190e5c89e4eSSatish Balay } 191e5c89e4eSSatish Balay 192e5c89e4eSSatish Balay /* 193e5c89e4eSSatish Balay Returns the max of the first entry owned by this processor and the 194e5c89e4eSSatish Balay sum of the second entry. 195b693b147SBarry Smith 19676f543a4SJed Brown The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero 197367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths 198b693b147SBarry Smith there would be no place to store the both needed results. 199e5c89e4eSSatish Balay */ 20076ec1555SBarry Smith PetscErrorCode PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum) 201e5c89e4eSSatish Balay { 202e5c89e4eSSatish Balay PetscFunctionBegin; 203d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 204d6e4c47cSJed Brown { 205d6e4c47cSJed Brown struct {PetscInt max,sum;} work; 2069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm)); 207d6e4c47cSJed Brown *max = work.max; 208d6e4c47cSJed Brown *sum = work.sum; 209d6e4c47cSJed Brown } 210d6e4c47cSJed Brown #else 211d6e4c47cSJed Brown { 212d6e4c47cSJed Brown PetscMPIInt size,rank; 213d6e4c47cSJed Brown struct {PetscInt max,sum;} *work; 2149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 2169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&work)); 2171c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm)); 2186ac3741eSJed Brown *max = work[rank].max; 2196ac3741eSJed Brown *sum = work[rank].sum; 2209566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 221d6e4c47cSJed Brown } 222d6e4c47cSJed Brown #endif 223e5c89e4eSSatish Balay PetscFunctionReturn(0); 224e5c89e4eSSatish Balay } 225e5c89e4eSSatish Balay 226e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/ 227e5c89e4eSSatish Balay 228de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 22906a205a8SBarry Smith MPI_Op MPIU_SUM = 0; 230e5c89e4eSSatish Balay 231092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 232e5c89e4eSSatish Balay { 233e5c89e4eSSatish Balay PetscInt i,count = *cnt; 234e5c89e4eSSatish Balay 235e5c89e4eSSatish Balay PetscFunctionBegin; 2367c2de775SJed Brown if (*datatype == MPIU_REAL) { 237e2e03761SBarry Smith PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 238a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] += xin[i]; 2397c2de775SJed Brown } 2407c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 241c5481aeeSJose E. Roman else if (*datatype == MPIU_COMPLEX) { 2427c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 243a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] += xin[i]; 2447c2de775SJed Brown } 2457c2de775SJed Brown #endif 2467c2de775SJed Brown else { 2477c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 24841e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 249e2e03761SBarry Smith } 250812af9f3SBarry Smith PetscFunctionReturnVoid(); 251e5c89e4eSSatish Balay } 252e5c89e4eSSatish Balay #endif 253e5c89e4eSSatish Balay 254570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 255d9822059SBarry Smith MPI_Op MPIU_MAX = 0; 256d9822059SBarry Smith MPI_Op MPIU_MIN = 0; 257d9822059SBarry Smith 258092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 259d9822059SBarry Smith { 260d9822059SBarry Smith PetscInt i,count = *cnt; 261d9822059SBarry Smith 262d9822059SBarry Smith PetscFunctionBegin; 2637c2de775SJed Brown if (*datatype == MPIU_REAL) { 2648c764dc5SJose Roman PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 265a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]); 2667c2de775SJed Brown } 2677c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2687c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2697c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 2707c2de775SJed Brown for (i=0; i<count; i++) { 2717c2de775SJed Brown xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 2727c2de775SJed Brown } 2737c2de775SJed Brown } 2747c2de775SJed Brown #endif 2757c2de775SJed Brown else { 2767c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 27741e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 2788c764dc5SJose Roman } 279d9822059SBarry Smith PetscFunctionReturnVoid(); 280d9822059SBarry Smith } 281d9822059SBarry Smith 282092991acSStefano Zampini PETSC_EXTERN void MPIAPI PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 283d9822059SBarry Smith { 284d9822059SBarry Smith PetscInt i,count = *cnt; 285d9822059SBarry Smith 286d9822059SBarry Smith PetscFunctionBegin; 2877c2de775SJed Brown if (*datatype == MPIU_REAL) { 2888c764dc5SJose Roman PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 289a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]); 2907c2de775SJed Brown } 2917c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2927c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2937c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 2947c2de775SJed Brown for (i=0; i<count; i++) { 2957c2de775SJed Brown xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 2967c2de775SJed Brown } 2977c2de775SJed Brown } 2987c2de775SJed Brown #endif 2997c2de775SJed Brown else { 3008c764dc5SJose Roman (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"); 30141e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 3028c764dc5SJose Roman } 303d9822059SBarry Smith PetscFunctionReturnVoid(); 304d9822059SBarry Smith } 305d9822059SBarry Smith #endif 306d9822059SBarry Smith 307480cf27aSJed Brown /* 308480cf27aSJed Brown Private routine to delete internal tag/name counter storage when a communicator is freed. 309480cf27aSJed Brown 310ff0e51ddSBarry 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. 311480cf27aSJed Brown 31212801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 313480cf27aSJed Brown 314480cf27aSJed Brown */ 31533779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state) 316480cf27aSJed Brown { 31705342407SJunchao Zhang PetscCommCounter *counter=(PetscCommCounter*)count_val; 31857f21012SBarry Smith struct PetscCommStash *comms = counter->comms, *pcomm; 319480cf27aSJed Brown 320480cf27aSJed Brown PetscFunctionBegin; 3219566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm)); 3229566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter->iflags)); 32357f21012SBarry Smith while (comms) { 3249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&comms->comm)); 32557f21012SBarry Smith pcomm = comms; 32657f21012SBarry Smith comms = comms->next; 3279566063dSJacob Faibussowitsch PetscCall(PetscFree(pcomm)); 32857f21012SBarry Smith } 3299566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter)); 330480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 331480cf27aSJed Brown } 332480cf27aSJed Brown 333480cf27aSJed Brown /* 33447435625SJed Brown This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user 335da3039f7SJed Brown calls MPI_Comm_free(). 336da3039f7SJed Brown 337da3039f7SJed Brown This is the only entry point for breaking the links between inner and outer comms. 338480cf27aSJed Brown 339ff0e51ddSBarry Smith This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator. 340480cf27aSJed Brown 34112801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 342480cf27aSJed Brown 343480cf27aSJed Brown */ 34433779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state) 345480cf27aSJed Brown { 34633779a13SJunchao Zhang union {MPI_Comm comm; void *ptr;} icomm; 347480cf27aSJed Brown 348480cf27aSJed Brown PetscFunctionBegin; 34912801b39SBarry Smith if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval"); 350ec4fadc2SJed Brown icomm.ptr = attr_val; 35176bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 35233779a13SJunchao Zhang /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */ 35333779a13SJunchao Zhang PetscMPIInt flg; 35433779a13SJunchao Zhang union {MPI_Comm comm; void *ptr;} ocomm; 3559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg)); 35633779a13SJunchao Zhang if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute"); 35733779a13SJunchao 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"); 35833779a13SJunchao Zhang } 3599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval)); 3609566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm)); 361da3039f7SJed Brown PetscFunctionReturn(MPI_SUCCESS); 362b89831e5SBarry Smith } 363da3039f7SJed Brown 364da3039f7SJed Brown /* 36533779a13SJunchao 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. 366da3039f7SJed Brown */ 36733779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state) 368da3039f7SJed Brown { 369da3039f7SJed Brown PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm)); 371480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 372480cf27aSJed Brown } 373480cf27aSJed Brown 37433779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *); 37542218b76SBarry Smith 376951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 3778cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*); 3788cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*); 3798cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*); 380e39fd77fSBarry Smith #endif 381e39fd77fSBarry Smith 3820f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE; 38312801b39SBarry Smith 384eb27c7c8SSatish Balay PETSC_INTERN int PetscGlobalArgc; 385eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs; 3866ae9a8a6SBarry Smith int PetscGlobalArgc = 0; 38702c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL; 388dff31646SBarry Smith PetscSegBuffer PetscCitationsList; 389e5c89e4eSSatish Balay 390dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void) 391051e4cf2SJed Brown { 392051e4cf2SJed Brown PetscFunctionBegin; 3939566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(1,10000,&PetscCitationsList)); 39430c35bf2SSatish Balay 39530c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\ 39630c35bf2SSatish Balay Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\ 39730c35bf2SSatish Balay and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\ 39830c35bf2SSatish Balay and Victor Eijkhout and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\ 39930c35bf2SSatish Balay and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\ 40030c35bf2SSatish Balay and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\ 40130c35bf2SSatish Balay and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\ 40230c35bf2SSatish Balay and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\ 40330c35bf2SSatish Balay Title = {{PETSc/TAO} Users Manual},\n\ 40430c35bf2SSatish Balay Number = {ANL-21/39 - Revision 3.17},\n\ 40530c35bf2SSatish Balay Institution = {Argonne National Laboratory},\n\ 40630c35bf2SSatish Balay Year = {2022}\n}\n",NULL)); 40730c35bf2SSatish Balay 40830c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\ 40930c35bf2SSatish Balay Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\ 41030c35bf2SSatish Balay Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\ 41130c35bf2SSatish Balay Booktitle = {Modern Software Tools in Scientific Computing},\n\ 41230c35bf2SSatish Balay Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\ 41330c35bf2SSatish Balay Pages = {163--202},\n\ 41430c35bf2SSatish Balay Publisher = {Birkh{\\\"{a}}user Press},\n\ 41530c35bf2SSatish Balay Year = {1997}\n}\n",NULL)); 41630c35bf2SSatish Balay 417051e4cf2SJed Brown PetscFunctionReturn(0); 418051e4cf2SJed Brown } 419e5c89e4eSSatish Balay 4202d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */ 4212d747510SLisandro Dalcin 4222d747510SLisandro Dalcin PetscErrorCode PetscSetProgramName(const char name[]) 4232d747510SLisandro Dalcin { 4242d747510SLisandro Dalcin PetscFunctionBegin; 4259566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(programname,name,sizeof(programname))); 4262d747510SLisandro Dalcin PetscFunctionReturn(0); 4272d747510SLisandro Dalcin } 4282d747510SLisandro Dalcin 4292d747510SLisandro Dalcin /*@C 4302d747510SLisandro Dalcin PetscGetProgramName - Gets the name of the running program. 4312d747510SLisandro Dalcin 4322d747510SLisandro Dalcin Not Collective 4332d747510SLisandro Dalcin 4342d747510SLisandro Dalcin Input Parameter: 4352d747510SLisandro Dalcin . len - length of the string name 4362d747510SLisandro Dalcin 4372d747510SLisandro Dalcin Output Parameter: 4382d747510SLisandro Dalcin . name - the name of the running program 4392d747510SLisandro Dalcin 4402d747510SLisandro Dalcin Level: advanced 4412d747510SLisandro Dalcin 4422d747510SLisandro Dalcin Notes: 4432d747510SLisandro Dalcin The name of the program is copied into the user-provided character 4442d747510SLisandro Dalcin array of length len. On some machines the program name includes 4452d747510SLisandro Dalcin its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN. 4462d747510SLisandro Dalcin @*/ 4472d747510SLisandro Dalcin PetscErrorCode PetscGetProgramName(char name[],size_t len) 4482d747510SLisandro Dalcin { 4492d747510SLisandro Dalcin PetscFunctionBegin; 4509566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(name,programname,len)); 4512d747510SLisandro Dalcin PetscFunctionReturn(0); 4522d747510SLisandro Dalcin } 4532d747510SLisandro Dalcin 454e5c89e4eSSatish Balay /*@C 455e5c89e4eSSatish Balay PetscGetArgs - Allows you to access the raw command line arguments anywhere 456e5c89e4eSSatish Balay after PetscInitialize() is called but before PetscFinalize(). 457e5c89e4eSSatish Balay 458e5c89e4eSSatish Balay Not Collective 459e5c89e4eSSatish Balay 460e5c89e4eSSatish Balay Output Parameters: 461e5c89e4eSSatish Balay + argc - count of number of command line arguments 462e5c89e4eSSatish Balay - args - the command line arguments 463e5c89e4eSSatish Balay 464e5c89e4eSSatish Balay Level: intermediate 465e5c89e4eSSatish Balay 466e5c89e4eSSatish Balay Notes: 467e5c89e4eSSatish Balay This is usually used to pass the command line arguments into other libraries 468e5c89e4eSSatish Balay that are called internally deep in PETSc or the application. 469e5c89e4eSSatish Balay 470f177e3b1SBarry Smith The first argument contains the program name as is normal for C arguments. 471f177e3b1SBarry Smith 472793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments() 473e5c89e4eSSatish Balay 474e5c89e4eSSatish Balay @*/ 4757087cfbeSBarry Smith PetscErrorCode PetscGetArgs(int *argc,char ***args) 476e5c89e4eSSatish Balay { 477e5c89e4eSSatish Balay PetscFunctionBegin; 47839a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()"); 479e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 480e5c89e4eSSatish Balay *args = PetscGlobalArgs; 481e5c89e4eSSatish Balay PetscFunctionReturn(0); 482e5c89e4eSSatish Balay } 483e5c89e4eSSatish Balay 484793721a6SBarry Smith /*@C 485793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 486793721a6SBarry Smith after PetscInitialize() is called but before PetscFinalize(). 487793721a6SBarry Smith 488793721a6SBarry Smith Not Collective 489793721a6SBarry Smith 490793721a6SBarry Smith Output Parameters: 491793721a6SBarry Smith . args - the command line arguments 492793721a6SBarry Smith 493793721a6SBarry Smith Level: intermediate 494793721a6SBarry Smith 495793721a6SBarry Smith Notes: 496793721a6SBarry Smith This does NOT start with the program name and IS null terminated (final arg is void) 497793721a6SBarry Smith 498793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments() 499793721a6SBarry Smith 500793721a6SBarry Smith @*/ 5017087cfbeSBarry Smith PetscErrorCode PetscGetArguments(char ***args) 502793721a6SBarry Smith { 503793721a6SBarry Smith PetscInt i,argc = PetscGlobalArgc; 504793721a6SBarry Smith 505793721a6SBarry Smith PetscFunctionBegin; 50639a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled,PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()"); 5072d747510SLisandro Dalcin if (!argc) {*args = NULL; PetscFunctionReturn(0);} 5089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(argc,args)); 5099566063dSJacob Faibussowitsch for (i=0; i<argc-1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i])); 5102d747510SLisandro Dalcin (*args)[argc-1] = NULL; 511793721a6SBarry Smith PetscFunctionReturn(0); 512793721a6SBarry Smith } 513793721a6SBarry Smith 514793721a6SBarry Smith /*@C 515793721a6SBarry Smith PetscFreeArguments - Frees the memory obtained with PetscGetArguments() 516793721a6SBarry Smith 517793721a6SBarry Smith Not Collective 518793721a6SBarry Smith 519793721a6SBarry Smith Output Parameters: 520793721a6SBarry Smith . args - the command line arguments 521793721a6SBarry Smith 522793721a6SBarry Smith Level: intermediate 523793721a6SBarry Smith 524793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments() 525793721a6SBarry Smith 526793721a6SBarry Smith @*/ 5277087cfbeSBarry Smith PetscErrorCode PetscFreeArguments(char **args) 528793721a6SBarry Smith { 52939a651e2SJacob Faibussowitsch PetscFunctionBegin; 53039a651e2SJacob Faibussowitsch if (args) { 531793721a6SBarry Smith PetscInt i = 0; 532793721a6SBarry Smith 5339566063dSJacob Faibussowitsch while (args[i]) PetscCall(PetscFree(args[i++])); 5349566063dSJacob Faibussowitsch PetscCall(PetscFree(args)); 53539a651e2SJacob Faibussowitsch } 536793721a6SBarry Smith PetscFunctionReturn(0); 537793721a6SBarry Smith } 538793721a6SBarry Smith 53927104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 54030befbd2SBarry Smith #include <petscconfiginfo.h> 54130befbd2SBarry Smith 54295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) 54311525c0dSBarry Smith { 54427104ee2SJacob Faibussowitsch PetscFunctionBegin; 54511525c0dSBarry Smith if (!PetscGlobalRank) { 54630befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64]; 54711525c0dSBarry Smith int port; 548ffbd1cfbSBarry Smith PetscBool flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE; 54911525c0dSBarry Smith size_t applinelen,introlen; 550ffbd1cfbSBarry Smith char sawsurl[256]; 55111525c0dSBarry Smith 5529566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-saws_log",&flg)); 55311525c0dSBarry Smith if (flg) { 55411525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 55511525c0dSBarry Smith 5569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL)); 55711525c0dSBarry Smith if (sawslog[0]) { 55811525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog)); 55911525c0dSBarry Smith } else { 56011525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL)); 56111525c0dSBarry Smith } 56211525c0dSBarry Smith } 5639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg)); 56411525c0dSBarry Smith if (flg) { 56511525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert)); 56611525c0dSBarry Smith } 5679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL)); 568ffbd1cfbSBarry Smith if (selectport) { 569ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Get_Available_Port,(&port)); 570ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 571ffbd1cfbSBarry Smith } else { 5729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg)); 57311525c0dSBarry Smith if (flg) { 57411525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 57511525c0dSBarry Smith } 576ffbd1cfbSBarry Smith } 5779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg)); 57811525c0dSBarry Smith if (flg) { 5792f613bf5SBarry Smith PetscStackCallSAWs(SAWs_Set_Document_Root,(root)); 5809566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(root,".",&rootlocal)); 5819c1e0ce8SBarry Smith } else { 5829566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-saws_options",&flg)); 5839c1e0ce8SBarry Smith if (flg) { 5849566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root))); 5852f613bf5SBarry Smith PetscStackCallSAWs(SAWs_Set_Document_Root,(root)); 5869c1e0ce8SBarry Smith } 58711525c0dSBarry Smith } 5889566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2)); 58911525c0dSBarry Smith if (flg2) { 59011525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 59128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option"); 5929566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root)); 5939566063dSJacob Faibussowitsch PetscCall(PetscTestDirectory(jsdir,'r',&flg)); 59428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory"); 5952f613bf5SBarry Smith PetscStackCallSAWs(SAWs_Push_Local_Header,()); 59611525c0dSBarry Smith } 5979566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(programname,sizeof(programname))); 5989566063dSJacob Faibussowitsch PetscCall(PetscStrlen(help,&applinelen)); 59911525c0dSBarry Smith introlen = 4096 + applinelen; 60030a8c9c0SSurtai Han applinelen += 1024; 6019566063dSJacob Faibussowitsch PetscCall(PetscMalloc(applinelen,&appline)); 6029566063dSJacob Faibussowitsch PetscCall(PetscMalloc(introlen,&intro)); 60311525c0dSBarry Smith 60411525c0dSBarry Smith if (rootlocal) { 6059566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline,applinelen,"%s.c.html",programname)); 6069566063dSJacob Faibussowitsch PetscCall(PetscTestFile(appline,'r',&rootlocal)); 60711525c0dSBarry Smith } 6089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetAll(NULL,&options)); 60911525c0dSBarry Smith if (rootlocal && help) { 6109566063dSJacob Faibussowitsch PetscCall(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)); 61111525c0dSBarry Smith } else if (help) { 6129566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help)); 61311525c0dSBarry Smith } else { 6149566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options)); 61511525c0dSBarry Smith } 6169566063dSJacob Faibussowitsch PetscCall(PetscFree(options)); 6179566063dSJacob Faibussowitsch PetscCall(PetscGetVersion(version,sizeof(version))); 6189566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(intro,introlen,"<body>\n" 619a17b96a8SKyle Gerard Felker "<center><h2> <a href=\"https://petsc.org/\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n" 620df62c222SBarry 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" 6215f80ce2aSJacob Faibussowitsch "%s",version,petscconfigureoptions,appline)); 62243da4ab2SBarry Smith PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro)); 6239566063dSJacob Faibussowitsch PetscCall(PetscFree(intro)); 6249566063dSJacob Faibussowitsch PetscCall(PetscFree(appline)); 625ffbd1cfbSBarry Smith if (selectport) { 626aa573868SBarry Smith PetscBool silent; 6277d812c46SBarry Smith 6287d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 62939a651e2SJacob Faibussowitsch while (SAWs_Initialize()) { 6307d812c46SBarry Smith PetscStackCallSAWs(SAWs_Get_Available_Port,(&port)); 6317d812c46SBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 6327d812c46SBarry Smith } 6337d812c46SBarry Smith 6349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL)); 635aa573868SBarry Smith if (!silent) { 636ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl)); 6379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl)); 638ffbd1cfbSBarry Smith } 6397d812c46SBarry Smith } else { 6407d812c46SBarry Smith PetscStackCallSAWs(SAWs_Initialize,()); 641aa573868SBarry Smith } 6429566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister("@TechReport{ saws,\n" 6430af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 6440af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 6450af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 6465f80ce2aSJacob Faibussowitsch " Year = 2013\n}\n",NULL)); 64711525c0dSBarry Smith } 64811525c0dSBarry Smith PetscFunctionReturn(0); 64911525c0dSBarry Smith } 65011525c0dSBarry Smith #endif 65111525c0dSBarry Smith 6524dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 6534dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) 6544dfee713SSatish Balay { 6554dfee713SSatish Balay PetscFunctionBegin; 6564dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 6574dfee713SSatish Balay /* see MPI.py for details on this bug */ 6584dfee713SSatish Balay (void) setenv("HWLOC_COMPONENTS","-x86",1); 6594dfee713SSatish Balay #endif 6604dfee713SSatish Balay PetscFunctionReturn(0); 6614dfee713SSatish Balay } 6624dfee713SSatish Balay 663a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS) 664a56f64adSBarry Smith #include <adios.h> 66522580e64SBarry Smith #include <adios_read.h> 6667b56e58cSSatish Balay int64_t Petsc_adios_group; 667a56f64adSBarry Smith #endif 668a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP) 669cd1b99a6SBarry Smith #include <omp.h> 670f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 671cd1b99a6SBarry Smith #endif 672a56f64adSBarry Smith 673a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_DEVICE) 674a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h> 675a4af0ceeSJacob Faibussowitsch # if PetscDefined(HAVE_CUDA) 676a4af0ceeSJacob Faibussowitsch // REMOVE ME 677a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL; 678a4af0ceeSJacob Faibussowitsch # endif 679a4af0ceeSJacob Faibussowitsch # if PetscDefined(HAVE_HIP) 680a4af0ceeSJacob Faibussowitsch // REMOVE ME 681a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL; 682a4af0ceeSJacob Faibussowitsch # endif 683a4af0ceeSJacob Faibussowitsch #endif 684a4af0ceeSJacob Faibussowitsch 68527104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H) 686bc8a24ecSBarry Smith #include <dlfcn.h> 687bc8a24ecSBarry Smith #endif 688a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 689a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void); 690a4af0ceeSJacob Faibussowitsch #endif 691a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 692a4af0ceeSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViennaCLInit(); 693a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE; 694a4af0ceeSJacob Faibussowitsch #endif 695bc8a24ecSBarry Smith 69685649d77SJunchao Zhang /* 69785649d77SJunchao Zhang PetscInitialize_Common - shared code between C and Fortran initialization 69885649d77SJunchao Zhang 69985649d77SJunchao Zhang prog: program name 70002101c96SSatish Balay file: optional PETSc database file name. Might be in Fortran string format when 'ftn' is true 70185649d77SJunchao Zhang help: program help message 70202101c96SSatish Balay ftn: is it called from Fortran initilization (petscinitializef_)? 70385649d77SJunchao Zhang readarguments,len: used when fortran is true 70485649d77SJunchao Zhang */ 70502101c96SSatish Balay PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char* prog,const char* file,const char *help,PetscBool ftn,PetscBool readarguments,PetscInt len) 70685649d77SJunchao Zhang { 70785649d77SJunchao Zhang PetscMPIInt size; 70885649d77SJunchao Zhang PetscBool flg = PETSC_TRUE; 70985649d77SJunchao Zhang char hostname[256]; 71085649d77SJunchao Zhang 71127104ee2SJacob Faibussowitsch PetscFunctionBegin; 71227104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(0); 71339a651e2SJacob Faibussowitsch /* these must be initialized in a routine, not as a constant declaration */ 71439a651e2SJacob Faibussowitsch PETSC_STDOUT = stdout; 71539a651e2SJacob Faibussowitsch PETSC_STDERR = stderr; 71639a651e2SJacob Faibussowitsch 7179566063dSJacob Faibussowitsch /* PetscCall can be used from now */ 71839a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_TRUE; 71939a651e2SJacob Faibussowitsch 72085649d77SJunchao Zhang /* 72185649d77SJunchao Zhang The checking over compatible runtime libraries is complicated by the MPI ABI initiative 72285649d77SJunchao Zhang https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 72385649d77SJunchao Zhang MPICH v3.1 (Released February 2014) 72485649d77SJunchao Zhang IBM MPI v2.1 (December 2014) 72585649d77SJunchao Zhang Intel MPI Library v5.0 (2014) 72685649d77SJunchao Zhang Cray MPT v7.0.0 (June 2014) 72785649d77SJunchao Zhang As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 72885649d77SJunchao Zhang listed above and since that time are compatible. 72985649d77SJunchao Zhang 73085649d77SJunchao Zhang Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 73185649d77SJunchao Zhang at compile time or runtime. Thus we will need to systematically track the allowed versions 73285649d77SJunchao Zhang and how they are represented in the mpi.h and MPI_Get_library_version() output in order 73385649d77SJunchao Zhang to perform the checking. 73485649d77SJunchao Zhang 73585649d77SJunchao Zhang Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 73685649d77SJunchao Zhang 73785649d77SJunchao Zhang Questions: 73885649d77SJunchao Zhang 73985649d77SJunchao Zhang Should the checks for ABI incompatibility be only on the major version number below? 74085649d77SJunchao Zhang Presumably the output to stderr will be removed before a release. 74185649d77SJunchao Zhang */ 74285649d77SJunchao Zhang 74385649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 74485649d77SJunchao Zhang { 74585649d77SJunchao Zhang char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 74685649d77SJunchao Zhang PetscMPIInt mpilibraryversionlength; 74739a651e2SJacob Faibussowitsch 7489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength)); 74985649d77SJunchao Zhang /* check for MPICH versions before MPI ABI initiative */ 75085649d77SJunchao Zhang #if defined(MPICH_VERSION) 75185649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000 75285649d77SJunchao Zhang { 75385649d77SJunchao Zhang char *ver,*lf; 75485649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 75539a651e2SJacob Faibussowitsch 7569566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion,"MPICH Version:",&ver)); 75739a651e2SJacob Faibussowitsch if (ver) { 7589566063dSJacob Faibussowitsch PetscCall(PetscStrchr(ver,'\n',&lf)); 75939a651e2SJacob Faibussowitsch if (lf) { 76085649d77SJunchao Zhang *lf = 0; 7619566063dSJacob Faibussowitsch PetscCall(PetscStrendswith(ver,MPICH_VERSION,&flg)); 76285649d77SJunchao Zhang } 76385649d77SJunchao Zhang } 76485649d77SJunchao Zhang if (!flg) { 7659566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n",mpilibraryversion,MPICH_VESION)); 76685649d77SJunchao Zhang flg = PETSC_TRUE; 76785649d77SJunchao Zhang } 76885649d77SJunchao Zhang } 76985649d77SJunchao Zhang #endif 77085649d77SJunchao Zhang /* check for OpenMPI version, it is not part of the MPI ABI initiative (is it part of another initiative that needs to be handled?) */ 77185649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION) 77285649d77SJunchao Zhang { 77385649d77SJunchao Zhang char *ver,bs[MPI_MAX_LIBRARY_VERSION_STRING],*bsf; 77485649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 77585649d77SJunchao Zhang #define PSTRSZ 2 77685649d77SJunchao Zhang char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI","FUJITSU MPI"}; 77785649d77SJunchao Zhang char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v","Library "}; 77885649d77SJunchao Zhang int i; 77985649d77SJunchao Zhang for (i=0; i<PSTRSZ; i++) { 7809566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion,ompistr1[i],&ver)); 78139a651e2SJacob Faibussowitsch if (ver) { 7829566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(bs,MPI_MAX_LIBRARY_VERSION_STRING,"%s%d.%d",ompistr2[i],OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION)); 7839566063dSJacob Faibussowitsch PetscCall(PetscStrstr(ver,bs,&bsf)); 78439a651e2SJacob Faibussowitsch if (bsf) flg = PETSC_TRUE; 78585649d77SJunchao Zhang break; 78685649d77SJunchao Zhang } 78785649d77SJunchao Zhang } 78885649d77SJunchao Zhang if (!flg) { 7897d3de750SJacob Faibussowitsch PetscInfo(NULL,"PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n",mpilibraryversion,OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION); 79085649d77SJunchao Zhang flg = PETSC_TRUE; 79185649d77SJunchao Zhang } 79285649d77SJunchao Zhang } 79385649d77SJunchao Zhang #endif 79485649d77SJunchao Zhang } 79585649d77SJunchao Zhang #endif 79685649d77SJunchao Zhang 797*6f5d4113SSatish Balay #if PetscDefined(HAVE_DLSYM) && defined(__USE_GNU) 79885649d77SJunchao Zhang /* 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 */ 79939a651e2SJacob Faibussowitsch PetscCheck(!dlsym(RTLD_DEFAULT,"ompi_mpi_init") || !dlsym(RTLD_DEFAULT,"MPID_Abort"),PETSC_COMM_SELF,PETSC_ERR_MPI_LIB_INCOMP,"Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly"); 80085649d77SJunchao Zhang #endif 80185649d77SJunchao Zhang 80285649d77SJunchao Zhang /* on Windows - set printf to default to printing 2 digit exponents */ 80385649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT) 80485649d77SJunchao Zhang _set_output_format(_TWO_DIGIT_EXPONENT); 80585649d77SJunchao Zhang #endif 80685649d77SJunchao Zhang 8079566063dSJacob Faibussowitsch PetscCall(PetscOptionsCreateDefault()); 80885649d77SJunchao Zhang 80985649d77SJunchao Zhang PetscFinalizeCalled = PETSC_FALSE; 81085649d77SJunchao Zhang 8119566063dSJacob Faibussowitsch PetscCall(PetscSetProgramName(prog)); 8129566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen)); 8139566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout)); 8149566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr)); 8159566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscCommSpinLock)); 81685649d77SJunchao Zhang 81785649d77SJunchao Zhang if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 8189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN)); 81985649d77SJunchao Zhang 82085649d77SJunchao Zhang if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) { 8219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS)); 8229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE)); 82385649d77SJunchao Zhang } 82485649d77SJunchao Zhang 82585649d77SJunchao Zhang /* Done after init due to a bug in MPICH-GM? */ 8269566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 82785649d77SJunchao Zhang 8289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank)); 8299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize)); 83085649d77SJunchao Zhang 83185649d77SJunchao Zhang MPIU_BOOL = MPI_INT; 83285649d77SJunchao Zhang MPIU_ENUM = MPI_INT; 83385649d77SJunchao Zhang MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64; 83485649d77SJunchao Zhang if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 83585649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 83685649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG) 83785649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 83885649d77SJunchao Zhang #endif 83939a651e2SJacob Faibussowitsch else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP_SYS,"Could not find MPI type for size_t"); 84085649d77SJunchao Zhang 84185649d77SJunchao Zhang /* 84285649d77SJunchao Zhang Initialized the global complex variable; this is because with 84385649d77SJunchao Zhang shared libraries the constructors for global variables 84485649d77SJunchao Zhang are not called; at least on IRIX. 84585649d77SJunchao Zhang */ 84685649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX) 84785649d77SJunchao Zhang { 84885649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 84985649d77SJunchao Zhang PetscComplex ic(0.0,1.0); 85085649d77SJunchao Zhang PETSC_i = ic; 85185649d77SJunchao Zhang #else 85285649d77SJunchao Zhang PETSC_i = _Complex_I; 85385649d77SJunchao Zhang #endif 85485649d77SJunchao Zhang } 85585649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */ 85685649d77SJunchao Zhang 85785649d77SJunchao Zhang /* 85885649d77SJunchao Zhang Create the PETSc MPI reduction operator that sums of the first 85985649d77SJunchao Zhang half of the entries and maxes the second half. 86085649d77SJunchao Zhang */ 8619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP)); 86285649d77SJunchao Zhang 86385649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) 8649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128)); 8659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128)); 86685649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX) 8679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128)); 8689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128)); 86985649d77SJunchao Zhang #endif 8709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscMax_Local,1,&MPIU_MAX)); 8719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscMin_Local,1,&MPIU_MIN)); 87285649d77SJunchao Zhang #elif defined(PETSC_USE_REAL___FP16) 8739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16)); 8749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FP16)); 8759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscMax_Local,1,&MPIU_MAX)); 8769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscMin_Local,1,&MPIU_MIN)); 87785649d77SJunchao Zhang #endif 87885649d77SJunchao Zhang 87985649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 8809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscSum_Local,1,&MPIU_SUM)); 88185649d77SJunchao Zhang #endif 88285649d77SJunchao Zhang 8839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR)); 8849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR)); 88585649d77SJunchao Zhang 88685649d77SJunchao Zhang /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */ 88785649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI) 88885649d77SJunchao Zhang { 88985649d77SJunchao Zhang struct PetscRealInt { PetscReal v; PetscInt i; }; 89085649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1,1}; 89185649d77SJunchao Zhang MPI_Aint blockOffsets[2] = {offsetof(struct PetscRealInt,v),offsetof(struct PetscRealInt,i)}; 89285649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_REAL,MPIU_INT}, tmpStruct; 89385649d77SJunchao Zhang 8949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct)); 8959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscRealInt),&MPIU_REAL_INT)); 8969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 8979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT)); 89885649d77SJunchao Zhang } 89985649d77SJunchao Zhang { 90085649d77SJunchao Zhang struct PetscScalarInt { PetscScalar v; PetscInt i; }; 90185649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1,1}; 90285649d77SJunchao Zhang MPI_Aint blockOffsets[2] = {offsetof(struct PetscScalarInt,v),offsetof(struct PetscScalarInt,i)}; 90385649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_SCALAR,MPIU_INT}, tmpStruct; 90485649d77SJunchao Zhang 9059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct)); 9069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscScalarInt),&MPIU_SCALAR_INT)); 9079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT)); 90985649d77SJunchao Zhang } 91085649d77SJunchao Zhang #endif 91185649d77SJunchao Zhang 91285649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES) 9139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT)); 9149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2INT)); 91585649d77SJunchao Zhang #endif 9169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4,MPI_INT,&MPI_4INT)); 9179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPI_4INT)); 9189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4,MPIU_INT,&MPIU_4INT)); 9199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_4INT)); 92085649d77SJunchao Zhang 92185649d77SJunchao Zhang /* 92285649d77SJunchao Zhang Attributes to be set on PETSc communicators 92385649d77SJunchao Zhang */ 9249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0)); 9259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0)); 9269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0)); 9279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0)); 92885649d77SJunchao Zhang 92985649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN) 9309566063dSJacob Faibussowitsch if (ftn) PetscCall(PetscInitFortran_Private(readarguments,file,len)); 93185649d77SJunchao Zhang else 93285649d77SJunchao Zhang #endif 9339566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL,&PetscGlobalArgc,&PetscGlobalArgs,file)); 93485649d77SJunchao Zhang 93585649d77SJunchao Zhang /* call a second time so it can look in the options database */ 9369566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 93785649d77SJunchao Zhang 93885649d77SJunchao Zhang /* 93985649d77SJunchao Zhang Check system options and print help 94085649d77SJunchao Zhang */ 9419566063dSJacob Faibussowitsch PetscCall(PetscOptionsCheckInitial_Private(help)); 94285649d77SJunchao Zhang 943a4af0ceeSJacob Faibussowitsch /* 944a4af0ceeSJacob Faibussowitsch Initialize PetscDevice and PetscDeviceContext 945a4af0ceeSJacob Faibussowitsch 946a4af0ceeSJacob Faibussowitsch Note to any future devs thinking of moving this, proper initialization requires: 947a4af0ceeSJacob Faibussowitsch 1. MPI initialized 948a4af0ceeSJacob Faibussowitsch 2. Options DB initialized 949a4af0ceeSJacob Faibussowitsch 3. Petsc error handling initialized, specifically signal handlers. This expects to set up its own SIGSEV handler via 950a4af0ceeSJacob Faibussowitsch the push/pop interface. 951a4af0ceeSJacob Faibussowitsch */ 952a2158755SJunchao Zhang #if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP) || PetscDefined(HAVE_SYCL)) 9539566063dSJacob Faibussowitsch PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD)); 954a4af0ceeSJacob Faibussowitsch #endif 955a4af0ceeSJacob Faibussowitsch 956a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 957a4af0ceeSJacob Faibussowitsch flg = PETSC_FALSE; 9589566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-log_summary",&flg)); 9599566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsHasName(NULL,NULL,"-log_view",&flg)); 9609566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsGetBool(NULL,NULL,"-viennacl_synchronize",&flg,NULL)); 961a4af0ceeSJacob Faibussowitsch PetscViennaCLSynchronize = flg; 9629566063dSJacob Faibussowitsch PetscCall(PetscViennaCLInit()); 963a4af0ceeSJacob Faibussowitsch #endif 964a4af0ceeSJacob Faibussowitsch 965a4af0ceeSJacob Faibussowitsch /* 966a4af0ceeSJacob Faibussowitsch Creates the logging data structures; this is enabled even if logging is not turned on 967a4af0ceeSJacob Faibussowitsch This is the last thing we do before returning to the user code to prevent having the 968a4af0ceeSJacob Faibussowitsch logging numbers contaminated by any startup time associated with MPI 969a4af0ceeSJacob Faibussowitsch */ 970a4af0ceeSJacob Faibussowitsch #if defined(PETSC_USE_LOG) 9719566063dSJacob Faibussowitsch PetscCall(PetscLogInitialize()); 972a4af0ceeSJacob Faibussowitsch #endif 973a4af0ceeSJacob Faibussowitsch 9749566063dSJacob Faibussowitsch PetscCall(PetscCitationsInitialize()); 97585649d77SJunchao Zhang 97685649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS) 9779566063dSJacob Faibussowitsch PetscCall(PetscInitializeSAWs(ftn ? NULL : help)); 97827104ee2SJacob Faibussowitsch flg = PETSC_FALSE; 9799566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-stack_view",&flg)); 9809566063dSJacob Faibussowitsch if (flg) PetscCall(PetscStackViewSAWs()); 98185649d77SJunchao Zhang #endif 98285649d77SJunchao Zhang 98385649d77SJunchao Zhang /* 98485649d77SJunchao Zhang Load the dynamic libraries (on machines that support them), this registers all 98585649d77SJunchao Zhang the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 98685649d77SJunchao Zhang */ 9879566063dSJacob Faibussowitsch PetscCall(PetscInitialize_DynamicLibraries()); 98885649d77SJunchao Zhang 9899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 9909566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"PETSc successfully started: number of processors = %d\n",size)); 9919566063dSJacob Faibussowitsch PetscCall(PetscGetHostName(hostname,256)); 9929566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Running on machine: %s\n",hostname)); 99385649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP) 99485649d77SJunchao Zhang { 99585649d77SJunchao Zhang PetscBool omp_view_flag; 99685649d77SJunchao Zhang char *threads = getenv("OMP_NUM_THREADS"); 9975f80ce2aSJacob Faibussowitsch PetscErrorCode ierr; 99885649d77SJunchao Zhang 99985649d77SJunchao Zhang if (threads) { 10009566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n",threads)); 100185649d77SJunchao Zhang (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads); 100285649d77SJunchao Zhang } else { 10032f840973SStefano Zampini PetscNumOMPThreads = (PetscInt) omp_get_max_threads(); 10049566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n",PetscNumOMPThreads)); 100585649d77SJunchao Zhang } 10069566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");PetscCall(ierr); 10079566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-omp_num_threads","Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS","None",PetscNumOMPThreads,&PetscNumOMPThreads,&flg)); 10089566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag)); 10099566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 101085649d77SJunchao Zhang if (flg) { 10119566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n",PetscNumOMPThreads)); 101285649d77SJunchao Zhang omp_set_num_threads((int)PetscNumOMPThreads); 101385649d77SJunchao Zhang } 101485649d77SJunchao Zhang if (omp_view_flag) { 10159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %" PetscInt_FMT "\n",PetscNumOMPThreads)); 101685649d77SJunchao Zhang } 101785649d77SJunchao Zhang } 101885649d77SJunchao Zhang #endif 101985649d77SJunchao Zhang 102085649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 102185649d77SJunchao Zhang /* 102285649d77SJunchao Zhang Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 102385649d77SJunchao Zhang 102485649d77SJunchao Zhang Currently not used because it is not supported by MPICH. 102585649d77SJunchao Zhang */ 10269566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL)); 102785649d77SJunchao Zhang #endif 102885649d77SJunchao Zhang 102985649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS) 10309566063dSJacob Faibussowitsch PetscCall(PetscFPTCreate(10000)); 103185649d77SJunchao Zhang #endif 103285649d77SJunchao Zhang 103385649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC) 103485649d77SJunchao Zhang { 103585649d77SJunchao Zhang PetscViewer viewer; 10369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg)); 103785649d77SJunchao Zhang if (flg) { 10389566063dSJacob Faibussowitsch PetscCall(PetscProcessPlacementView(viewer)); 10399566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 104085649d77SJunchao Zhang } 104185649d77SJunchao Zhang } 104285649d77SJunchao Zhang #endif 104385649d77SJunchao Zhang 104485649d77SJunchao Zhang flg = PETSC_TRUE; 10459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL)); 10469566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE)); 104785649d77SJunchao Zhang 104885649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS) 10499566063dSJacob Faibussowitsch PetscCall(adios_init_noxml(PETSC_COMM_WORLD)); 10509566063dSJacob Faibussowitsch PetscCall(adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default)); 10519566063dSJacob Faibussowitsch PetscCall(adios_select_method(Petsc_adios_group,"MPI","","")); 10529566063dSJacob Faibussowitsch PetscCall(adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"")); 105385649d77SJunchao Zhang #endif 105485649d77SJunchao Zhang 105585649d77SJunchao Zhang #if defined(__VALGRIND_H) 105685649d77SJunchao Zhang PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND? PETSC_TRUE: PETSC_FALSE; 105785649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE) 10589566063dSJacob Faibussowitsch if (PETSC_RUNNING_ON_VALGRIND) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING: Running valgrind with the MacOS native BLAS and LAPACK can fail. If it fails suggest configuring with --download-fblaslapack or --download-f2cblaslapack")); 105985649d77SJunchao Zhang #endif 106085649d77SJunchao Zhang #endif 106185649d77SJunchao Zhang /* 106285649d77SJunchao Zhang Set flag that we are completely initialized 106385649d77SJunchao Zhang */ 106485649d77SJunchao Zhang PetscInitializeCalled = PETSC_TRUE; 106585649d77SJunchao Zhang 10669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-python",&flg)); 10679566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPythonInitialize(NULL,NULL)); 106827104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 106985649d77SJunchao Zhang } 107085649d77SJunchao Zhang 1071e5c89e4eSSatish Balay /*@C 1072e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 1073e5c89e4eSSatish Balay PetscInitialize() calls MPI_Init() if that has yet to be called, 1074e5c89e4eSSatish Balay so this routine should always be called near the beginning of 1075e5c89e4eSSatish Balay your program -- usually the very first line! 1076e5c89e4eSSatish Balay 1077e5c89e4eSSatish Balay Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set 1078e5c89e4eSSatish Balay 1079e5c89e4eSSatish Balay Input Parameters: 1080e5c89e4eSSatish Balay + argc - count of number of command line arguments 1081e5c89e4eSSatish Balay . args - the command line arguments 1082be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format. 1083be10d61cSLisandro Dalcin Use NULL or empty string to not check for code specific file. 1084be10d61cSLisandro Dalcin Also checks ~/.petscrc, .petscrc and petscrc. 1085c5b5d8d5SVaclav Hapla Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files. 10860298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message 1087e5c89e4eSSatish Balay 108805827820SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that 108905827820SBarry Smith communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a 109005827820SBarry Smith four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not, 109105827820SBarry Smith then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even 109205827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 1093e5c89e4eSSatish Balay 1094e5c89e4eSSatish Balay Options Database Keys: 10957ca660e7SBarry Smith + -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message 10967ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 1097e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 10988a690491SBarry Smith . -on_error_emacs <machinename> - causes emacsclient to jump to error file 10998a690491SBarry Smith . -on_error_abort - calls abort() when error detected (no traceback) 11008a690491SBarry Smith . -on_error_mpiabort - calls MPI_abort() when error detected 11018a690491SBarry Smith . -error_output_stderr - prints error messages to stderr instead of the default stdout 11028a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 1103bf4d2887SBarry Smith . -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger 1104e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds) - Pauses debugger 1105e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 1106e5c89e4eSSatish Balay process and wait (-debugger_pause) seconds for attachment 110779dccf82SBarry Smith . -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug) 110879dccf82SBarry Smith . -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no) 110979dccf82SBarry Smith . -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug() 1110aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 111192f119d6SBarry 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 111292f119d6SBarry Smith . -malloc_view - show a list of all allocated memory during PetscFinalize() 111392f119d6SBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view 1114608c71bfSMatthew G. Knepley . -malloc_requested_size - malloc logging will record the requested size rather than size after alignment 111592f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 1116e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 1117e5c89e4eSSatish Balay . -shared_tmp - indicates /tmp directory is shared by all processors 1118e5c89e4eSSatish Balay . -not_shared_tmp - each processor has own /tmp 1119e5c89e4eSSatish Balay . -tmp - alternative name of /tmp directory 1120e5c89e4eSSatish Balay . -get_total_flops - returns total flops done by all processors 11210841954dSBarry Smith - -memory_view - Print memory usage at end of run 1122e5c89e4eSSatish Balay 1123c5b5d8d5SVaclav Hapla Options Database Keys for Option Database: 1124c5b5d8d5SVaclav Hapla + -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc 1125c5b5d8d5SVaclav Hapla . -options_monitor - monitor all set options to standard output for the whole program run 1126c5b5d8d5SVaclav Hapla - -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet() 1127c5b5d8d5SVaclav Hapla 1128c5b5d8d5SVaclav Hapla Options -options_monitor_{all,cancel} are 1129c5b5d8d5SVaclav Hapla position-independent and apply to all options set since the PETSc start. 1130c5b5d8d5SVaclav Hapla They can be used also in option files. 1131c5b5d8d5SVaclav Hapla 1132c5b5d8d5SVaclav Hapla See PetscOptionsMonitorSet() to do monitoring programmatically. 1133c5b5d8d5SVaclav Hapla 1134e5c89e4eSSatish Balay Options Database Keys for Profiling: 1135a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 1136fe9b927eSVaclav Hapla + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo(). 1137217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 1138217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 1139495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 1140e5c89e4eSSatish Balay hangs without running in the debugger). See PetscLogTraceBegin(). 11419a9a5d4cSBarry Smith . -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView(). 114279dccf82SBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView(). 11439a9a5d4cSBarry Smith . -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the 1144495fc317SBarry Smith summary is written to the file. See PetscLogView(). 1145f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 1146495fc317SBarry Smith . -log_all [filename] - Logs extensive profiling information See PetscLogDump(). 1147495fc317SBarry Smith . -log [filename] - Logs basic profiline information See PetscLogDump(). 1148c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 114987c3beb6SLisandro Dalcin . -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off 1150c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 1151495fc317SBarry Smith 11526a77f485SPierre Jolivet Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time 1153e5c89e4eSSatish Balay 1154ffbd1cfbSBarry Smith Options Database Keys for SAWs: 1155ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 1156ffbd1cfbSBarry 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 1157ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 1158ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 1159ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 1160ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 1161ffbd1cfbSBarry Smith 1162e5c89e4eSSatish Balay Environmental Variables: 1163e5c89e4eSSatish Balay + PETSC_TMP - alternative tmp directory 1164e5c89e4eSSatish Balay . PETSC_SHARED_TMP - tmp is shared by all processes 1165e5c89e4eSSatish Balay . PETSC_NOT_SHARED_TMP - each process has its own private tmp 11664a971ea4SToby Isaac . PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs 11674a971ea4SToby Isaac . PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document 1168e5c89e4eSSatish Balay . PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer 1169e5c89e4eSSatish Balay - PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to 1170e5c89e4eSSatish Balay 1171e5c89e4eSSatish Balay Level: beginner 1172e5c89e4eSSatish Balay 1173e5c89e4eSSatish Balay Notes: 1174e5c89e4eSSatish Balay If for some reason you must call MPI_Init() separately, call 1175e5c89e4eSSatish Balay it before PetscInitialize(). 1176e5c89e4eSSatish Balay 1177e5c89e4eSSatish Balay Fortran Version: 1178e5c89e4eSSatish Balay In Fortran this routine has the format 1179e5c89e4eSSatish Balay $ call PetscInitialize(file,ierr) 1180e5c89e4eSSatish Balay 1181e5c89e4eSSatish Balay + ierr - error return code 1182c5b5d8d5SVaclav Hapla - file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc. 1183c5b5d8d5SVaclav Hapla Use PETSC_NULL_CHARACTER to not check for code specific file. 1184c5b5d8d5SVaclav Hapla Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files. 1185e5c89e4eSSatish Balay 1186e5c89e4eSSatish Balay Important Fortran Note: 11870eb4c9c0SBarry Smith In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a 11880298fd71SBarry Smith null character string; you CANNOT just use NULL as 1189a7f22e61SSatish Balay in the C version. See Users-Manual: ch_fortran for details. 1190e5c89e4eSSatish Balay 119101cb0274SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after 119201cb0274SBarry Smith calling PetscInitialize(). 1193e5c89e4eSSatish Balay 119401cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments() 1195e5c89e4eSSatish Balay 1196e5c89e4eSSatish Balay @*/ 11977087cfbeSBarry Smith PetscErrorCode PetscInitialize(int *argc,char ***args,const char file[],const char help[]) 1198e5c89e4eSSatish Balay { 119985649d77SJunchao Zhang PetscMPIInt flag; 120085649d77SJunchao Zhang const char *prog = "Unknown Name"; 1201e5c89e4eSSatish Balay 120227104ee2SJacob Faibussowitsch PetscFunctionBegin; 120327104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(0); 12049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Initialized(&flag)); 1205e5c89e4eSSatish Balay if (!flag) { 120639a651e2SJacob Faibussowitsch PetscCheck(PETSC_COMM_WORLD == MPI_COMM_NULL,PETSC_COMM_SELF,PETSC_ERR_SUP,"You cannot set PETSC_COMM_WORLD if you have not initialized MPI first"); 12079566063dSJacob Faibussowitsch PetscCall(PetscPreMPIInit_Private()); 12085e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 12095e765c61SJed Brown { 121039a651e2SJacob Faibussowitsch PetscMPIInt PETSC_UNUSED provided; 12119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided)); 12125e765c61SJed Brown } 12135e765c61SJed Brown #else 12149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init(argc,args)); 12155e765c61SJed Brown #endif 1216e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 1217e5c89e4eSSatish Balay } 12184dfee713SSatish Balay 121985649d77SJunchao Zhang if (argc && *argc) prog = **args; 1220e5c89e4eSSatish Balay if (argc && args) { 1221e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 1222e5c89e4eSSatish Balay PetscGlobalArgs = *args; 1223e5c89e4eSSatish Balay } 12249566063dSJacob Faibussowitsch PetscCall(PetscInitialize_Common(prog,file,help,PETSC_FALSE/*C*/,PETSC_FALSE,0)); 122527104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 1226e5c89e4eSSatish Balay } 1227e5c89e4eSSatish Balay 1228a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 122995c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1230ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1231ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 123295c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 12334097062eSBarry Smith #endif 1234e5c89e4eSSatish Balay 1235008a6e76SBarry Smith /* 1236008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1237008a6e76SBarry Smith */ 1238008a6e76SBarry Smith PetscErrorCode PetscFreeMPIResources(void) 1239008a6e76SBarry Smith { 1240008a6e76SBarry Smith PetscFunctionBegin; 1241008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 12429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128)); 1243008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX) 12449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128)); 1245008a6e76SBarry Smith #endif 12469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAX)); 12479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MIN)); 1248008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16) 12499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FP16)); 12509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAX)); 12519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MIN)); 1252008a6e76SBarry Smith #endif 1253008a6e76SBarry Smith 1254de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 12559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_SUM)); 1256008a6e76SBarry Smith #endif 1257008a6e76SBarry Smith 12589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR)); 12599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT)); 12609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT)); 126140df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 12629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2INT)); 1263008a6e76SBarry Smith #endif 12649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPI_4INT)); 12659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_4INT)); 12669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP)); 1267008a6e76SBarry Smith PetscFunctionReturn(0); 1268008a6e76SBarry Smith } 1269008a6e76SBarry Smith 1270a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 1271a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void); 1272a4af0ceeSJacob Faibussowitsch #endif 1273a4af0ceeSJacob Faibussowitsch 1274e5c89e4eSSatish Balay /*@C 1275e5c89e4eSSatish Balay PetscFinalize - Checks for options to be called at the conclusion 1276e5c89e4eSSatish Balay of the program. MPI_Finalize() is called only if the user had not 1277e5c89e4eSSatish Balay called MPI_Init() before calling PetscInitialize(). 1278e5c89e4eSSatish Balay 1279e5c89e4eSSatish Balay Collective on PETSC_COMM_WORLD 1280e5c89e4eSSatish Balay 1281e5c89e4eSSatish Balay Options Database Keys: 128226a7e8d4SBarry Smith + -options_view - Calls PetscOptionsView() 1283e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 12847eb1d149SBarry 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 1285e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 128692f119d6SBarry Smith . -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed 1287e5c89e4eSSatish Balay . -malloc_info - Prints total memory usage 128892f119d6SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and where 1289e5c89e4eSSatish Balay 1290e5c89e4eSSatish Balay Level: beginner 1291e5c89e4eSSatish Balay 1292e5c89e4eSSatish Balay Note: 1293e5c89e4eSSatish Balay See PetscInitialize() for more general runtime options. 1294e5c89e4eSSatish Balay 129588c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd() 1296e5c89e4eSSatish Balay @*/ 12977087cfbeSBarry Smith PetscErrorCode PetscFinalize(void) 1298e5c89e4eSSatish Balay { 12994bb5149bSJed Brown PetscMPIInt rank; 1300a8d2bbe5SBarry Smith PetscInt nopt; 13012bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE; 1302dff31646SBarry Smith PetscBool flg; 130310463e74SBarry Smith #if defined(PETSC_USE_LOG) 130410463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 130510463e74SBarry Smith #endif 1306e5c89e4eSSatish Balay 13073cbbc5ffSBarry Smith PetscFunctionBegin; 130839a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscInitialize() must be called before PetscFinalize()"); 13099566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"PetscFinalize() called\n")); 1310b022a5c1SBarry Smith 13119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 1312a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 13139566063dSJacob Faibussowitsch PetscCall(adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE)); 13149566063dSJacob Faibussowitsch PetscCall(adios_finalize(rank)); 1315a56f64adSBarry Smith #endif 13169566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-citations",&flg)); 1317dff31646SBarry Smith if (flg) { 13181f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 13191f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 13201f817a21SBarry Smith 13219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL)); 13221f817a21SBarry Smith if (filename[0]) { 13239566063dSJacob Faibussowitsch PetscCall(PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd)); 13241f817a21SBarry Smith } 13259566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(PetscCitationsList,1,&cits)); 1326dff31646SBarry Smith cits[0] = 0; 13279566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList,&cits)); 13289566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n")); 13299566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n")); 13309566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits)); 13319566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n")); 13329566063dSJacob Faibussowitsch PetscCall(PetscFClose(PETSC_COMM_WORLD,fd)); 13339566063dSJacob Faibussowitsch PetscCall(PetscFree(cits)); 1334dff31646SBarry Smith } 13359566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&PetscCitationsList)); 1336dff31646SBarry Smith 1337c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER) 133804102261SBarry Smith /* TextBelt is run for testing purposes only, please do not use this feature often */ 133904102261SBarry Smith { 134004102261SBarry Smith PetscInt nmax = 2; 134104102261SBarry Smith char **buffs; 13429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2,&buffs)); 13439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1)); 134404102261SBarry Smith if (flg1) { 134528b400f6SJacob Faibussowitsch PetscCheck(nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\""); 134604102261SBarry Smith if (nmax == 1) { 13479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(128,&buffs[1])); 13489566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(buffs[1],32)); 13499566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buffs[1]," has completed")); 135004102261SBarry Smith } 13519566063dSJacob Faibussowitsch PetscCall(PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL)); 13529566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[0])); 13539566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[1])); 135404102261SBarry Smith } 13559566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs)); 135604102261SBarry Smith } 1357763ec7b1SBarry Smith { 1358763ec7b1SBarry Smith PetscInt nmax = 2; 1359763ec7b1SBarry Smith char **buffs; 13609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2,&buffs)); 13619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1)); 1362763ec7b1SBarry Smith if (flg1) { 136328b400f6SJacob Faibussowitsch PetscCheck(nmax,PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\""); 1364763ec7b1SBarry Smith if (nmax == 1) { 13659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(128,&buffs[1])); 13669566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(buffs[1],32)); 13679566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buffs[1]," has completed")); 1368763ec7b1SBarry Smith } 13699566063dSJacob Faibussowitsch PetscCall(PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL)); 13709566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[0])); 13719566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[1])); 1372763ec7b1SBarry Smith } 13739566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs)); 1374763ec7b1SBarry Smith } 137504102261SBarry Smith #endif 137604102261SBarry Smith 13772d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 13789566063dSJacob Faibussowitsch PetscCall(PetscFPTDestroy()); 13792d53ad75SBarry Smith #endif 13802d53ad75SBarry Smith 1381e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1382dff31646SBarry Smith flg = PETSC_FALSE; 13839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL)); 1384d5649816SBarry Smith if (flg) { 13859566063dSJacob Faibussowitsch PetscCall(PetscOptionsSAWsDestroy()); 1386d5649816SBarry Smith } 1387d5649816SBarry Smith #endif 1388d5649816SBarry Smith 1389681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1390681455b2SBarry Smith flg1 = PETSC_FALSE; 13919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL)); 1392681455b2SBarry Smith if (flg1) { 1393681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 13949566063dSJacob Faibussowitsch PetscCall(PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL)); 1395681455b2SBarry Smith } 1396681455b2SBarry Smith #endif 1397681455b2SBarry Smith 139867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 13999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL)); 1400e5c89e4eSSatish Balay if (!flg2) { 140190d69ab7SBarry Smith flg2 = PETSC_FALSE; 14029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL)); 1403e5c89e4eSSatish Balay } 1404e5c89e4eSSatish Balay if (flg2) { 14059566063dSJacob Faibussowitsch PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n")); 1406e5c89e4eSSatish Balay } 140767584ceeSBarry Smith #endif 1408e5c89e4eSSatish Balay 1409e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 141090d69ab7SBarry Smith flg1 = PETSC_FALSE; 14119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL)); 1412e5c89e4eSSatish Balay if (flg1) { 1413e5c89e4eSSatish Balay PetscLogDouble flops = 0; 14149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD)); 14159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops)); 1416e5c89e4eSSatish Balay } 1417e5c89e4eSSatish Balay #endif 1418e5c89e4eSSatish Balay 1419e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 1420e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE) 1421e5c89e4eSSatish Balay mname[0] = 0; 14229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1)); 1423e5c89e4eSSatish Balay if (flg1) { 14249566063dSJacob Faibussowitsch if (mname[0]) PetscCall(PetscLogMPEDump(mname)); 14259566063dSJacob Faibussowitsch else PetscCall(PetscLogMPEDump(0)); 1426e5c89e4eSSatish Balay } 1427e5c89e4eSSatish Balay #endif 1428356e5837SBarry Smith #endif 1429a297a907SKarl Rupp 1430dd710f27SBarry Smith /* 1431dd710f27SBarry Smith Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 1432dd710f27SBarry Smith */ 14339566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1434dd710f27SBarry Smith 1435356e5837SBarry Smith #if defined(PETSC_USE_LOG) 14369566063dSJacob Faibussowitsch PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE)); 14379566063dSJacob Faibussowitsch PetscCall(PetscLogViewFromOptions()); 14389566063dSJacob Faibussowitsch PetscCall(PetscOptionsPopGetViewerOff()); 143987c3beb6SLisandro Dalcin 1440356e5837SBarry Smith mname[0] = 0; 14419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1)); 1442e5c89e4eSSatish Balay if (flg1) { 144391eabc43SBarry Smith PetscViewer viewer; 14449566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING: -log_summary is being deprecated; switch to -log_view\n\n\n")); 144591eabc43SBarry Smith if (mname[0]) { 14469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer)); 14479566063dSJacob Faibussowitsch PetscCall(PetscLogView(viewer)); 14489566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 144933f85c2fSBarry Smith } else { 145033f85c2fSBarry Smith viewer = PETSC_VIEWER_STDOUT_WORLD; 14519566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT)); 14529566063dSJacob Faibussowitsch PetscCall(PetscLogView(viewer)); 14539566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 145433f85c2fSBarry Smith } 1455e5c89e4eSSatish Balay } 1456a297a907SKarl Rupp 1457dd710f27SBarry Smith /* 1458dd710f27SBarry Smith Free any objects created by the last block of code. 1459dd710f27SBarry Smith */ 14609566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1461dd710f27SBarry Smith 1462dd710f27SBarry Smith mname[0] = 0; 14639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1)); 14649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2)); 14659566063dSJacob Faibussowitsch if (flg1 || flg2) PetscCall(PetscLogDump(mname)); 1466e5c89e4eSSatish Balay #endif 146710463e74SBarry Smith 146890d69ab7SBarry Smith flg1 = PETSC_FALSE; 14699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL)); 14709566063dSJacob Faibussowitsch if (!flg1) PetscCall(PetscPopSignalHandler()); 147190d69ab7SBarry Smith flg1 = PETSC_FALSE; 14729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL)); 1473e5c89e4eSSatish Balay if (flg1) { 14749566063dSJacob Faibussowitsch PetscCall(PetscMPIDump(stdout)); 1475e5c89e4eSSatish Balay } 147690d69ab7SBarry Smith flg1 = PETSC_FALSE; 147790d69ab7SBarry Smith flg2 = PETSC_FALSE; 14788bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 14799566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1)); 14809566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1)); 14819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL)); 1482e4c476e2SSatish Balay 1483e5c89e4eSSatish Balay if (flg2) { 1484be56827dSJed Brown PetscViewer viewer; 14859566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 14869566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer,PETSCVIEWERASCII)); 14879566063dSJacob Faibussowitsch PetscCall(PetscOptionsView(NULL,viewer)); 14889566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1489e5c89e4eSSatish Balay } 1490e5c89e4eSSatish Balay 1491e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 14929566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-nox",&flg1)); 14939566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1)); 1494e5c89e4eSSatish Balay 149533fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 14969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1)); 14979245e749SBarry Smith if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE; 1498e5c89e4eSSatish Balay if (flg3) { 14993de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 1500be56827dSJed Brown PetscViewer viewer; 15019566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 15029566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer,PETSCVIEWERASCII)); 15039566063dSJacob Faibussowitsch PetscCall(PetscOptionsView(NULL,viewer)); 15049566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1505e5c89e4eSSatish Balay } 15069566063dSJacob Faibussowitsch PetscCall(PetscOptionsAllUsed(NULL,&nopt)); 15073de2bfdfSBarry Smith if (nopt) { 15089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n")); 15099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n")); 15103de2bfdfSBarry Smith if (nopt == 1) { 15119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n")); 1512e5c89e4eSSatish Balay } else { 15139566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"There are %" PetscInt_FMT " unused database options. They are:\n",nopt)); 1514e5c89e4eSSatish Balay } 15153de2bfdfSBarry Smith } else if (flg3 && flg1) { 15169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n")); 1517df12ba86SBarry Smith } 15189566063dSJacob Faibussowitsch PetscCall(PetscOptionsLeft(NULL)); 1519e5c89e4eSSatish Balay } 1520e5c89e4eSSatish Balay 1521e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1522d45a07a7SBarry Smith if (!PetscGlobalRank) { 15239566063dSJacob Faibussowitsch PetscCall(PetscStackSAWsViewOff()); 152416ad0300SBarry Smith PetscStackCallSAWs(SAWs_Finalize,()); 1525d45a07a7SBarry Smith } 1526ec957eceSBarry Smith #endif 1527ec957eceSBarry Smith 15284097062eSBarry Smith #if defined(PETSC_USE_LOG) 152910463e74SBarry Smith /* 1530dbc8283eSBarry Smith List all objects the user may have forgot to free 15312eff7a51SBarry Smith */ 153205df10baSBarry Smith if (PetscObjectsLog) { 15339566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1)); 1534a64a8e02SBarry Smith if (flg1) { 1535a64a8e02SBarry Smith MPI_Comm local_comm; 15367eb1d149SBarry Smith char string[64]; 1537a64a8e02SBarry Smith 15389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL)); 15399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(MPI_COMM_WORLD,&local_comm)); 15409566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm,1)); 15419566063dSJacob Faibussowitsch PetscCall(PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE)); 15429566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm,1)); 15439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 15440a1571b3SBarry Smith } 154505df10baSBarry Smith } 15464097062eSBarry Smith #endif 15474097062eSBarry Smith 15484097062eSBarry Smith #if defined(PETSC_USE_LOG) 1549dbc8283eSBarry Smith PetscObjectsCounts = 0; 1550dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 15519566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscObjects)); 15524097062eSBarry Smith #endif 15532eff7a51SBarry Smith 155433f85c2fSBarry Smith /* 155533f85c2fSBarry Smith Destroy any packages that registered a finalize 155633f85c2fSBarry Smith */ 15579566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalizeAll()); 155833f85c2fSBarry Smith 1559101409b8SToby Isaac #if defined(PETSC_USE_LOG) 15609566063dSJacob Faibussowitsch PetscCall(PetscLogFinalize()); 1561101409b8SToby Isaac #endif 1562101409b8SToby Isaac 156333f85c2fSBarry Smith /* 156448dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 156548dd1dffSBarry Smith 15669566063dSJacob Faibussowitsch PetscCall(PetscFunctionListPrintAll()); 156748dd1dffSBarry Smith */ 156837e93019SBarry Smith 15694028d114SSatish Balay if (petsc_history) { 15709566063dSJacob Faibussowitsch PetscCall(PetscCloseHistoryFile(&petsc_history)); 157102c9f0b5SLisandro Dalcin petsc_history = NULL; 1572e5c89e4eSSatish Balay } 15739566063dSJacob Faibussowitsch PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton)); 15749566063dSJacob Faibussowitsch PetscCall(PetscInfoDestroy()); 1575e5c89e4eSSatish Balay 157667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 157792f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1578e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 157992f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1580e5c89e4eSSatish Balay FILE *fd; 1581ed9cf6e9SBarry Smith int err; 1582e5c89e4eSSatish Balay 1583dc92acbaSJed Brown flg2 = PETSC_FALSE; 158492f119d6SBarry Smith flg3 = PETSC_FALSE; 15859566063dSJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL)); 15869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL)); 158792f119d6SBarry Smith fname[0] = 0; 15889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1)); 1589e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1590e5c89e4eSSatish Balay 1591589a23caSBarry Smith PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank); 159228b400f6SJacob Faibussowitsch fd = fopen(sname,"w"); PetscCheck(fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname); 15939566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(fd)); 1594ed9cf6e9SBarry Smith err = fclose(fd); 159528b400f6SJacob Faibussowitsch PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 159692f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1597e5c89e4eSSatish Balay MPI_Comm local_comm; 1598e5c89e4eSSatish Balay 15999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(MPI_COMM_WORLD,&local_comm)); 16009566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm,1)); 16019566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(stdout)); 16029566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm,1)); 16039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1604e5c89e4eSSatish Balay } 1605e5c89e4eSSatish Balay fname[0] = 0; 16069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1)); 1607e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1608e5c89e4eSSatish Balay 1609589a23caSBarry Smith PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank); 161028b400f6SJacob Faibussowitsch fd = fopen(sname,"w"); PetscCheck(fd,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname); 16119566063dSJacob Faibussowitsch PetscCall(PetscMallocView(fd)); 1612ed9cf6e9SBarry Smith err = fclose(fd); 161328b400f6SJacob Faibussowitsch PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 161492f119d6SBarry Smith } else if (flg1) { 161592f119d6SBarry Smith MPI_Comm local_comm; 161692f119d6SBarry Smith 16179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(MPI_COMM_WORLD,&local_comm)); 16189566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm,1)); 16199566063dSJacob Faibussowitsch PetscCall(PetscMallocView(stdout)); 16209566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm,1)); 16219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1622e5c89e4eSSatish Balay } 1623e5c89e4eSSatish Balay } 162467584ceeSBarry Smith #endif 162520e2c332SMatthew G. Knepley 16265486ca60SMatthew G. Knepley /* 16275486ca60SMatthew G. Knepley Close any open dynamic libraries 16285486ca60SMatthew G. Knepley */ 16299566063dSJacob Faibussowitsch PetscCall(PetscFinalize_DynamicLibraries()); 16305486ca60SMatthew G. Knepley 1631e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 16329566063dSJacob Faibussowitsch PetscCall(PetscOptionsDestroyDefault()); 1633e5c89e4eSSatish Balay 1634e5c89e4eSSatish Balay PetscGlobalArgc = 0; 163502c9f0b5SLisandro Dalcin PetscGlobalArgs = NULL; 1636e5c89e4eSSatish Balay 1637c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 1638c2b86a48SJunchao Zhang if (PetscBeganKokkos) { 16399566063dSJacob Faibussowitsch PetscCall(PetscKokkosFinalize_Private()); 1640c2b86a48SJunchao Zhang PetscBeganKokkos = PETSC_FALSE; 164145639126SStefano Zampini PetscKokkosInitialized = PETSC_FALSE; 1642c2b86a48SJunchao Zhang } 1643c2b86a48SJunchao Zhang #endif 1644c2b86a48SJunchao Zhang 164571438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM) 164671438e86SJunchao Zhang if (PetscBeganNvshmem) { 16479566063dSJacob Faibussowitsch PetscCall(PetscNvshmemFinalize()); 164871438e86SJunchao Zhang PetscBeganNvshmem = PETSC_FALSE; 164971438e86SJunchao Zhang } 165071438e86SJunchao Zhang #endif 1651a0e72f99SJunchao Zhang 16529566063dSJacob Faibussowitsch PetscCall(PetscFreeMPIResources()); 1653e5c89e4eSSatish Balay 1654dbc8283eSBarry Smith /* 1655efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1656efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1657efb80d3cSBarry Smith 1658efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1659efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1660dbc8283eSBarry Smith */ 1661b770b1f6SSatish Balay { 1662dbc8283eSBarry Smith PetscCommCounter *counter; 1663dbc8283eSBarry Smith PetscMPIInt flg; 1664dbc8283eSBarry Smith MPI_Comm icomm; 1665265f3f35SJed Brown union {MPI_Comm comm; void *ptr;} ucomm; 16669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg)); 1667dbc8283eSBarry Smith if (flg) { 1668265f3f35SJed Brown icomm = ucomm.comm; 16699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg)); 167028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1671dbc8283eSBarry Smith 16729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval)); 16739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval)); 16749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1675dbc8283eSBarry Smith } 16769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg)); 1677dbc8283eSBarry Smith if (flg) { 1678265f3f35SJed Brown icomm = ucomm.comm; 16799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg)); 168028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1681dbc8283eSBarry Smith 16829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval)); 16839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval)); 16849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1685dbc8283eSBarry Smith } 1686b770b1f6SSatish Balay } 1687dbc8283eSBarry Smith 16889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval)); 16899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval)); 16909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval)); 16919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval)); 1692480cf27aSJed Brown 16939566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen)); 16949566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout)); 16959566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr)); 16969566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock)); 1697ef19f930SBarry Smith 1698e5c89e4eSSatish Balay if (PetscBeganMPI) { 169999b1327fSBarry Smith PetscMPIInt flag; 17009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalized(&flag)); 170139a651e2SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 170239a651e2SJacob Faibussowitsch /* wait until the very last moment to disable error handling */ 170339a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_FALSE; 17049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalize()); 170539a651e2SJacob Faibussowitsch } else PetscErrorHandlingInitialized = PETSC_FALSE; 170639a651e2SJacob Faibussowitsch 1707e5c89e4eSSatish Balay /* 1708e5c89e4eSSatish Balay 1709e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1710e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1711e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1712e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1713e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 17140e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1715e5c89e4eSSatish Balay memory was not freed. 1716e5c89e4eSSatish Balay 1717e5c89e4eSSatish Balay */ 17189566063dSJacob Faibussowitsch PetscCall(PetscMallocClear()); 17199566063dSJacob Faibussowitsch PetscCall(PetscStackReset()); 1720a297a907SKarl Rupp 1721e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1722e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 172356883f60SBarry Smith #if defined(PETSC_USE_GCOV) 172456883f60SBarry Smith /* 172556883f60SBarry Smith flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the 172656883f60SBarry Smith gcov files are still being added to the directories as git tries to remove the directories. 172756883f60SBarry Smith */ 172856883f60SBarry Smith __gcov_flush(); 172956883f60SBarry Smith #endif 17301724198aSStefano Zampini /* To match PetscFunctionBegin() at the beginning of this function */ 17311724198aSStefano Zampini PetscStackClearTop; 173227104ee2SJacob Faibussowitsch return 0; 1733e5c89e4eSSatish Balay } 1734e5c89e4eSSatish Balay 173543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 17368cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b) 173743db4dbbSBarry Smith { 173843db4dbbSBarry Smith if (*a == *b) return 1; 173943db4dbbSBarry Smith if (*a + 32 == *b) return 1; 174043db4dbbSBarry Smith if (*a - 32 == *b) return 1; 174143db4dbbSBarry Smith return 0; 174243db4dbbSBarry Smith } 1743a70650f6SBarry Smith #endif 174443db4dbbSBarry Smith 174543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 17468cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b) 174743db4dbbSBarry Smith { 174843db4dbbSBarry Smith if (*a == *b) return 1; 174943db4dbbSBarry Smith if (*a + 32 == *b) return 1; 175043db4dbbSBarry Smith if (*a - 32 == *b) return 1; 175143db4dbbSBarry Smith return 0; 175243db4dbbSBarry Smith } 175343db4dbbSBarry Smith #endif 1754