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*/ 6473903fcSJunchao Zhang #include <petsc/private/logimpl.h> 7665c2dedSJed Brown #include <petscviewer.h> 862e5d2d2SJDBetteridge #include <petsc/private/garbagecollector.h> 9a0e72f99SJunchao Zhang 106edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS) 116edef35eSSatish Balay #include <petsc/private/valgrind/valgrind.h> 126edef35eSSatish Balay #endif 136edef35eSSatish Balay 14fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS) 1585649d77SJunchao Zhang #include <petsc/private/fortranimpl.h> 1685649d77SJunchao Zhang #endif 1785649d77SJunchao Zhang 187ce81a4bSJacob Faibussowitsch #if PetscDefined(USE_COVERAGE) 1956883f60SBarry Smith EXTERN_C_BEGIN 20aaf3846cSSatish Balay #if defined(PETSC_HAVE___GCOV_DUMP) 21aaf3846cSSatish Balay #define __gcov_flush(x) __gcov_dump(x) 22aaf3846cSSatish Balay #endif 2356883f60SBarry Smith void __gcov_flush(void); 2456883f60SBarry Smith EXTERN_C_END 2556883f60SBarry Smith #endif 268101f56cSMatthew Knepley 272d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 2895c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData; 292d53ad75SBarry Smith PetscFPT PetscFPTData = 0; 302d53ad75SBarry Smith #endif 312d53ad75SBarry Smith 3227104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 3316ad0300SBarry Smith #include <petscviewersaws.h> 34a6790183SBarry Smith #endif 35eb02082bSJunchao Zhang 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 PetscSequentialPhaseBegin_Private(MPI_Comm, int); 4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm, int); 4295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE **); 430069ddf5SShri Abhyankar 446de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */ 45e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL; 4627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD) 47835d5ba2SJunchao Zhang PetscMPIInt PETSC_MPI_THREAD_REQUIRED = PETSC_DECIDE; 486de5d289SStefano Zampini #else 49835d5ba2SJunchao Zhang PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_SINGLE; 506de5d289SStefano Zampini #endif 51e5c89e4eSSatish Balay 52480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval = MPI_KEYVAL_INVALID; 53480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID; 54480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID; 555f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval = MPI_KEYVAL_INVALID; 5662e5d2d2SJDBetteridge PetscMPIInt Petsc_CreationIdx_keyval = MPI_KEYVAL_INVALID; 5762e5d2d2SJDBetteridge PetscMPIInt Petsc_Garbage_HMap_keyval = MPI_KEYVAL_INVALID; 58480cf27aSJed Brown 595ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedWD_keyval = MPI_KEYVAL_INVALID; 605ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedTmp_keyval = MPI_KEYVAL_INVALID; 615ea2b939SDuncan Campbell 62e5c89e4eSSatish Balay /* 63e5c89e4eSSatish Balay Declare and set all the string names of the PETSc enums 64e5c89e4eSSatish Balay */ 6502c9f0b5SLisandro Dalcin const char *const PetscBools[] = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL}; 6602c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL}; 67e5c89e4eSSatish Balay 68ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE; 69ace3abfcSBarry Smith PetscBool PetscPreLoadingOn = PETSC_FALSE; 700f8e0872SSatish Balay 71a2f94806SJed Brown PetscInt PetscHotRegionDepth; 72a2f94806SJed Brown 736edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE; 746edef35eSSatish Balay 75b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY) 76b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen; 77b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout; 78b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr; 79b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock; 80b22622e2STadeu Manoel #endif 81b22622e2STadeu Manoel 82046ed546SBarry Smith extern PetscInt PetscNumBLASThreads; 83046ed546SBarry Smith 84aec76313SJacob Faibussowitsch /*@C 85945d1669SBarry Smith PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args 8672a42c3cSBarry Smith 8772a42c3cSBarry Smith Collective 8872a42c3cSBarry Smith 8910450e9eSJacob Faibussowitsch Input Parameters: 9010450e9eSJacob Faibussowitsch + argc - number of args 9110450e9eSJacob Faibussowitsch . args - array of command line arguments 9210450e9eSJacob Faibussowitsch . filename - optional name of the program file, pass `NULL` to ignore 9310450e9eSJacob Faibussowitsch - help - optional help, pass `NULL` to ignore 9410450e9eSJacob Faibussowitsch 9572a42c3cSBarry Smith Level: advanced 9672a42c3cSBarry Smith 9795452b02SPatrick Sanan Notes: 98a56f64adSBarry Smith this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to 99a3b724e8SBarry Smith indicate that it did NOT start MPI so that the `PetscFinalize()` does not end MPI, thus allowing `PetscInitialize()` to 100a56f64adSBarry Smith be called multiple times from Julia without the problem of trying to initialize MPI more than once. 1010f11a792SBarry Smith 10210450e9eSJacob Faibussowitsch Developer Notes: 10310450e9eSJacob Faibussowitsch Turns off PETSc signal handling to allow Julia to manage signals 1041ea5a559SBarry Smith 105db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()` 1060f11a792SBarry Smith */ 107d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoPointers(int argc, char **args, const char *filename, const char *help) 108d71ae5a4SJacob Faibussowitsch { 10972a42c3cSBarry Smith int myargc = argc; 11072a42c3cSBarry Smith char **myargs = args; 11172a42c3cSBarry Smith 11272a42c3cSBarry Smith PetscFunctionBegin; 1139566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&myargc, &myargs, filename, help)); 1149566063dSJacob Faibussowitsch PetscCall(PetscPopSignalHandler()); 115df413903SBarry Smith PetscBeganMPI = PETSC_FALSE; 1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11772a42c3cSBarry Smith } 11872a42c3cSBarry Smith 119e5c89e4eSSatish Balay /*@C 120811af0c4SBarry Smith PetscInitializeNoArguments - Calls `PetscInitialize()` from C/C++ without 121e5c89e4eSSatish Balay the command line arguments. 122e5c89e4eSSatish Balay 123e5c89e4eSSatish Balay Collective 124e5c89e4eSSatish Balay 125e5c89e4eSSatish Balay Level: advanced 126e5c89e4eSSatish Balay 127db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()` 128e5c89e4eSSatish Balay @*/ 129d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoArguments(void) 130d71ae5a4SJacob Faibussowitsch { 131e5c89e4eSSatish Balay int argc = 0; 13202c9f0b5SLisandro Dalcin char **args = NULL; 133e5c89e4eSSatish Balay 134e5c89e4eSSatish Balay PetscFunctionBegin; 1359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, NULL)); 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137e5c89e4eSSatish Balay } 138e5c89e4eSSatish Balay 139e5c89e4eSSatish Balay /*@ 140e5c89e4eSSatish Balay PetscInitialized - Determine whether PETSc is initialized. 141e5c89e4eSSatish Balay 14210450e9eSJacob Faibussowitsch Output Parameter: 14310450e9eSJacob Faibussowitsch . isInitialized - `PETSC_TRUE` if PETSc is initialized, `PETSC_FALSE` otherwise 14410450e9eSJacob Faibussowitsch 14593b6d2d1SJed Brown Level: beginner 146e5c89e4eSSatish Balay 147db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()` 148e5c89e4eSSatish Balay @*/ 149d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialized(PetscBool *isInitialized) 150d71ae5a4SJacob Faibussowitsch { 15127104ee2SJacob Faibussowitsch PetscFunctionBegin; 1524f572ea9SToby Isaac if (PetscInitializeCalled) PetscAssertPointer(isInitialized, 1); 153e5c89e4eSSatish Balay *isInitialized = PetscInitializeCalled; 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155e5c89e4eSSatish Balay } 156e5c89e4eSSatish Balay 157e5c89e4eSSatish Balay /*@ 158811af0c4SBarry Smith PetscFinalized - Determine whether `PetscFinalize()` has been called yet 159e5c89e4eSSatish Balay 16010450e9eSJacob Faibussowitsch Output Parameter: 16110450e9eSJacob Faibussowitsch . isFinalized - `PETSC_TRUE` if PETSc is finalized, `PETSC_FALSE` otherwise 16210450e9eSJacob Faibussowitsch 163e5c89e4eSSatish Balay Level: developer 164e5c89e4eSSatish Balay 165db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()` 166e5c89e4eSSatish Balay @*/ 167d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalized(PetscBool *isFinalized) 168d71ae5a4SJacob Faibussowitsch { 16927104ee2SJacob Faibussowitsch PetscFunctionBegin; 1704f572ea9SToby Isaac if (!PetscFinalizeCalled) PetscAssertPointer(isFinalized, 1); 171e5c89e4eSSatish Balay *isFinalized = PetscFinalizeCalled; 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173e5c89e4eSSatish Balay } 174e5c89e4eSSatish Balay 17557171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]); 176e5c89e4eSSatish Balay 177e5c89e4eSSatish Balay /* 178e5c89e4eSSatish Balay This function is the MPI reduction operation used to compute the sum of the 179e5c89e4eSSatish Balay first half of the datatype and the max of the second half. 180e5c89e4eSSatish Balay */ 181367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0; 18262e5d2d2SJDBetteridge MPI_Op Petsc_Garbage_SetIntersectOp = 0; 183e5c89e4eSSatish Balay 184d71ae5a4SJacob Faibussowitsch PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, int *cnt, MPI_Datatype *datatype) 185d71ae5a4SJacob Faibussowitsch { 186e5c89e4eSSatish Balay PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out, i, count = *cnt; 187e5c89e4eSSatish Balay 188e5c89e4eSSatish Balay PetscFunctionBegin; 189e5c89e4eSSatish Balay if (*datatype != MPIU_2INT) { 1903ba16761SJacob Faibussowitsch PetscErrorCode ierr = (*PetscErrorPrintf)("Can only handle MPIU_2INT data types"); 1913ba16761SJacob Faibussowitsch (void)ierr; 19241e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 193e5c89e4eSSatish Balay } 194e5c89e4eSSatish Balay 195e5c89e4eSSatish Balay for (i = 0; i < count; i++) { 196e5c89e4eSSatish Balay xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]); 197e5c89e4eSSatish Balay xout[2 * i + 1] += xin[2 * i + 1]; 198e5c89e4eSSatish Balay } 199812af9f3SBarry Smith PetscFunctionReturnVoid(); 200e5c89e4eSSatish Balay } 201e5c89e4eSSatish Balay 202e5c89e4eSSatish Balay /* 203e5c89e4eSSatish Balay Returns the max of the first entry owned by this processor and the 204e5c89e4eSSatish Balay sum of the second entry. 205b693b147SBarry Smith 20676f543a4SJed Brown The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero 207367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths 208b693b147SBarry Smith there would be no place to store the both needed results. 209e5c89e4eSSatish Balay */ 210d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt sizes[], PetscInt *max, PetscInt *sum) 211d71ae5a4SJacob Faibussowitsch { 212e5c89e4eSSatish Balay PetscFunctionBegin; 213d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 214d6e4c47cSJed Brown { 2159371c9d4SSatish Balay struct { 2169371c9d4SSatish Balay PetscInt max, sum; 2179371c9d4SSatish Balay } work; 2189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce_scatter_block((void *)sizes, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm)); 219d6e4c47cSJed Brown *max = work.max; 220d6e4c47cSJed Brown *sum = work.sum; 221d6e4c47cSJed Brown } 222d6e4c47cSJed Brown #else 223d6e4c47cSJed Brown { 224d6e4c47cSJed Brown PetscMPIInt size, rank; 2259371c9d4SSatish Balay struct { 2269371c9d4SSatish Balay PetscInt max, sum; 2279371c9d4SSatish Balay } *work; 2289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &work)); 2311c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((void *)sizes, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm)); 2326ac3741eSJed Brown *max = work[rank].max; 2336ac3741eSJed Brown *sum = work[rank].sum; 2349566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 235d6e4c47cSJed Brown } 236d6e4c47cSJed Brown #endif 2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 238e5c89e4eSSatish Balay } 239e5c89e4eSSatish Balay 240*ef384fdeSPierre Jolivet #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16)) 241*ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 242613bf2b2SPierre Jolivet #include <quadmath.h> 243613bf2b2SPierre Jolivet #endif 2449e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0; 245de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 24606a205a8SBarry Smith MPI_Op MPIU_SUM = 0; 247613bf2b2SPierre Jolivet #endif 248e5c89e4eSSatish Balay 249d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 250d71ae5a4SJacob Faibussowitsch { 251e5c89e4eSSatish Balay PetscInt i, count = *cnt; 252e5c89e4eSSatish Balay 253e5c89e4eSSatish Balay PetscFunctionBegin; 2547c2de775SJed Brown if (*datatype == MPIU_REAL) { 255e2e03761SBarry Smith PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 256a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] += xin[i]; 2577c2de775SJed Brown } 2587c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 259c5481aeeSJose E. Roman else if (*datatype == MPIU_COMPLEX) { 2607c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 261a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] += xin[i]; 2627c2de775SJed Brown } 2637c2de775SJed Brown #endif 264*ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 265613bf2b2SPierre Jolivet else if (*datatype == MPIU___FLOAT128) { 266613bf2b2SPierre Jolivet __float128 *xin = (__float128 *)in, *xout = (__float128 *)out; 267613bf2b2SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 26874c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 2699371c9d4SSatish Balay } else if (*datatype == MPIU___COMPLEX128) { 270613bf2b2SPierre Jolivet __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out; 271613bf2b2SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 27274c01ffaSSatish Balay #endif 273613bf2b2SPierre Jolivet } 274613bf2b2SPierre Jolivet #endif 275*ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16) 2769e517322SPierre Jolivet else if (*datatype == MPIU___FP16) { 2779e517322SPierre Jolivet __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out; 2789e517322SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 2799e517322SPierre Jolivet } 2809e517322SPierre Jolivet #endif 2817c2de775SJed Brown else { 282*ef384fdeSPierre Jolivet #if (!defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_SKIP_REAL___FLOAT128)) && (!defined(PETSC_HAVE_REAL___FP16) || defined(PETSC_SKIP_REAL___FP16)) 283*ef384fdeSPierre Jolivet PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types")); 284*ef384fdeSPierre Jolivet #elif !defined(PETSC_HAVE_REAL___FP16) || defined(PETSC_SKIP_REAL___FP16) 2853ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types")); 286*ef384fdeSPierre Jolivet #elif !defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_SKIP_REAL___FLOAT128) 2873ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types")); 2889e517322SPierre Jolivet #else 2893ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types")); 290613bf2b2SPierre Jolivet #endif 29141e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 292e2e03761SBarry Smith } 293812af9f3SBarry Smith PetscFunctionReturnVoid(); 294e5c89e4eSSatish Balay } 295e5c89e4eSSatish Balay #endif 296e5c89e4eSSatish Balay 297570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 298d9822059SBarry Smith MPI_Op MPIU_MAX = 0; 299d9822059SBarry Smith MPI_Op MPIU_MIN = 0; 300d9822059SBarry Smith 301d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 302d71ae5a4SJacob Faibussowitsch { 303d9822059SBarry Smith PetscInt i, count = *cnt; 304d9822059SBarry Smith 305d9822059SBarry Smith PetscFunctionBegin; 3067c2de775SJed Brown if (*datatype == MPIU_REAL) { 3078c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 308a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]); 3097c2de775SJed Brown } 3107c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 3117c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 3127c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 313ad540459SPierre Jolivet for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 3147c2de775SJed Brown } 3157c2de775SJed Brown #endif 3167c2de775SJed Brown else { 3173ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types")); 31841e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 3198c764dc5SJose Roman } 320d9822059SBarry Smith PetscFunctionReturnVoid(); 321d9822059SBarry Smith } 322d9822059SBarry Smith 323d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 324d71ae5a4SJacob Faibussowitsch { 325d9822059SBarry Smith PetscInt i, count = *cnt; 326d9822059SBarry Smith 327d9822059SBarry Smith PetscFunctionBegin; 3287c2de775SJed Brown if (*datatype == MPIU_REAL) { 3298c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 330a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]); 3317c2de775SJed Brown } 3327c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 3337c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 3347c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 335ad540459SPierre Jolivet for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 3367c2de775SJed Brown } 3377c2de775SJed Brown #endif 3387c2de775SJed Brown else { 3393ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types")); 34041e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 3418c764dc5SJose Roman } 342d9822059SBarry Smith PetscFunctionReturnVoid(); 343d9822059SBarry Smith } 344d9822059SBarry Smith #endif 345d9822059SBarry Smith 346480cf27aSJed Brown /* 347480cf27aSJed Brown Private routine to delete internal tag/name counter storage when a communicator is freed. 348480cf27aSJed Brown 349ff0e51ddSBarry 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. 350480cf27aSJed Brown 35112801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 352480cf27aSJed Brown 353480cf27aSJed Brown */ 3548434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_DeleteFn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state) 355d71ae5a4SJacob Faibussowitsch { 35605342407SJunchao Zhang PetscCommCounter *counter = (PetscCommCounter *)count_val; 35757f21012SBarry Smith struct PetscCommStash *comms = counter->comms, *pcomm; 358480cf27aSJed Brown 359480cf27aSJed Brown PetscFunctionBegin; 3609566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm)); 3619566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter->iflags)); 36257f21012SBarry Smith while (comms) { 3639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&comms->comm)); 36457f21012SBarry Smith pcomm = comms; 36557f21012SBarry Smith comms = comms->next; 3669566063dSJacob Faibussowitsch PetscCall(PetscFree(pcomm)); 36757f21012SBarry Smith } 3689566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter)); 369480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 370480cf27aSJed Brown } 371480cf27aSJed Brown 372480cf27aSJed Brown /* 37347435625SJed Brown This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user 374da3039f7SJed Brown calls MPI_Comm_free(). 375da3039f7SJed Brown 376da3039f7SJed Brown This is the only entry point for breaking the links between inner and outer comms. 377480cf27aSJed Brown 378ff0e51ddSBarry Smith This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator. 379480cf27aSJed Brown 38012801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 381480cf27aSJed Brown 382480cf27aSJed Brown */ 3838434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_DeleteFn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) 384d71ae5a4SJacob Faibussowitsch { 3859371c9d4SSatish Balay union 386480cf27aSJed Brown { 3879371c9d4SSatish Balay MPI_Comm comm; 3889371c9d4SSatish Balay void *ptr; 3899371c9d4SSatish Balay } icomm; 390480cf27aSJed Brown 391480cf27aSJed Brown PetscFunctionBegin; 39212801b39SBarry Smith if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval"); 393ec4fadc2SJed Brown icomm.ptr = attr_val; 39476bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 39533779a13SJunchao Zhang /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */ 39633779a13SJunchao Zhang PetscMPIInt flg; 3979371c9d4SSatish Balay union 3989371c9d4SSatish Balay { 3999371c9d4SSatish Balay MPI_Comm comm; 4009371c9d4SSatish Balay void *ptr; 4019371c9d4SSatish Balay } ocomm; 4029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg)); 40333779a13SJunchao Zhang if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute"); 40433779a13SJunchao 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"); 40533779a13SJunchao Zhang } 4069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval)); 4079566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm)); 408da3039f7SJed Brown PetscFunctionReturn(MPI_SUCCESS); 409b89831e5SBarry Smith } 410da3039f7SJed Brown 411da3039f7SJed Brown /* 4128434afd1SBarry Smith * This is invoked on the inner comm when Petsc_InnerComm_Attr_DeleteFn calls MPI_Comm_delete_attr(). It should not be reached any other way. 413da3039f7SJed Brown */ 4148434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_DeleteFn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) 415d71ae5a4SJacob Faibussowitsch { 416da3039f7SJed Brown PetscFunctionBegin; 4179566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm)); 418480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 419480cf27aSJed Brown } 420480cf27aSJed Brown 4218434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_DeleteFn(MPI_Comm, PetscMPIInt, void *, void *); 42242218b76SBarry Smith 423951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 4248cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *); 4258cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 4268cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 427e39fd77fSBarry Smith #endif 428e39fd77fSBarry Smith 4290f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE; 43012801b39SBarry Smith 431eb27c7c8SSatish Balay PETSC_INTERN int PetscGlobalArgc; 432eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs; 4336ae9a8a6SBarry Smith int PetscGlobalArgc = 0; 43402c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL; 435dff31646SBarry Smith PetscSegBuffer PetscCitationsList; 436e5c89e4eSSatish Balay 437d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void) 438d71ae5a4SJacob Faibussowitsch { 439051e4cf2SJed Brown PetscFunctionBegin; 4409566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList)); 44130c35bf2SSatish Balay 44230c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\ 44330c35bf2SSatish Balay Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\ 44430c35bf2SSatish Balay and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\ 4453f81df79SBarry Smith and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\ 44630c35bf2SSatish Balay and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\ 44730c35bf2SSatish Balay and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\ 44830c35bf2SSatish Balay and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\ 44930c35bf2SSatish Balay and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\ 45030c35bf2SSatish Balay Title = {{PETSc/TAO} Users Manual},\n\ 45130694c2dSSatish Balay Number = {ANL-21/39 - Revision 3.21},\n\ 452f5a9e7c7SSatish Balay Doi = {10.2172/2205494},\n\ 45330c35bf2SSatish Balay Institution = {Argonne National Laboratory},\n\ 45430694c2dSSatish Balay Year = {2024}\n}\n", 4559371c9d4SSatish Balay NULL)); 45630c35bf2SSatish Balay 45730c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\ 45830c35bf2SSatish Balay Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\ 45930c35bf2SSatish Balay Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\ 46030c35bf2SSatish Balay Booktitle = {Modern Software Tools in Scientific Computing},\n\ 46130c35bf2SSatish Balay Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\ 46230c35bf2SSatish Balay Pages = {163--202},\n\ 46330c35bf2SSatish Balay Publisher = {Birkh{\\\"{a}}user Press},\n\ 4649371c9d4SSatish Balay Year = {1997}\n}\n", 4659371c9d4SSatish Balay NULL)); 4663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 467051e4cf2SJed Brown } 468e5c89e4eSSatish Balay 4692d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */ 4702d747510SLisandro Dalcin 471d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[]) 472d71ae5a4SJacob Faibussowitsch { 4732d747510SLisandro Dalcin PetscFunctionBegin; 4749566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(programname, name, sizeof(programname))); 4753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4762d747510SLisandro Dalcin } 4772d747510SLisandro Dalcin 4782d747510SLisandro Dalcin /*@C 4792d747510SLisandro Dalcin PetscGetProgramName - Gets the name of the running program. 4802d747510SLisandro Dalcin 4812d747510SLisandro Dalcin Not Collective 4822d747510SLisandro Dalcin 4832d747510SLisandro Dalcin Input Parameter: 4842d747510SLisandro Dalcin . len - length of the string name 4852d747510SLisandro Dalcin 4862d747510SLisandro Dalcin Output Parameter: 487811af0c4SBarry Smith . name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN` 4882d747510SLisandro Dalcin 4892d747510SLisandro Dalcin Level: advanced 4902d747510SLisandro Dalcin 49121532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()` 4922d747510SLisandro Dalcin @*/ 493d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProgramName(char name[], size_t len) 494d71ae5a4SJacob Faibussowitsch { 4952d747510SLisandro Dalcin PetscFunctionBegin; 4969566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(name, programname, len)); 4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4982d747510SLisandro Dalcin } 4992d747510SLisandro Dalcin 500e5c89e4eSSatish Balay /*@C 501e5c89e4eSSatish Balay PetscGetArgs - Allows you to access the raw command line arguments anywhere 502811af0c4SBarry Smith after PetscInitialize() is called but before `PetscFinalize()`. 503e5c89e4eSSatish Balay 504e5c89e4eSSatish Balay Not Collective 505e5c89e4eSSatish Balay 506e5c89e4eSSatish Balay Output Parameters: 507e5c89e4eSSatish Balay + argc - count of number of command line arguments 508e5c89e4eSSatish Balay - args - the command line arguments 509e5c89e4eSSatish Balay 510e5c89e4eSSatish Balay Level: intermediate 511e5c89e4eSSatish Balay 512e5c89e4eSSatish Balay Notes: 513e5c89e4eSSatish Balay This is usually used to pass the command line arguments into other libraries 514e5c89e4eSSatish Balay that are called internally deep in PETSc or the application. 515e5c89e4eSSatish Balay 516f177e3b1SBarry Smith The first argument contains the program name as is normal for C arguments. 517f177e3b1SBarry Smith 51821532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()` 519e5c89e4eSSatish Balay @*/ 520d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args) 521d71ae5a4SJacob Faibussowitsch { 522e5c89e4eSSatish Balay PetscFunctionBegin; 52339a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 524e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 525e5c89e4eSSatish Balay *args = PetscGlobalArgs; 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 527e5c89e4eSSatish Balay } 528e5c89e4eSSatish Balay 529793721a6SBarry Smith /*@C 530793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 531811af0c4SBarry Smith after `PetscInitialize()` is called but before `PetscFinalize()`. 532793721a6SBarry Smith 533793721a6SBarry Smith Not Collective 534793721a6SBarry Smith 5352fe279fdSBarry Smith Output Parameter: 536793721a6SBarry Smith . args - the command line arguments 537793721a6SBarry Smith 538793721a6SBarry Smith Level: intermediate 539793721a6SBarry Smith 54021532e8aSBarry Smith Note: 54121532e8aSBarry Smith This does NOT start with the program name and IS `NULL` terminated (final arg is void) 542793721a6SBarry Smith 54321532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`, `PetscInitialize()` 544793721a6SBarry Smith @*/ 545d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args) 546d71ae5a4SJacob Faibussowitsch { 547793721a6SBarry Smith PetscInt i, argc = PetscGlobalArgc; 548793721a6SBarry Smith 549793721a6SBarry Smith PetscFunctionBegin; 55039a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 5519371c9d4SSatish Balay if (!argc) { 5529371c9d4SSatish Balay *args = NULL; 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5549371c9d4SSatish Balay } 5559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(argc, args)); 5569566063dSJacob Faibussowitsch for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i])); 5572d747510SLisandro Dalcin (*args)[argc - 1] = NULL; 5583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 559793721a6SBarry Smith } 560793721a6SBarry Smith 561793721a6SBarry Smith /*@C 562811af0c4SBarry Smith PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()` 563793721a6SBarry Smith 564793721a6SBarry Smith Not Collective 565793721a6SBarry Smith 5662fe279fdSBarry Smith Output Parameter: 567793721a6SBarry Smith . args - the command line arguments 568793721a6SBarry Smith 569793721a6SBarry Smith Level: intermediate 570793721a6SBarry Smith 571db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()` 572793721a6SBarry Smith @*/ 573d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args) 574d71ae5a4SJacob Faibussowitsch { 57539a651e2SJacob Faibussowitsch PetscFunctionBegin; 57639a651e2SJacob Faibussowitsch if (args) { 577793721a6SBarry Smith PetscInt i = 0; 578793721a6SBarry Smith 5799566063dSJacob Faibussowitsch while (args[i]) PetscCall(PetscFree(args[i++])); 5809566063dSJacob Faibussowitsch PetscCall(PetscFree(args)); 58139a651e2SJacob Faibussowitsch } 5823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 583793721a6SBarry Smith } 584793721a6SBarry Smith 58527104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 58630befbd2SBarry Smith #include <petscconfiginfo.h> 58730befbd2SBarry Smith 588d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) 589d71ae5a4SJacob Faibussowitsch { 59027104ee2SJacob Faibussowitsch PetscFunctionBegin; 59111525c0dSBarry Smith if (!PetscGlobalRank) { 59230befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64]; 59311525c0dSBarry Smith int port; 594ffbd1cfbSBarry Smith PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE; 59511525c0dSBarry Smith size_t applinelen, introlen; 596ffbd1cfbSBarry Smith char sawsurl[256]; 59711525c0dSBarry Smith 5989566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg)); 59911525c0dSBarry Smith if (flg) { 60011525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 60111525c0dSBarry Smith 6029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL)); 60311525c0dSBarry Smith if (sawslog[0]) { 604792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog)); 60511525c0dSBarry Smith } else { 606792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL)); 60711525c0dSBarry Smith } 60811525c0dSBarry Smith } 6099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg)); 61048a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert)); 6119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL)); 612ffbd1cfbSBarry Smith if (selectport) { 613792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 614792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 615ffbd1cfbSBarry Smith } else { 6169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg)); 61748a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Port, (port)); 618ffbd1cfbSBarry Smith } 6199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg)); 62011525c0dSBarry Smith if (flg) { 621792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 6229566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(root, ".", &rootlocal)); 6239c1e0ce8SBarry Smith } else { 6249566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg)); 6259c1e0ce8SBarry Smith if (flg) { 6269566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root))); 627792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 6289c1e0ce8SBarry Smith } 62911525c0dSBarry Smith } 6309566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2)); 63111525c0dSBarry Smith if (flg2) { 63211525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 63328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option"); 6349566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root)); 6359566063dSJacob Faibussowitsch PetscCall(PetscTestDirectory(jsdir, 'r', &flg)); 63628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory"); 637792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Local_Header, ()); 63811525c0dSBarry Smith } 6399566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(programname, sizeof(programname))); 6409566063dSJacob Faibussowitsch PetscCall(PetscStrlen(help, &applinelen)); 64111525c0dSBarry Smith introlen = 4096 + applinelen; 64230a8c9c0SSurtai Han applinelen += 1024; 6439566063dSJacob Faibussowitsch PetscCall(PetscMalloc(applinelen, &appline)); 6449566063dSJacob Faibussowitsch PetscCall(PetscMalloc(introlen, &intro)); 64511525c0dSBarry Smith 64611525c0dSBarry Smith if (rootlocal) { 6479566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname)); 6489566063dSJacob Faibussowitsch PetscCall(PetscTestFile(appline, 'r', &rootlocal)); 64911525c0dSBarry Smith } 6509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetAll(NULL, &options)); 65111525c0dSBarry Smith if (rootlocal && help) { 6529566063dSJacob 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)); 65311525c0dSBarry Smith } else if (help) { 6549566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help)); 65511525c0dSBarry Smith } else { 6569566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options)); 65711525c0dSBarry Smith } 6589566063dSJacob Faibussowitsch PetscCall(PetscFree(options)); 6599566063dSJacob Faibussowitsch PetscCall(PetscGetVersion(version, sizeof(version))); 6609371c9d4SSatish Balay PetscCall(PetscSNPrintf(intro, introlen, 6619371c9d4SSatish Balay "<body>\n" 662a17b96a8SKyle 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" 663df62c222SBarry 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" 6649371c9d4SSatish Balay "%s", 6659371c9d4SSatish Balay version, petscconfigureoptions, appline)); 666792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro)); 6679566063dSJacob Faibussowitsch PetscCall(PetscFree(intro)); 6689566063dSJacob Faibussowitsch PetscCall(PetscFree(appline)); 669ffbd1cfbSBarry Smith if (selectport) { 670aa573868SBarry Smith PetscBool silent; 6717d812c46SBarry Smith 6727d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 67339a651e2SJacob Faibussowitsch while (SAWs_Initialize()) { 674792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 675792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 6767d812c46SBarry Smith } 6777d812c46SBarry Smith 6789566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL)); 679aa573868SBarry Smith if (!silent) { 680792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl)); 6819566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl)); 682ffbd1cfbSBarry Smith } 6837d812c46SBarry Smith } else { 684792fecdfSBarry Smith PetscCallSAWs(SAWs_Initialize, ()); 685aa573868SBarry Smith } 6869566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister("@TechReport{ saws,\n" 6870af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 6880af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 6890af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 6909371c9d4SSatish Balay " Year = 2013\n}\n", 6919371c9d4SSatish Balay NULL)); 69211525c0dSBarry Smith } 6933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69411525c0dSBarry Smith } 69511525c0dSBarry Smith #endif 69611525c0dSBarry Smith 6974dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 698d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) 699d71ae5a4SJacob Faibussowitsch { 7004dfee713SSatish Balay PetscFunctionBegin; 7014dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 7024dfee713SSatish Balay /* see MPI.py for details on this bug */ 7034dfee713SSatish Balay (void)setenv("HWLOC_COMPONENTS", "-x86", 1); 7044dfee713SSatish Balay #endif 7053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7064dfee713SSatish Balay } 7074dfee713SSatish Balay 708a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS) 709a56f64adSBarry Smith #include <adios.h> 71022580e64SBarry Smith #include <adios_read.h> 7117b56e58cSSatish Balay int64_t Petsc_adios_group; 712a56f64adSBarry Smith #endif 713a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP) 714cd1b99a6SBarry Smith #include <omp.h> 715f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 716cd1b99a6SBarry Smith #endif 717a56f64adSBarry Smith 718a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h> 719a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA) 7200e6b6b59SJacob Faibussowitsch #include <petscdevice_cuda.h> 721a4af0ceeSJacob Faibussowitsch // REMOVE ME 722a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL; 723a4af0ceeSJacob Faibussowitsch #endif 724a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP) 7250e6b6b59SJacob Faibussowitsch #include <petscdevice_hip.h> 726a4af0ceeSJacob Faibussowitsch // REMOVE ME 727a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL; 728a4af0ceeSJacob Faibussowitsch #endif 729a4af0ceeSJacob Faibussowitsch 73027104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H) 731bc8a24ecSBarry Smith #include <dlfcn.h> 732bc8a24ecSBarry Smith #endif 733a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void); 734a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 7353274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscViennaCLInit(void); 736a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE; 737a4af0ceeSJacob Faibussowitsch #endif 738bc8a24ecSBarry Smith 739660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE; 740660278c0SBarry Smith 74185649d77SJunchao Zhang /* 74285649d77SJunchao Zhang PetscInitialize_Common - shared code between C and Fortran initialization 74385649d77SJunchao Zhang 74485649d77SJunchao Zhang prog: program name 74502101c96SSatish Balay file: optional PETSc database file name. Might be in Fortran string format when 'ftn' is true 74685649d77SJunchao Zhang help: program help message 747da81f932SPierre Jolivet ftn: is it called from Fortran initialization (petscinitializef_)? 74885649d77SJunchao Zhang readarguments,len: used when fortran is true 74985649d77SJunchao Zhang */ 750d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len) 751d71ae5a4SJacob Faibussowitsch { 75285649d77SJunchao Zhang PetscMPIInt size; 75385649d77SJunchao Zhang PetscBool flg = PETSC_TRUE; 75485649d77SJunchao Zhang char hostname[256]; 755046ed546SBarry Smith PetscBool blas_view_flag = PETSC_FALSE; 75685649d77SJunchao Zhang 75727104ee2SJacob Faibussowitsch PetscFunctionBegin; 7583ba16761SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS); 75939a651e2SJacob Faibussowitsch /* these must be initialized in a routine, not as a constant declaration */ 76039a651e2SJacob Faibussowitsch PETSC_STDOUT = stdout; 76139a651e2SJacob Faibussowitsch PETSC_STDERR = stderr; 76239a651e2SJacob Faibussowitsch 7639566063dSJacob Faibussowitsch /* PetscCall can be used from now */ 76439a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_TRUE; 76539a651e2SJacob Faibussowitsch 76685649d77SJunchao Zhang /* 76785649d77SJunchao Zhang The checking over compatible runtime libraries is complicated by the MPI ABI initiative 76885649d77SJunchao Zhang https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 76985649d77SJunchao Zhang MPICH v3.1 (Released February 2014) 77085649d77SJunchao Zhang IBM MPI v2.1 (December 2014) 77185649d77SJunchao Zhang Intel MPI Library v5.0 (2014) 77285649d77SJunchao Zhang Cray MPT v7.0.0 (June 2014) 77385649d77SJunchao Zhang As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 77485649d77SJunchao Zhang listed above and since that time are compatible. 77585649d77SJunchao Zhang 77685649d77SJunchao Zhang Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 77785649d77SJunchao Zhang at compile time or runtime. Thus we will need to systematically track the allowed versions 77885649d77SJunchao Zhang and how they are represented in the mpi.h and MPI_Get_library_version() output in order 77985649d77SJunchao Zhang to perform the checking. 78085649d77SJunchao Zhang 78185649d77SJunchao Zhang Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 78285649d77SJunchao Zhang 78385649d77SJunchao Zhang Questions: 78485649d77SJunchao Zhang 78585649d77SJunchao Zhang Should the checks for ABI incompatibility be only on the major version number below? 78685649d77SJunchao Zhang Presumably the output to stderr will be removed before a release. 78785649d77SJunchao Zhang */ 78885649d77SJunchao Zhang 78985649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 79085649d77SJunchao Zhang { 79185649d77SJunchao Zhang char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 79285649d77SJunchao Zhang PetscMPIInt mpilibraryversionlength; 79339a651e2SJacob Faibussowitsch 7949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength)); 79585649d77SJunchao Zhang /* check for MPICH versions before MPI ABI initiative */ 79685649d77SJunchao Zhang #if defined(MPICH_VERSION) 79785649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000 79885649d77SJunchao Zhang { 79985649d77SJunchao Zhang char *ver, *lf; 80085649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 80139a651e2SJacob Faibussowitsch 8029566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver)); 80339a651e2SJacob Faibussowitsch if (ver) { 8049566063dSJacob Faibussowitsch PetscCall(PetscStrchr(ver, '\n', &lf)); 80539a651e2SJacob Faibussowitsch if (lf) { 80685649d77SJunchao Zhang *lf = 0; 8079566063dSJacob Faibussowitsch PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg)); 80885649d77SJunchao Zhang } 80985649d77SJunchao Zhang } 81085649d77SJunchao Zhang if (!flg) { 811d5b396e8SSatish Balay PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION)); 81285649d77SJunchao Zhang flg = PETSC_TRUE; 81385649d77SJunchao Zhang } 81485649d77SJunchao Zhang } 81585649d77SJunchao Zhang #endif 81685649d77SJunchao Zhang /* check for Open MPI version, it is not part of the MPI ABI initiative (is it part of another initiative that needs to be handled?) */ 81785649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION) 81885649d77SJunchao Zhang { 81985649d77SJunchao Zhang char *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf; 82085649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 82185649d77SJunchao Zhang #define PSTRSZ 2 82285649d77SJunchao Zhang char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"}; 82385649d77SJunchao Zhang char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "}; 82485649d77SJunchao Zhang int i; 82585649d77SJunchao Zhang for (i = 0; i < PSTRSZ; i++) { 8269566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver)); 82739a651e2SJacob Faibussowitsch if (ver) { 8289566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION)); 8299566063dSJacob Faibussowitsch PetscCall(PetscStrstr(ver, bs, &bsf)); 83039a651e2SJacob Faibussowitsch if (bsf) flg = PETSC_TRUE; 83185649d77SJunchao Zhang break; 83285649d77SJunchao Zhang } 83385649d77SJunchao Zhang } 83485649d77SJunchao Zhang if (!flg) { 8353ba16761SJacob Faibussowitsch PetscCall(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)); 83685649d77SJunchao Zhang flg = PETSC_TRUE; 83785649d77SJunchao Zhang } 83885649d77SJunchao Zhang } 83985649d77SJunchao Zhang #endif 84085649d77SJunchao Zhang } 84185649d77SJunchao Zhang #endif 84285649d77SJunchao Zhang 843e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__)) 84485649d77SJunchao Zhang /* These symbols are currently in the Open MPI and MPICH libraries; they may not always be, in that case the test will simply not detect the problem */ 84539a651e2SJacob 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 Open MPI and MPICH based MPI libraries and will not run correctly"); 84685649d77SJunchao Zhang #endif 84785649d77SJunchao Zhang 84885649d77SJunchao Zhang /* on Windows - set printf to default to printing 2 digit exponents */ 84985649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT) 85085649d77SJunchao Zhang _set_output_format(_TWO_DIGIT_EXPONENT); 85185649d77SJunchao Zhang #endif 85285649d77SJunchao Zhang 8539566063dSJacob Faibussowitsch PetscCall(PetscOptionsCreateDefault()); 85485649d77SJunchao Zhang 85585649d77SJunchao Zhang PetscFinalizeCalled = PETSC_FALSE; 85685649d77SJunchao Zhang 8579566063dSJacob Faibussowitsch PetscCall(PetscSetProgramName(prog)); 8589566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen)); 8599566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout)); 8609566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr)); 8619566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscCommSpinLock)); 86285649d77SJunchao Zhang 86385649d77SJunchao Zhang if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 8649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN)); 86585649d77SJunchao Zhang 86685649d77SJunchao Zhang if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) { 8679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS)); 8689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE)); 86985649d77SJunchao Zhang } 87085649d77SJunchao Zhang 87185649d77SJunchao Zhang /* Done after init due to a bug in MPICH-GM? */ 8729566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 87385649d77SJunchao Zhang 8749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank)); 8759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize)); 87685649d77SJunchao Zhang 87785649d77SJunchao Zhang MPIU_BOOL = MPI_INT; 87885649d77SJunchao Zhang MPIU_ENUM = MPI_INT; 87985649d77SJunchao Zhang MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64; 88085649d77SJunchao Zhang if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 88185649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 88285649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG) 88385649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 88485649d77SJunchao Zhang #endif 88539a651e2SJacob Faibussowitsch else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t"); 88685649d77SJunchao Zhang 88785649d77SJunchao Zhang /* 88885649d77SJunchao Zhang Initialized the global complex variable; this is because with 88985649d77SJunchao Zhang shared libraries the constructors for global variables 89085649d77SJunchao Zhang are not called; at least on IRIX. 89185649d77SJunchao Zhang */ 89285649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX) 89385649d77SJunchao Zhang { 89485649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 89585649d77SJunchao Zhang PetscComplex ic(0.0, 1.0); 89685649d77SJunchao Zhang PETSC_i = ic; 89785649d77SJunchao Zhang #else 89885649d77SJunchao Zhang PETSC_i = _Complex_I; 89985649d77SJunchao Zhang #endif 90085649d77SJunchao Zhang } 90185649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */ 90285649d77SJunchao Zhang 90385649d77SJunchao Zhang /* 90485649d77SJunchao Zhang Create the PETSc MPI reduction operator that sums of the first 90585649d77SJunchao Zhang half of the entries and maxes the second half. 90685649d77SJunchao Zhang */ 9079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP)); 90885649d77SJunchao Zhang 909*ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 9109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128)); 9119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128)); 91274c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 9139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128)); 9149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128)); 91585649d77SJunchao Zhang #endif 91674c01ffaSSatish Balay #endif 917*ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16) 9189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16)); 9199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FP16)); 92085649d77SJunchao Zhang #endif 92185649d77SJunchao Zhang 92285649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 9239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM)); 924613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX)); 925613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN)); 926*ef384fdeSPierre Jolivet #elif (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16)) 9279e517322SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128)); 92885649d77SJunchao Zhang #endif 92985649d77SJunchao Zhang 9309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR)); 93162e5d2d2SJDBetteridge PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp)); 9329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR)); 93385649d77SJunchao Zhang 93485649d77SJunchao Zhang /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */ 93585649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI) 93685649d77SJunchao Zhang { 93785649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 93893d501b3SJacob Faibussowitsch MPI_Aint blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)}; 93985649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_REAL, MPIU_INT}, tmpStruct; 94085649d77SJunchao Zhang 9419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 94293d501b3SJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT)); 9439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT)); 94585649d77SJunchao Zhang } 94685649d77SJunchao Zhang { 94785649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 94893d501b3SJacob Faibussowitsch MPI_Aint blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)}; 94985649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_SCALAR, MPIU_INT}, tmpStruct; 95085649d77SJunchao Zhang 9519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 95293d501b3SJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT)); 9539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT)); 95585649d77SJunchao Zhang } 95685649d77SJunchao Zhang #endif 95785649d77SJunchao Zhang 95885649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES) 9599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT)); 9609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2INT)); 96185649d77SJunchao Zhang #endif 9629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT)); 9639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPI_4INT)); 9649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT)); 9659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_4INT)); 96685649d77SJunchao Zhang 96785649d77SJunchao Zhang /* 96885649d77SJunchao Zhang Attributes to be set on PETSc communicators 96985649d77SJunchao Zhang */ 9708434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_DeleteFn, &Petsc_Counter_keyval, (void *)0)); 9718434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_DeleteFn, &Petsc_InnerComm_keyval, (void *)0)); 9728434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_DeleteFn, &Petsc_OuterComm_keyval, (void *)0)); 9738434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_DeleteFn, &Petsc_ShmComm_keyval, (void *)0)); 97462e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, (void *)0)); 97562e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, (void *)0)); 97685649d77SJunchao Zhang 977fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS) 9789566063dSJacob Faibussowitsch if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len)); 97985649d77SJunchao Zhang else 98085649d77SJunchao Zhang #endif 9819566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file)); 98285649d77SJunchao Zhang 98385649d77SJunchao Zhang /* call a second time so it can look in the options database */ 9849566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 98585649d77SJunchao Zhang 98685649d77SJunchao Zhang /* 98785649d77SJunchao Zhang Check system options and print help 98885649d77SJunchao Zhang */ 9899566063dSJacob Faibussowitsch PetscCall(PetscOptionsCheckInitial_Private(help)); 99085649d77SJunchao Zhang 991a4af0ceeSJacob Faibussowitsch /* 9920e6b6b59SJacob Faibussowitsch Creates the logging data structures; this is enabled even if logging is not turned on 9930e6b6b59SJacob Faibussowitsch This is the last thing we do before returning to the user code to prevent having the 9940e6b6b59SJacob Faibussowitsch logging numbers contaminated by any startup time associated with MPI 9950e6b6b59SJacob Faibussowitsch */ 9960e6b6b59SJacob Faibussowitsch PetscCall(PetscLogInitialize()); 9970e6b6b59SJacob Faibussowitsch 9980e6b6b59SJacob Faibussowitsch /* 999a4af0ceeSJacob Faibussowitsch Initialize PetscDevice and PetscDeviceContext 1000a4af0ceeSJacob Faibussowitsch 1001a4af0ceeSJacob Faibussowitsch Note to any future devs thinking of moving this, proper initialization requires: 1002a4af0ceeSJacob Faibussowitsch 1. MPI initialized 1003a4af0ceeSJacob Faibussowitsch 2. Options DB initialized 10040e6b6b59SJacob Faibussowitsch 3. Petsc error handling initialized, specifically signal handlers. This expects to set up 10050e6b6b59SJacob Faibussowitsch its own SIGSEV handler via the push/pop interface. 10060e6b6b59SJacob Faibussowitsch 4. Logging initialized 1007a4af0ceeSJacob Faibussowitsch */ 10089566063dSJacob Faibussowitsch PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD)); 1009a4af0ceeSJacob Faibussowitsch 1010a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 1011a4af0ceeSJacob Faibussowitsch flg = PETSC_FALSE; 1012d75802c7SJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg)); 10139566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL)); 1014a4af0ceeSJacob Faibussowitsch PetscViennaCLSynchronize = flg; 10159566063dSJacob Faibussowitsch PetscCall(PetscViennaCLInit()); 1016a4af0ceeSJacob Faibussowitsch #endif 1017a4af0ceeSJacob Faibussowitsch 10189566063dSJacob Faibussowitsch PetscCall(PetscCitationsInitialize()); 101985649d77SJunchao Zhang 102085649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS) 10219566063dSJacob Faibussowitsch PetscCall(PetscInitializeSAWs(ftn ? NULL : help)); 102227104ee2SJacob Faibussowitsch flg = PETSC_FALSE; 10239566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg)); 10249566063dSJacob Faibussowitsch if (flg) PetscCall(PetscStackViewSAWs()); 102585649d77SJunchao Zhang #endif 102685649d77SJunchao Zhang 102785649d77SJunchao Zhang /* 102885649d77SJunchao Zhang Load the dynamic libraries (on machines that support them), this registers all 102985649d77SJunchao Zhang the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 103085649d77SJunchao Zhang */ 10319566063dSJacob Faibussowitsch PetscCall(PetscInitialize_DynamicLibraries()); 103285649d77SJunchao Zhang 10339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 10349566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size)); 103596dcfd6cSBarry Smith PetscCall(PetscGetHostName(hostname, sizeof(hostname))); 10369566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname)); 103785649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP) 103885649d77SJunchao Zhang { 103985649d77SJunchao Zhang PetscBool omp_view_flag; 104085649d77SJunchao Zhang char *threads = getenv("OMP_NUM_THREADS"); 104185649d77SJunchao Zhang 104285649d77SJunchao Zhang if (threads) { 10439566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads)); 104485649d77SJunchao Zhang (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads); 104585649d77SJunchao Zhang } else { 10462f840973SStefano Zampini PetscNumOMPThreads = (PetscInt)omp_get_max_threads(); 10479566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads)); 104885649d77SJunchao Zhang } 1049d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys"); 10509566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg)); 10519566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag)); 1052d0609cedSBarry Smith PetscOptionsEnd(); 105385649d77SJunchao Zhang if (flg) { 10541a0ee928SPierre Jolivet PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads)); 105585649d77SJunchao Zhang omp_set_num_threads((int)PetscNumOMPThreads); 105685649d77SJunchao Zhang } 105748a46eb9SPierre Jolivet if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads)); 105885649d77SJunchao Zhang } 105985649d77SJunchao Zhang #endif 106085649d77SJunchao Zhang 1061046ed546SBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "BLAS options", "Sys"); 1062046ed546SBarry Smith PetscCall(PetscOptionsName("-blas_view", "Display number of threads to use for BLAS operations", NULL, &blas_view_flag)); 1063046ed546SBarry Smith #if defined(PETSC_HAVE_BLI_THREAD_SET_NUM_THREADS) || defined(PETSC_HAVE_MKL_SET_NUM_THREADS) || defined(PETSC_HAVE_OPENBLAS_SET_NUM_THREADS) 1064046ed546SBarry Smith { 1065046ed546SBarry Smith char *threads = NULL; 1066046ed546SBarry Smith 1067046ed546SBarry Smith /* determine any default number of threads requested in the environment; TODO: Apple libraries? */ 1068046ed546SBarry Smith #if defined(PETSC_HAVE_BLI_THREAD_SET_NUM_THREADS) 1069046ed546SBarry Smith threads = getenv("BLIS_NUM_THREADS"); 1070046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of BLIS threads %s given by BLIS_NUM_THREADS\n", threads)); 1071046ed546SBarry Smith if (!threads) { 1072046ed546SBarry Smith threads = getenv("OMP_NUM_THREADS"); 1073046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of BLIS threads %s given by OMP_NUM_THREADS\n", threads)); 1074046ed546SBarry Smith } 1075046ed546SBarry Smith #elif defined(PETSC_HAVE_MKL_SET_NUM_THREADS) 1076046ed546SBarry Smith threads = getenv("MKL_NUM_THREADS"); 1077046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of MKL threads %s given by MKL_NUM_THREADS\n", threads)); 1078046ed546SBarry Smith #elif defined(PETSC_HAVE_OPENBLAS_SET_NUM_THREADS) 1079046ed546SBarry Smith threads = getenv("OPENBLAS_NUM_THREADS"); 1080046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of OpenBLAS threads %s given by OPENBLAS_NUM_THREADS\n", threads)); 1081046ed546SBarry Smith if (!threads) { 1082046ed546SBarry Smith threads = getenv("OMP_NUM_THREADS"); 1083046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of OpenBLAS threads %s given by OMP_NUM_THREADS\n", threads)); 1084046ed546SBarry Smith } 1085046ed546SBarry Smith #endif 1086046ed546SBarry Smith if (threads) (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumBLASThreads); 1087046ed546SBarry Smith PetscCall(PetscOptionsInt("-blas_num_threads", "Number of threads to use for BLAS operations", "None", PetscNumBLASThreads, &PetscNumBLASThreads, &flg)); 1088046ed546SBarry Smith PetscCall(PetscBLASSetNumThreads(PetscNumBLASThreads)); 1089046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: number of threads %" PetscInt_FMT "\n", PetscNumBLASThreads)); 1090046ed546SBarry Smith } 1091046ed546SBarry Smith #elif defined(PETSC_HAVE_APPLE_ACCELERATE) 1092046ed546SBarry Smith PetscCall(PetscInfo(NULL, "BLAS: Apple Accelerate library, thread support with no user control\n")); 1093046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: Apple Accelerate library, thread support with no user control\n")); 1094046ed546SBarry Smith #else 1095046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: no thread support\n")); 1096046ed546SBarry Smith #endif 1097046ed546SBarry Smith PetscOptionsEnd(); 1098046ed546SBarry Smith 109985649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 110085649d77SJunchao Zhang /* 110185649d77SJunchao Zhang Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 110285649d77SJunchao Zhang 110385649d77SJunchao Zhang Currently not used because it is not supported by MPICH. 110485649d77SJunchao Zhang */ 11059566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL)); 110685649d77SJunchao Zhang #endif 110785649d77SJunchao Zhang 110885649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS) 11099566063dSJacob Faibussowitsch PetscCall(PetscFPTCreate(10000)); 111085649d77SJunchao Zhang #endif 111185649d77SJunchao Zhang 111285649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC) 111385649d77SJunchao Zhang { 111485649d77SJunchao Zhang PetscViewer viewer; 11159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg)); 111685649d77SJunchao Zhang if (flg) { 11179566063dSJacob Faibussowitsch PetscCall(PetscProcessPlacementView(viewer)); 1118cd791dc2SBarry Smith PetscCall(PetscOptionsRestoreViewer(&viewer)); 111985649d77SJunchao Zhang } 112085649d77SJunchao Zhang } 112185649d77SJunchao Zhang #endif 112285649d77SJunchao Zhang 112385649d77SJunchao Zhang flg = PETSC_TRUE; 11249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL)); 11259566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE)); 112685649d77SJunchao Zhang 112785649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS) 11283ba16761SJacob Faibussowitsch PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD); 11293ba16761SJacob Faibussowitsch PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default); 11303ba16761SJacob Faibussowitsch PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", ""); 11313ba16761SJacob Faibussowitsch PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, ""); 113285649d77SJunchao Zhang #endif 113385649d77SJunchao Zhang 113485649d77SJunchao Zhang #if defined(__VALGRIND_H) 113585649d77SJunchao Zhang PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE; 113685649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE) 1137337bb527SBarry Smith 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, try configuring with --download-fblaslapack or --download-f2cblaslapack")); 113885649d77SJunchao Zhang #endif 113985649d77SJunchao Zhang #endif 114085649d77SJunchao Zhang /* 114185649d77SJunchao Zhang Set flag that we are completely initialized 114285649d77SJunchao Zhang */ 114385649d77SJunchao Zhang PetscInitializeCalled = PETSC_TRUE; 114485649d77SJunchao Zhang 11459566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg)); 11469566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPythonInitialize(NULL, NULL)); 1147f1f2ae84SBarry Smith 1148f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1149f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin()); 1150f1f2ae84SBarry Smith else PetscCheck(!flg, PETSC_COMM_WORLD, PETSC_ERR_SUP, "PETSc configured using -with-single-library=0; -mpi_linear_solver_server not supported in that case"); 11513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115285649d77SJunchao Zhang } 115385649d77SJunchao Zhang 115410450e9eSJacob Faibussowitsch // "Unknown section 'Environmental Variables'" 115510450e9eSJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-section-header-unknown 1156e5c89e4eSSatish Balay /*@C 1157e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 1158811af0c4SBarry Smith `PetscInitialize()` calls MPI_Init() if that has yet to be called, 1159e5c89e4eSSatish Balay so this routine should always be called near the beginning of 1160e5c89e4eSSatish Balay your program -- usually the very first line! 1161e5c89e4eSSatish Balay 1162811af0c4SBarry Smith Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set 1163e5c89e4eSSatish Balay 1164e5c89e4eSSatish Balay Input Parameters: 1165e5c89e4eSSatish Balay + argc - count of number of command line arguments 1166e5c89e4eSSatish Balay . args - the command line arguments 1167be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format. 1168be10d61cSLisandro Dalcin Use NULL or empty string to not check for code specific file. 1169be10d61cSLisandro Dalcin Also checks ~/.petscrc, .petscrc and petscrc. 1170c5b5d8d5SVaclav Hapla Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files. 11710298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message 1172e5c89e4eSSatish Balay 1173811af0c4SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that 1174811af0c4SBarry Smith communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a 1175811af0c4SBarry Smith four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not, 1176811af0c4SBarry Smith then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even 117705827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 1178e5c89e4eSSatish Balay 1179e5c89e4eSSatish Balay Options Database Keys: 11807ca660e7SBarry Smith + -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message 11817ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 1182e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 11838a690491SBarry Smith . -on_error_emacs <machinename> - causes emacsclient to jump to error file 1184811af0c4SBarry Smith . -on_error_abort - calls `abort()` when error detected (no traceback) 1185811af0c4SBarry Smith . -on_error_mpiabort - calls `MPI_abort()` when error detected 11861af3601dSBarry Smith . -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr 11878a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 1188bf4d2887SBarry Smith . -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger 1189e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds) - Pauses debugger 1190e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 1191e5c89e4eSSatish Balay process and wait (-debugger_pause) seconds for attachment 1192aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 119392f119d6SBarry 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 1194811af0c4SBarry Smith . -malloc_view - show a list of all allocated memory during `PetscFinalize()` 119592f119d6SBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view 1196608c71bfSMatthew G. Knepley . -malloc_requested_size - malloc logging will record the requested size rather than size after alignment 119792f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 1198e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 1199e5c89e4eSSatish Balay . -shared_tmp - indicates /tmp directory is shared by all processors 1200e5c89e4eSSatish Balay . -not_shared_tmp - each processor has own /tmp 1201e5c89e4eSSatish Balay . -tmp - alternative name of /tmp directory 1202e5c89e4eSSatish Balay . -get_total_flops - returns total flops done by all processors 12030841954dSBarry Smith - -memory_view - Print memory usage at end of run 1204e5c89e4eSSatish Balay 1205c5b5d8d5SVaclav Hapla Options Database Keys for Option Database: 1206c5b5d8d5SVaclav Hapla + -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc 1207c5b5d8d5SVaclav Hapla . -options_monitor - monitor all set options to standard output for the whole program run 1208811af0c4SBarry Smith - -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()` 1209c5b5d8d5SVaclav Hapla 1210c5b5d8d5SVaclav Hapla Options -options_monitor_{all,cancel} are 1211c5b5d8d5SVaclav Hapla position-independent and apply to all options set since the PETSc start. 1212c5b5d8d5SVaclav Hapla They can be used also in option files. 1213c5b5d8d5SVaclav Hapla 1214811af0c4SBarry Smith See `PetscOptionsMonitorSet()` to do monitoring programmatically. 1215c5b5d8d5SVaclav Hapla 1216e5c89e4eSSatish Balay Options Database Keys for Profiling: 1217a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 1218811af0c4SBarry Smith + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`. 1219217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 1220217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 1221495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 1222811af0c4SBarry Smith hangs without running in the debugger). See `PetscLogTraceBegin()`. 1223ad2e3d55SToby Isaac . -log_view [:filename:format][,[:filename:format]...] - Prints summary of flop and timing information to screen or file, see `PetscLogView()` (up to 4 viewers) 1224811af0c4SBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`. 1225811af0c4SBarry Smith . -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView(). 1226f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 1227b665b14eSToby Isaac . -log [filename] - Logs profiling information in a dump file, see `PetscLogDump()`. 1228b665b14eSToby Isaac . -log_all [filename] - Same as `-log`. 1229c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 12302611ad71SToby Isaac . -log_perfstubs - Starts a log handler with the perfstubs interface (which is used by TAU) 123161cc7448SToby Isaac . -log_nvtx - Starts an nvtx log handler for use with Nsight 1232811af0c4SBarry Smith . -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off 1233c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 1234495fc317SBarry Smith 1235ffbd1cfbSBarry Smith Options Database Keys for SAWs: 1236ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 1237ffbd1cfbSBarry 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 1238ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 1239ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 1240ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 1241ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 1242ffbd1cfbSBarry Smith 1243e5c89e4eSSatish Balay Environmental Variables: 1244811af0c4SBarry Smith + `PETSC_TMP` - alternative tmp directory 1245811af0c4SBarry Smith . `PETSC_SHARED_TMP` - tmp is shared by all processes 1246811af0c4SBarry Smith . `PETSC_NOT_SHARED_TMP` - each process has its own private tmp 1247811af0c4SBarry Smith . `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs 1248811af0c4SBarry Smith . `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document 1249811af0c4SBarry Smith . `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer 1250811af0c4SBarry Smith - `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to 1251e5c89e4eSSatish Balay 1252e5c89e4eSSatish Balay Level: beginner 1253e5c89e4eSSatish Balay 1254811af0c4SBarry Smith Note: 1255811af0c4SBarry Smith If for some reason you must call `MPI_Init()` separately, call 1256811af0c4SBarry Smith it before `PetscInitialize()`. 1257e5c89e4eSSatish Balay 1258811af0c4SBarry Smith Fortran Notes: 1259811af0c4SBarry Smith In Fortran this routine can be called with 1260811af0c4SBarry Smith .vb 1261811af0c4SBarry Smith call PetscInitialize(ierr) 1262811af0c4SBarry Smith call PetscInitialize(file,ierr) or 1263811af0c4SBarry Smith call PetscInitialize(file,help,ierr) 1264811af0c4SBarry Smith .ve 1265e5c89e4eSSatish Balay 1266811af0c4SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after 1267811af0c4SBarry Smith calling `PetscInitialize()`. 1268e5c89e4eSSatish Balay 12695dedd997SBarry Smith Options Database Key for Developers: 12705dedd997SBarry Smith . -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form: 12715dedd997SBarry Smith "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py) 12725dedd997SBarry Smith 1273db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()` 1274e5c89e4eSSatish Balay @*/ 1275d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[]) 1276d71ae5a4SJacob Faibussowitsch { 127785649d77SJunchao Zhang PetscMPIInt flag; 127821ef0414SBarry Smith const char *prog = "Unknown Name", *mpienv; 1279e5c89e4eSSatish Balay 128027104ee2SJacob Faibussowitsch PetscFunctionBegin; 12813ba16761SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS); 12829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Initialized(&flag)); 1283e5c89e4eSSatish Balay if (!flag) { 128439a651e2SJacob 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"); 12859566063dSJacob Faibussowitsch PetscCall(PetscPreMPIInit_Private()); 12865e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 12875e765c61SJed Brown { 1288835d5ba2SJunchao Zhang PetscMPIInt provided; 1289835d5ba2SJunchao Zhang PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE ? MPI_THREAD_FUNNELED : PETSC_MPI_THREAD_REQUIRED, &provided)); 1290835d5ba2SJunchao Zhang PetscCheck(PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE || provided >= PETSC_MPI_THREAD_REQUIRED, PETSC_COMM_SELF, PETSC_ERR_MPI, "The MPI implementation's provided thread level is less than what you required"); 1291835d5ba2SJunchao Zhang if (PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE) PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED; // assign it a valid value after check-up 12925e765c61SJed Brown } 12935e765c61SJed Brown #else 12949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init(argc, args)); 12955e765c61SJed Brown #endif 129621ef0414SBarry Smith if (PetscDefined(HAVE_MPIUNI)) { 129721ef0414SBarry Smith mpienv = getenv("PMI_SIZE"); 129821ef0414SBarry Smith if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE"); 129921ef0414SBarry Smith if (mpienv) { 130021ef0414SBarry Smith PetscInt isize; 130121ef0414SBarry Smith PetscCall(PetscOptionsStringToInt(mpienv, &isize)); 130221ef0414SBarry Smith if (isize != 1) printf("You are using an MPI-uni (sequential) install of PETSc but trying to launch parallel jobs; you need full MPI version of PETSc\n"); 130321ef0414SBarry Smith PetscCheck(isize == 1, MPI_COMM_SELF, PETSC_ERR_MPI, "You are using an MPI-uni (sequential) install of PETSc but trying to launch parallel jobs; you need full MPI version of PETSc"); 130421ef0414SBarry Smith } 130521ef0414SBarry Smith } 1306e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 1307e5c89e4eSSatish Balay } 13084dfee713SSatish Balay 130985649d77SJunchao Zhang if (argc && *argc) prog = **args; 1310e5c89e4eSSatish Balay if (argc && args) { 1311e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 1312e5c89e4eSSatish Balay PetscGlobalArgs = *args; 1313e5c89e4eSSatish Balay } 1314811af0c4SBarry Smith PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0)); 13153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1316e5c89e4eSSatish Balay } 1317e5c89e4eSSatish Balay 131895c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1319ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1320ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 132195c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 1322e5c89e4eSSatish Balay 1323008a6e76SBarry Smith /* 1324008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1325008a6e76SBarry Smith */ 1326d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void) 1327d71ae5a4SJacob Faibussowitsch { 1328008a6e76SBarry Smith PetscFunctionBegin; 1329*ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 13309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128)); 133174c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 13329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128)); 1333008a6e76SBarry Smith #endif 133474c01ffaSSatish Balay #endif 1335*ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16) 13369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FP16)); 1337008a6e76SBarry Smith #endif 1338008a6e76SBarry Smith 1339de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 13409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_SUM)); 1341613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MAX)); 1342613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MIN)); 1343*ef384fdeSPierre Jolivet #elif (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16)) 13449e517322SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128)); 1345008a6e76SBarry Smith #endif 1346008a6e76SBarry Smith 13479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR)); 13489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT)); 13499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT)); 135040df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 13519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2INT)); 1352008a6e76SBarry Smith #endif 13539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPI_4INT)); 13549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_4INT)); 13559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP)); 135662e5d2d2SJDBetteridge PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp)); 13573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1358008a6e76SBarry Smith } 1359008a6e76SBarry Smith 1360a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void); 1361a4af0ceeSJacob Faibussowitsch 1362e5c89e4eSSatish Balay /*@C 1363e5c89e4eSSatish Balay PetscFinalize - Checks for options to be called at the conclusion 1364811af0c4SBarry Smith of the program. `MPI_Finalize()` is called only if the user had not 1365811af0c4SBarry Smith called `MPI_Init()` before calling `PetscInitialize()`. 1366e5c89e4eSSatish Balay 1367811af0c4SBarry Smith Collective on `PETSC_COMM_WORLD` 1368e5c89e4eSSatish Balay 1369e5c89e4eSSatish Balay Options Database Keys: 1370811af0c4SBarry Smith + -options_view - Calls `PetscOptionsView()` 1371e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 13727eb1d149SBarry 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 1373e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 1374811af0c4SBarry Smith . -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed 13754bd3d7f8SBarry Smith . -memory_view - Prints total memory usage 13764bd3d7f8SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and in what functions 1377e5c89e4eSSatish Balay 1378e5c89e4eSSatish Balay Level: beginner 1379e5c89e4eSSatish Balay 1380e5c89e4eSSatish Balay Note: 1381811af0c4SBarry Smith See `PetscInitialize()` for other runtime options. 1382e5c89e4eSSatish Balay 1383db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()` 1384e5c89e4eSSatish Balay @*/ 1385d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void) 1386d71ae5a4SJacob Faibussowitsch { 13874bb5149bSJed Brown PetscMPIInt rank; 1388a8d2bbe5SBarry Smith PetscInt nopt; 13892bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE; 1390dff31646SBarry Smith PetscBool flg; 139110463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 1392e5c89e4eSSatish Balay 13933cbbc5ffSBarry Smith PetscFunctionBegin; 139439a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()"); 13959566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PetscFinalize() called\n")); 1396b022a5c1SBarry Smith 1397f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1398f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd()); 1399f1f2ae84SBarry Smith 140062e5d2d2SJDBetteridge /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */ 140162e5d2d2SJDBetteridge { 140262e5d2d2SJDBetteridge union 140362e5d2d2SJDBetteridge { 140462e5d2d2SJDBetteridge MPI_Comm comm; 140562e5d2d2SJDBetteridge void *ptr; 140662e5d2d2SJDBetteridge } ucomm; 140762e5d2d2SJDBetteridge PetscMPIInt flg; 140862e5d2d2SJDBetteridge void *tmp; 140962e5d2d2SJDBetteridge 141062e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 141162e5d2d2SJDBetteridge if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg)); 141262e5d2d2SJDBetteridge if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF)); 141362e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 141462e5d2d2SJDBetteridge if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg)); 141562e5d2d2SJDBetteridge if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD)); 141662e5d2d2SJDBetteridge } 141762e5d2d2SJDBetteridge 14189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1419a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 14203ba16761SJacob Faibussowitsch PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE); 14213ba16761SJacob Faibussowitsch PetscCallExternal(adios_finalize, rank); 1422a56f64adSBarry Smith #endif 14239566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg)); 1424dff31646SBarry Smith if (flg) { 14251f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 14261f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 14271f817a21SBarry Smith 14289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL)); 142948a46eb9SPierre Jolivet if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd)); 14309566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits)); 1431dff31646SBarry Smith cits[0] = 0; 14329566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits)); 14339566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n")); 14349566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 14359566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits)); 14369566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 14379566063dSJacob Faibussowitsch PetscCall(PetscFClose(PETSC_COMM_WORLD, fd)); 14389566063dSJacob Faibussowitsch PetscCall(PetscFree(cits)); 1439dff31646SBarry Smith } 14409566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&PetscCitationsList)); 1441dff31646SBarry Smith 14422d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 14439566063dSJacob Faibussowitsch PetscCall(PetscFPTDestroy()); 14442d53ad75SBarry Smith #endif 14452d53ad75SBarry Smith 1446e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1447dff31646SBarry Smith flg = PETSC_FALSE; 14489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL)); 14491baa6e33SBarry Smith if (flg) PetscCall(PetscOptionsSAWsDestroy()); 1450d5649816SBarry Smith #endif 1451d5649816SBarry Smith 1452681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1453681455b2SBarry Smith flg1 = PETSC_FALSE; 14549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL)); 1455681455b2SBarry Smith if (flg1) { 1456681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 1457a5af8288SJose E. Roman PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -15 Xvfb", "r", NULL)); 1458681455b2SBarry Smith } 1459681455b2SBarry Smith #endif 1460681455b2SBarry Smith 146167584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 14629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL)); 146348a46eb9SPierre Jolivet if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n")); 146467584ceeSBarry Smith #endif 1465e5c89e4eSSatish Balay 14662611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 146790d69ab7SBarry Smith flg1 = PETSC_FALSE; 14689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL)); 1469e5c89e4eSSatish Balay if (flg1) { 1470e5c89e4eSSatish Balay PetscLogDouble flops = 0; 14719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD)); 14729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops)); 1473e5c89e4eSSatish Balay } 14742611ad71SToby Isaac } 1475e5c89e4eSSatish Balay 14762611ad71SToby Isaac if (PetscDefined(USE_LOG) && PetscDefined(HAVE_MPE)) { 1477e5c89e4eSSatish Balay mname[0] = 0; 14789566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1)); 14792611ad71SToby Isaac if (flg1) PetscCall(PetscLogMPEDump(mname[0] ? mname : NULL)); 1480e5c89e4eSSatish Balay } 1481a297a907SKarl Rupp 1482c9903f8fSJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 1483c9903f8fSJunchao Zhang // Free petsc/kokkos stuff before the potentially non-null petsc default gpu stream is destroyed by PetscObjectRegisterDestroyAll 1484c9903f8fSJunchao Zhang if (PetscKokkosInitialized) { 1485c9903f8fSJunchao Zhang PetscCall(PetscKokkosFinalize_Private()); 1486c9903f8fSJunchao Zhang PetscKokkosInitialized = PETSC_FALSE; 1487c9903f8fSJunchao Zhang } 1488c9903f8fSJunchao Zhang #endif 1489c9903f8fSJunchao Zhang 14903048253cSJacob Faibussowitsch // Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 14919566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1492dd710f27SBarry Smith 14932611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 14949566063dSJacob Faibussowitsch PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE)); 14959566063dSJacob Faibussowitsch PetscCall(PetscLogViewFromOptions()); 14969566063dSJacob Faibussowitsch PetscCall(PetscOptionsPopGetViewerOff()); 1497473903fcSJunchao Zhang // It should be turned on with PetscLogGpuTime() and never turned off except in this place 1498473903fcSJunchao Zhang PetscLogGpuTimeFlag = PETSC_FALSE; 149987c3beb6SLisandro Dalcin 15003048253cSJacob Faibussowitsch // Free any objects created by the last block of code. 15019566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1502dd710f27SBarry Smith 1503dd710f27SBarry Smith mname[0] = 0; 15049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1)); 15059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2)); 15069566063dSJacob Faibussowitsch if (flg1 || flg2) PetscCall(PetscLogDump(mname)); 15072611ad71SToby Isaac } 150810463e74SBarry Smith 150990d69ab7SBarry Smith flg1 = PETSC_FALSE; 15109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL)); 15119566063dSJacob Faibussowitsch if (!flg1) PetscCall(PetscPopSignalHandler()); 151290d69ab7SBarry Smith flg1 = PETSC_FALSE; 15139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL)); 15141baa6e33SBarry Smith if (flg1) PetscCall(PetscMPIDump(stdout)); 151590d69ab7SBarry Smith flg1 = PETSC_FALSE; 151690d69ab7SBarry Smith flg2 = PETSC_FALSE; 15178bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 15189566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1)); 15199566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 15209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL)); 1521e4c476e2SSatish Balay 152257b1f488SBarry Smith if (flg2) { PetscCall(PetscOptionsView(NULL, PETSC_VIEWER_STDOUT_WORLD)); } 1523e5c89e4eSSatish Balay 1524e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 15259566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1)); 15269566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1)); 1527e5c89e4eSSatish Balay 152833fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 15299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1)); 153059f199a7SJed Brown if (!flg1) flg3 = PETSC_TRUE; 1531e5c89e4eSSatish Balay if (flg3) { 15323de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 153357b1f488SBarry Smith PetscCall(PetscOptionsView(NULL, PETSC_VIEWER_STDOUT_WORLD)); 1534e5c89e4eSSatish Balay } 15359566063dSJacob Faibussowitsch PetscCall(PetscOptionsAllUsed(NULL, &nopt)); 15363de2bfdfSBarry Smith if (nopt) { 15379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n")); 15389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n")); 15393de2bfdfSBarry Smith if (nopt == 1) { 15409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n")); 1541e5c89e4eSSatish Balay } else { 15429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt)); 1543e5c89e4eSSatish Balay } 15443de2bfdfSBarry Smith } else if (flg3 && flg1) { 15459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n")); 1546df12ba86SBarry Smith } 15479566063dSJacob Faibussowitsch PetscCall(PetscOptionsLeft(NULL)); 1548e5c89e4eSSatish Balay } 1549e5c89e4eSSatish Balay 1550e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1551d45a07a7SBarry Smith if (!PetscGlobalRank) { 15529566063dSJacob Faibussowitsch PetscCall(PetscStackSAWsViewOff()); 1553792fecdfSBarry Smith PetscCallSAWs(SAWs_Finalize, ()); 1554d45a07a7SBarry Smith } 1555ec957eceSBarry Smith #endif 1556ec957eceSBarry Smith 155710463e74SBarry Smith /* 1558dbc8283eSBarry Smith List all objects the user may have forgot to free 15592eff7a51SBarry Smith */ 15602611ad71SToby Isaac if (PetscDefined(USE_LOG) && PetscObjectsLog) { 15619566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 1562a64a8e02SBarry Smith if (flg1) { 1563a64a8e02SBarry Smith MPI_Comm local_comm; 15647eb1d149SBarry Smith char string[64]; 1565a64a8e02SBarry Smith 15669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL)); 1567afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 15689566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 15699566063dSJacob Faibussowitsch PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE)); 15709566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 15719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 15720a1571b3SBarry Smith } 157305df10baSBarry Smith } 15744097062eSBarry Smith 1575dbc8283eSBarry Smith PetscObjectsCounts = 0; 1576dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 15779566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscObjects)); 15782eff7a51SBarry Smith 157933f85c2fSBarry Smith /* 158033f85c2fSBarry Smith Destroy any packages that registered a finalize 158133f85c2fSBarry Smith */ 15829566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalizeAll()); 158333f85c2fSBarry Smith 15849566063dSJacob Faibussowitsch PetscCall(PetscLogFinalize()); 1585101409b8SToby Isaac 158633f85c2fSBarry Smith /* 158748dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 158848dd1dffSBarry Smith */ 15892e956fe4SStefano Zampini if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll()); 159037e93019SBarry Smith 15914028d114SSatish Balay if (petsc_history) { 15929566063dSJacob Faibussowitsch PetscCall(PetscCloseHistoryFile(&petsc_history)); 159302c9f0b5SLisandro Dalcin petsc_history = NULL; 1594e5c89e4eSSatish Balay } 15959566063dSJacob Faibussowitsch PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton)); 15969566063dSJacob Faibussowitsch PetscCall(PetscInfoDestroy()); 1597e5c89e4eSSatish Balay 159867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 159992f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1600e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 160192f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1602e5c89e4eSSatish Balay FILE *fd; 1603ed9cf6e9SBarry Smith int err; 1604e5c89e4eSSatish Balay 1605dc92acbaSJed Brown flg2 = PETSC_FALSE; 160692f119d6SBarry Smith flg3 = PETSC_FALSE; 16079566063dSJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL)); 16089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL)); 160992f119d6SBarry Smith fname[0] = 0; 16109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1)); 1611e5c89e4eSSatish Balay if (flg1 && fname[0]) { 16123ba16761SJacob Faibussowitsch PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank)); 16139371c9d4SSatish Balay fd = fopen(sname, "w"); 16149371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16159566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(fd)); 1616ed9cf6e9SBarry Smith err = fclose(fd); 161728b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 161892f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1619e5c89e4eSSatish Balay MPI_Comm local_comm; 1620e5c89e4eSSatish Balay 1621afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16229566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16239566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(stdout)); 16249566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1626e5c89e4eSSatish Balay } 1627e5c89e4eSSatish Balay fname[0] = 0; 16289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1)); 1629e5c89e4eSSatish Balay if (flg1 && fname[0]) { 16303ba16761SJacob Faibussowitsch PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank)); 16319371c9d4SSatish Balay fd = fopen(sname, "w"); 16329371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16339566063dSJacob Faibussowitsch PetscCall(PetscMallocView(fd)); 1634ed9cf6e9SBarry Smith err = fclose(fd); 163528b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 163692f119d6SBarry Smith } else if (flg1) { 163792f119d6SBarry Smith MPI_Comm local_comm; 163892f119d6SBarry Smith 1639afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16409566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16419566063dSJacob Faibussowitsch PetscCall(PetscMallocView(stdout)); 16429566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1644e5c89e4eSSatish Balay } 1645e5c89e4eSSatish Balay } 164667584ceeSBarry Smith #endif 164720e2c332SMatthew G. Knepley 16485486ca60SMatthew G. Knepley /* 16495486ca60SMatthew G. Knepley Close any open dynamic libraries 16505486ca60SMatthew G. Knepley */ 16519566063dSJacob Faibussowitsch PetscCall(PetscFinalize_DynamicLibraries()); 16525486ca60SMatthew G. Knepley 1653e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 16549566063dSJacob Faibussowitsch PetscCall(PetscOptionsDestroyDefault()); 1655e5c89e4eSSatish Balay 1656e5c89e4eSSatish Balay PetscGlobalArgc = 0; 165702c9f0b5SLisandro Dalcin PetscGlobalArgs = NULL; 1658e5c89e4eSSatish Balay 165971438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM) 166071438e86SJunchao Zhang if (PetscBeganNvshmem) { 16619566063dSJacob Faibussowitsch PetscCall(PetscNvshmemFinalize()); 166271438e86SJunchao Zhang PetscBeganNvshmem = PETSC_FALSE; 166371438e86SJunchao Zhang } 166471438e86SJunchao Zhang #endif 1665a0e72f99SJunchao Zhang 16669566063dSJacob Faibussowitsch PetscCall(PetscFreeMPIResources()); 1667e5c89e4eSSatish Balay 1668dbc8283eSBarry Smith /* 1669efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1670efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1671efb80d3cSBarry Smith 1672efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1673efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1674dbc8283eSBarry Smith */ 1675b770b1f6SSatish Balay { 1676dbc8283eSBarry Smith PetscCommCounter *counter; 1677dbc8283eSBarry Smith PetscMPIInt flg; 1678dbc8283eSBarry Smith MPI_Comm icomm; 16799371c9d4SSatish Balay union 16809371c9d4SSatish Balay { 16819371c9d4SSatish Balay MPI_Comm comm; 16829371c9d4SSatish Balay void *ptr; 16839371c9d4SSatish Balay } ucomm; 16849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 1685dbc8283eSBarry Smith if (flg) { 1686265f3f35SJed Brown icomm = ucomm.comm; 16879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 168828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1689dbc8283eSBarry Smith 16909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval)); 16919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 16929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1693dbc8283eSBarry Smith } 16949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 1695dbc8283eSBarry Smith if (flg) { 1696265f3f35SJed Brown icomm = ucomm.comm; 16979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 169828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1699dbc8283eSBarry Smith 17009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval)); 17019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 17029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1703dbc8283eSBarry Smith } 1704b770b1f6SSatish Balay } 1705dbc8283eSBarry Smith 17069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval)); 17079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval)); 17089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval)); 17099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval)); 171062e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval)); 171162e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval)); 1712480cf27aSJed Brown 17135ea2b939SDuncan Campbell // Free keyvals which may be silently created by some routines 17145ea2b939SDuncan Campbell if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval)); 17155ea2b939SDuncan Campbell if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval)); 17165ea2b939SDuncan Campbell 17179566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen)); 17189566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout)); 17199566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr)); 17209566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock)); 1721ef19f930SBarry Smith 1722e5c89e4eSSatish Balay if (PetscBeganMPI) { 172399b1327fSBarry Smith PetscMPIInt flag; 17249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalized(&flag)); 172539a651e2SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 172639a651e2SJacob Faibussowitsch /* wait until the very last moment to disable error handling */ 172739a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_FALSE; 17289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalize()); 172939a651e2SJacob Faibussowitsch } else PetscErrorHandlingInitialized = PETSC_FALSE; 173039a651e2SJacob Faibussowitsch 1731e5c89e4eSSatish Balay /* 1732e5c89e4eSSatish Balay 1733e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1734e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1735e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1736e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1737e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 17380e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1739e5c89e4eSSatish Balay memory was not freed. 1740e5c89e4eSSatish Balay 1741e5c89e4eSSatish Balay */ 17429566063dSJacob Faibussowitsch PetscCall(PetscMallocClear()); 17439566063dSJacob Faibussowitsch PetscCall(PetscStackReset()); 1744a297a907SKarl Rupp 1745e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1746e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 17477ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE) 174856883f60SBarry Smith /* 174956883f60SBarry Smith flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the 175056883f60SBarry Smith gcov files are still being added to the directories as git tries to remove the directories. 175156883f60SBarry Smith */ 175256883f60SBarry Smith __gcov_flush(); 175356883f60SBarry Smith #endif 17541724198aSStefano Zampini /* To match PetscFunctionBegin() at the beginning of this function */ 17551724198aSStefano Zampini PetscStackClearTop; 17563ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 1757e5c89e4eSSatish Balay } 1758e5c89e4eSSatish Balay 175943db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 1760d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b) 1761d71ae5a4SJacob Faibussowitsch { 176243db4dbbSBarry Smith if (*a == *b) return 1; 176343db4dbbSBarry Smith if (*a + 32 == *b) return 1; 176443db4dbbSBarry Smith if (*a - 32 == *b) return 1; 176543db4dbbSBarry Smith return 0; 176643db4dbbSBarry Smith } 1767a70650f6SBarry Smith #endif 176843db4dbbSBarry Smith 176943db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 1770d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b) 1771d71ae5a4SJacob Faibussowitsch { 177243db4dbbSBarry Smith if (*a == *b) return 1; 177343db4dbbSBarry Smith if (*a + 32 == *b) return 1; 177443db4dbbSBarry Smith if (*a - 32 == *b) return 1; 177543db4dbbSBarry Smith return 0; 177643db4dbbSBarry Smith } 177743db4dbbSBarry Smith #endif 1778