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) 156dd63270SBarry Smith #include <petsc/private/ftnimpl.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 87cc4c1da9SBarry Smith Collective, No Fortran Support 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 @*/ 129ce78bad3SBarry Smith PetscErrorCode PetscInitializeNoArguments(void) PeNS 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 1846497c311SBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 185d71ae5a4SJacob Faibussowitsch { 186e5c89e4eSSatish Balay PetscFunctionBegin; 1876497c311SBarry Smith if (*datatype == MPIU_INT_MPIINT && PetscDefined(USE_64BIT_INDICES)) { 1886497c311SBarry Smith #if defined(PETSC_USE_64BIT_INDICES) 1896497c311SBarry Smith struct petsc_mpiu_int_mpiint *xin = (struct petsc_mpiu_int_mpiint *)in, *xout = (struct petsc_mpiu_int_mpiint *)out; 1906497c311SBarry Smith PetscMPIInt count = *cnt; 191e5c89e4eSSatish Balay 1926497c311SBarry Smith for (PetscMPIInt i = 0; i < count; i++) { 1936497c311SBarry Smith xout[i].a = PetscMax(xout[i].a, xin[i].a); 1946497c311SBarry Smith xout[i].b += xin[i].b; 1956497c311SBarry Smith } 1966497c311SBarry Smith #endif 1976497c311SBarry Smith } else if (*datatype == MPIU_2INT || *datatype == MPIU_INT_MPIINT) { 1986497c311SBarry Smith PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out; 1996497c311SBarry Smith PetscMPIInt count = *cnt; 2006497c311SBarry Smith 2016497c311SBarry Smith for (PetscMPIInt i = 0; i < count; i++) { 202e5c89e4eSSatish Balay xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]); 203e5c89e4eSSatish Balay xout[2 * i + 1] += xin[2 * i + 1]; 204e5c89e4eSSatish Balay } 2056497c311SBarry Smith } else { 2066497c311SBarry Smith PetscErrorCode ierr = (*PetscErrorPrintf)("Can only handle MPIU_2INT and MPIU_INT_MPIINT data types"); 2076497c311SBarry Smith (void)ierr; 2086497c311SBarry Smith PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 2096497c311SBarry Smith } 210812af9f3SBarry Smith PetscFunctionReturnVoid(); 211e5c89e4eSSatish Balay } 212e5c89e4eSSatish Balay 2136497c311SBarry Smith /*@ 2146497c311SBarry Smith PetscMaxSum - Returns the max of the first entry over all MPI processes and the sum of the second entry. 215b693b147SBarry Smith 2166497c311SBarry Smith Collective 2176497c311SBarry Smith 2186497c311SBarry Smith Input Parameters: 2196497c311SBarry Smith + comm - the communicator 2206497c311SBarry Smith - array - an arry of length 2 times `size`, the number of MPI processes 2216497c311SBarry Smith 2226497c311SBarry Smith Output Parameters: 2236497c311SBarry Smith + max - the maximum of `array[2*rank]` over all MPI processes 2246497c311SBarry Smith - sum - the sum of the `array[2*rank + 1]` over all MPI processes 2256497c311SBarry Smith 2266497c311SBarry Smith Level: developer 2276497c311SBarry Smith 2286497c311SBarry Smith .seealso: `PetscInitialize()` 2296497c311SBarry Smith @*/ 2306497c311SBarry Smith PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt array[], PetscInt *max, PetscInt *sum) 231d71ae5a4SJacob Faibussowitsch { 232e5c89e4eSSatish Balay PetscFunctionBegin; 233d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 234d6e4c47cSJed Brown { 2359371c9d4SSatish Balay struct { 2369371c9d4SSatish Balay PetscInt max, sum; 2379371c9d4SSatish Balay } work; 2386497c311SBarry Smith PetscCallMPI(MPI_Reduce_scatter_block((void *)array, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm)); 239d6e4c47cSJed Brown *max = work.max; 240d6e4c47cSJed Brown *sum = work.sum; 241d6e4c47cSJed Brown } 242d6e4c47cSJed Brown #else 243d6e4c47cSJed Brown { 244d6e4c47cSJed Brown PetscMPIInt size, rank; 2459371c9d4SSatish Balay struct { 2469371c9d4SSatish Balay PetscInt max, sum; 2479371c9d4SSatish Balay } *work; 2489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &work)); 251462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce((void *)array, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm)); 2526ac3741eSJed Brown *max = work[rank].max; 2536ac3741eSJed Brown *sum = work[rank].sum; 2549566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 255d6e4c47cSJed Brown } 256d6e4c47cSJed Brown #endif 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 258e5c89e4eSSatish Balay } 259e5c89e4eSSatish Balay 260ef384fdeSPierre Jolivet #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16)) 261ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 262613bf2b2SPierre Jolivet #include <quadmath.h> 263613bf2b2SPierre Jolivet #endif 2649e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0; 265de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 26606a205a8SBarry Smith MPI_Op MPIU_SUM = 0; 267613bf2b2SPierre Jolivet #endif 268e5c89e4eSSatish Balay 269d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 270d71ae5a4SJacob Faibussowitsch { 2716497c311SBarry Smith PetscMPIInt i, count = *cnt; 272e5c89e4eSSatish Balay 273e5c89e4eSSatish Balay PetscFunctionBegin; 2747c2de775SJed Brown if (*datatype == MPIU_REAL) { 275e2e03761SBarry Smith PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 276a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] += xin[i]; 2777c2de775SJed Brown } 2787c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 279c5481aeeSJose E. Roman else if (*datatype == MPIU_COMPLEX) { 2807c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 281a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] += xin[i]; 2827c2de775SJed Brown } 2837c2de775SJed Brown #endif 284ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 285613bf2b2SPierre Jolivet else if (*datatype == MPIU___FLOAT128) { 286613bf2b2SPierre Jolivet __float128 *xin = (__float128 *)in, *xout = (__float128 *)out; 287613bf2b2SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 28874c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 2899371c9d4SSatish Balay } else if (*datatype == MPIU___COMPLEX128) { 290613bf2b2SPierre Jolivet __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out; 291613bf2b2SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 29274c01ffaSSatish Balay #endif 293613bf2b2SPierre Jolivet } 294613bf2b2SPierre Jolivet #endif 295ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16) 2969e517322SPierre Jolivet else if (*datatype == MPIU___FP16) { 2979e517322SPierre Jolivet __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out; 2986497c311SBarry Smith for (i = 0; i < count; i++) xout[i] = (__fp16)(xin[i] + xout[i]); 2999e517322SPierre Jolivet } 3009e517322SPierre Jolivet #endif 3017c2de775SJed Brown else { 302ef384fdeSPierre Jolivet #if (!defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_SKIP_REAL___FLOAT128)) && (!defined(PETSC_HAVE_REAL___FP16) || defined(PETSC_SKIP_REAL___FP16)) 303ef384fdeSPierre Jolivet PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types")); 304ef384fdeSPierre Jolivet #elif !defined(PETSC_HAVE_REAL___FP16) || defined(PETSC_SKIP_REAL___FP16) 3053ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types")); 306ef384fdeSPierre Jolivet #elif !defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_SKIP_REAL___FLOAT128) 3073ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types")); 3089e517322SPierre Jolivet #else 3093ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types")); 310613bf2b2SPierre Jolivet #endif 31141e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 312e2e03761SBarry Smith } 313812af9f3SBarry Smith PetscFunctionReturnVoid(); 314e5c89e4eSSatish Balay } 315e5c89e4eSSatish Balay #endif 316e5c89e4eSSatish Balay 317570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 318d9822059SBarry Smith MPI_Op MPIU_MAX = 0; 319d9822059SBarry Smith MPI_Op MPIU_MIN = 0; 320d9822059SBarry Smith 321d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 322d71ae5a4SJacob Faibussowitsch { 323d9822059SBarry Smith PetscInt i, count = *cnt; 324d9822059SBarry Smith 325d9822059SBarry Smith PetscFunctionBegin; 3267c2de775SJed Brown if (*datatype == MPIU_REAL) { 3278c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 328a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]); 3297c2de775SJed Brown } 3307c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 3317c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 3327c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 333ad540459SPierre Jolivet for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 3347c2de775SJed Brown } 3357c2de775SJed Brown #endif 3367c2de775SJed Brown else { 3373ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types")); 33841e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 3398c764dc5SJose Roman } 340d9822059SBarry Smith PetscFunctionReturnVoid(); 341d9822059SBarry Smith } 342d9822059SBarry Smith 343d71ae5a4SJacob Faibussowitsch PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 344d71ae5a4SJacob Faibussowitsch { 345d9822059SBarry Smith PetscInt i, count = *cnt; 346d9822059SBarry Smith 347d9822059SBarry Smith PetscFunctionBegin; 3487c2de775SJed Brown if (*datatype == MPIU_REAL) { 3498c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 350a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]); 3517c2de775SJed Brown } 3527c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 3537c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 3547c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 355ad540459SPierre Jolivet for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 3567c2de775SJed Brown } 3577c2de775SJed Brown #endif 3587c2de775SJed Brown else { 3593ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types")); 36041e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 3618c764dc5SJose Roman } 362d9822059SBarry Smith PetscFunctionReturnVoid(); 363d9822059SBarry Smith } 364d9822059SBarry Smith #endif 365d9822059SBarry Smith 366480cf27aSJed Brown /* 367480cf27aSJed Brown Private routine to delete internal tag/name counter storage when a communicator is freed. 368480cf27aSJed Brown 369ff0e51ddSBarry 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. 370480cf27aSJed Brown 37112801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 372480cf27aSJed Brown 373480cf27aSJed Brown */ 3748434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_DeleteFn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state) 375d71ae5a4SJacob Faibussowitsch { 37605342407SJunchao Zhang PetscCommCounter *counter = (PetscCommCounter *)count_val; 37757f21012SBarry Smith struct PetscCommStash *comms = counter->comms, *pcomm; 378480cf27aSJed Brown 379480cf27aSJed Brown PetscFunctionBegin; 3807c5b2466SBarry Smith PetscCallReturnMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm)); 3817c5b2466SBarry Smith PetscCallReturnMPI(PetscFree(counter->iflags)); 38257f21012SBarry Smith while (comms) { 3837c5b2466SBarry Smith PetscCallMPIReturnMPI(MPI_Comm_free(&comms->comm)); 38457f21012SBarry Smith pcomm = comms; 38557f21012SBarry Smith comms = comms->next; 3867c5b2466SBarry Smith PetscCallReturnMPI(PetscFree(pcomm)); 38757f21012SBarry Smith } 3887c5b2466SBarry Smith PetscCallReturnMPI(PetscFree(counter)); 389480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 390480cf27aSJed Brown } 391480cf27aSJed Brown 392480cf27aSJed Brown /* 39347435625SJed Brown This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user 394da3039f7SJed Brown calls MPI_Comm_free(). 395da3039f7SJed Brown 396da3039f7SJed Brown This is the only entry point for breaking the links between inner and outer comms. 397480cf27aSJed Brown 398ff0e51ddSBarry Smith This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator. 399480cf27aSJed Brown 40012801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 401480cf27aSJed Brown 402480cf27aSJed Brown */ 4038434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_DeleteFn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) 404d71ae5a4SJacob Faibussowitsch { 4059371c9d4SSatish Balay union 406480cf27aSJed Brown { 4079371c9d4SSatish Balay MPI_Comm comm; 4089371c9d4SSatish Balay void *ptr; 4099371c9d4SSatish Balay } icomm; 410480cf27aSJed Brown 411480cf27aSJed Brown PetscFunctionBegin; 4127c5b2466SBarry Smith PetscCheckReturnMPI(keyval == Petsc_InnerComm_keyval, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval"); 413ec4fadc2SJed Brown icomm.ptr = attr_val; 41476bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 41533779a13SJunchao Zhang /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */ 41633779a13SJunchao Zhang PetscMPIInt flg; 4179371c9d4SSatish Balay union 4189371c9d4SSatish Balay { 4199371c9d4SSatish Balay MPI_Comm comm; 4209371c9d4SSatish Balay void *ptr; 4219371c9d4SSatish Balay } ocomm; 4227c5b2466SBarry Smith PetscCallMPIReturnMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg)); 4237c5b2466SBarry Smith PetscCheckReturnMPI(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute"); 4247c5b2466SBarry Smith PetscCheckReturnMPI(ocomm.comm == comm, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm's OuterComm attribute does not point to outer PETSc comm"); 42533779a13SJunchao Zhang } 4267c5b2466SBarry Smith PetscCallMPIReturnMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval)); 4277c5b2466SBarry Smith PetscCallReturnMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm)); 428da3039f7SJed Brown PetscFunctionReturn(MPI_SUCCESS); 429b89831e5SBarry Smith } 430da3039f7SJed Brown 431da3039f7SJed Brown /* 4328434afd1SBarry 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. 433da3039f7SJed Brown */ 4348434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_DeleteFn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) 435d71ae5a4SJacob Faibussowitsch { 436da3039f7SJed Brown PetscFunctionBegin; 4377c5b2466SBarry Smith PetscCallReturnMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm)); 438480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 439480cf27aSJed Brown } 440480cf27aSJed Brown 4418434afd1SBarry Smith PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_DeleteFn(MPI_Comm, PetscMPIInt, void *, void *); 44242218b76SBarry Smith 443951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 4448cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *); 4458cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 4468cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 447e39fd77fSBarry Smith #endif 448e39fd77fSBarry Smith 4490f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE; 45012801b39SBarry Smith 451eb27c7c8SSatish Balay PETSC_INTERN int PetscGlobalArgc; 4529f0612e4SBarry Smith PETSC_INTERN char **PetscGlobalArgs, **PetscGlobalArgsFortran; 4536ae9a8a6SBarry Smith int PetscGlobalArgc = 0; 45402c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL; 4559f0612e4SBarry Smith char **PetscGlobalArgsFortran = NULL; 456dff31646SBarry Smith PetscSegBuffer PetscCitationsList; 457e5c89e4eSSatish Balay 458d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCitationsInitialize(void) 459d71ae5a4SJacob Faibussowitsch { 460051e4cf2SJed Brown PetscFunctionBegin; 4619566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList)); 46230c35bf2SSatish Balay 46330c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\ 46430c35bf2SSatish Balay Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\ 46530c35bf2SSatish Balay and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\ 4663f81df79SBarry Smith and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\ 46730c35bf2SSatish Balay and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\ 46830c35bf2SSatish Balay and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\ 4699281ddf3SSatish Balay and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith and Hansol Suh\n\ 47030c35bf2SSatish Balay and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\ 47130c35bf2SSatish Balay Title = {{PETSc/TAO} Users Manual},\n\ 47294492ad7SSatish Balay Number = {ANL-21/39 - Revision 3.23},\n\ 473fef8b845SSatish Balay Doi = {10.2172/2565610},\n\ 47430c35bf2SSatish Balay Institution = {Argonne National Laboratory},\n\ 47594492ad7SSatish Balay Year = {2025}\n}\n", 4769371c9d4SSatish Balay NULL)); 47730c35bf2SSatish Balay 47830c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\ 47930c35bf2SSatish Balay Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\ 48030c35bf2SSatish Balay Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\ 48130c35bf2SSatish Balay Booktitle = {Modern Software Tools in Scientific Computing},\n\ 48230c35bf2SSatish Balay Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\ 48330c35bf2SSatish Balay Pages = {163--202},\n\ 48430c35bf2SSatish Balay Publisher = {Birkh{\\\"{a}}user Press},\n\ 4859371c9d4SSatish Balay Year = {1997}\n}\n", 4869371c9d4SSatish Balay NULL)); 4873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 488051e4cf2SJed Brown } 489e5c89e4eSSatish Balay 4902d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */ 4912d747510SLisandro Dalcin 492d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSetProgramName(const char name[]) 493d71ae5a4SJacob Faibussowitsch { 4942d747510SLisandro Dalcin PetscFunctionBegin; 4959566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(programname, name, sizeof(programname))); 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4972d747510SLisandro Dalcin } 4982d747510SLisandro Dalcin 4992d747510SLisandro Dalcin /*@C 5002d747510SLisandro Dalcin PetscGetProgramName - Gets the name of the running program. 5012d747510SLisandro Dalcin 5022d747510SLisandro Dalcin Not Collective 5032d747510SLisandro Dalcin 5042d747510SLisandro Dalcin Input Parameter: 5052d747510SLisandro Dalcin . len - length of the string name 5062d747510SLisandro Dalcin 5072d747510SLisandro Dalcin Output Parameter: 508811af0c4SBarry Smith . name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN` 5092d747510SLisandro Dalcin 5102d747510SLisandro Dalcin Level: advanced 5112d747510SLisandro Dalcin 51221532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()` 5132d747510SLisandro Dalcin @*/ 514d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProgramName(char name[], size_t len) 515d71ae5a4SJacob Faibussowitsch { 5162d747510SLisandro Dalcin PetscFunctionBegin; 5179566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(name, programname, len)); 5183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5192d747510SLisandro Dalcin } 5202d747510SLisandro Dalcin 521e5c89e4eSSatish Balay /*@C 522e5c89e4eSSatish Balay PetscGetArgs - Allows you to access the raw command line arguments anywhere 5230b4b7b1cSBarry Smith after `PetscInitialize()` is called but before `PetscFinalize()`. 524e5c89e4eSSatish Balay 525cc4c1da9SBarry Smith Not Collective, No Fortran Support 526e5c89e4eSSatish Balay 527e5c89e4eSSatish Balay Output Parameters: 5280b4b7b1cSBarry Smith + argc - count of the number of command line arguments 529e5c89e4eSSatish Balay - args - the command line arguments 530e5c89e4eSSatish Balay 531e5c89e4eSSatish Balay Level: intermediate 532e5c89e4eSSatish Balay 533e5c89e4eSSatish Balay Notes: 534e5c89e4eSSatish Balay This is usually used to pass the command line arguments into other libraries 535e5c89e4eSSatish Balay that are called internally deep in PETSc or the application. 536e5c89e4eSSatish Balay 5379f0612e4SBarry Smith The first argument contains the program name as is normal for C programs. 538f177e3b1SBarry Smith 5390b4b7b1cSBarry Smith See `PetscGetArguments()` for a variant of this routine. 5400b4b7b1cSBarry Smith 54121532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()` 542e5c89e4eSSatish Balay @*/ 543d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args) 544d71ae5a4SJacob Faibussowitsch { 545e5c89e4eSSatish Balay PetscFunctionBegin; 54639a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 547e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 548e5c89e4eSSatish Balay *args = PetscGlobalArgs; 5493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 550e5c89e4eSSatish Balay } 551e5c89e4eSSatish Balay 552793721a6SBarry Smith /*@C 553793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 554811af0c4SBarry Smith after `PetscInitialize()` is called but before `PetscFinalize()`. 555793721a6SBarry Smith 556cc4c1da9SBarry Smith Not Collective, No Fortran Support 557793721a6SBarry Smith 5582fe279fdSBarry Smith Output Parameter: 559793721a6SBarry Smith . args - the command line arguments 560793721a6SBarry Smith 561793721a6SBarry Smith Level: intermediate 562793721a6SBarry Smith 56321532e8aSBarry Smith Note: 5640b4b7b1cSBarry Smith This does NOT start with the program name and IS `NULL` terminated (the final argument is void) 5650b4b7b1cSBarry Smith 5660b4b7b1cSBarry Smith Use `PetscFreeArguments()` to return the memory used by the arguments. 5670b4b7b1cSBarry Smith 5680b4b7b1cSBarry Smith This makes a copy of the arguments and the array of arguments, while `PetscGetArgs()` does not make a copy, 5690b4b7b1cSBarry Smith it returns the array of arguments that was passed into the main program. 570793721a6SBarry Smith 57121532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`, `PetscInitialize()` 572793721a6SBarry Smith @*/ 573d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args) 574d71ae5a4SJacob Faibussowitsch { 575793721a6SBarry Smith PetscInt i, argc = PetscGlobalArgc; 576793721a6SBarry Smith 577793721a6SBarry Smith PetscFunctionBegin; 57839a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 5799371c9d4SSatish Balay if (!argc) { 5809371c9d4SSatish Balay *args = NULL; 5813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5829371c9d4SSatish Balay } 5839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(argc, args)); 5849566063dSJacob Faibussowitsch for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i])); 5852d747510SLisandro Dalcin (*args)[argc - 1] = NULL; 5863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 587793721a6SBarry Smith } 588793721a6SBarry Smith 589793721a6SBarry Smith /*@C 590811af0c4SBarry Smith PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()` 591793721a6SBarry Smith 592cc4c1da9SBarry Smith Not Collective, No Fortran Support 593793721a6SBarry Smith 5942fe279fdSBarry Smith Output Parameter: 595793721a6SBarry Smith . args - the command line arguments 596793721a6SBarry Smith 597793721a6SBarry Smith Level: intermediate 598793721a6SBarry Smith 5990b4b7b1cSBarry Smith Developer Note: 6000b4b7b1cSBarry Smith This should be PetscRestoreArguments() 6010b4b7b1cSBarry Smith 602db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()` 603793721a6SBarry Smith @*/ 604d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args) 605d71ae5a4SJacob Faibussowitsch { 60639a651e2SJacob Faibussowitsch PetscFunctionBegin; 60739a651e2SJacob Faibussowitsch if (args) { 608793721a6SBarry Smith PetscInt i = 0; 609793721a6SBarry Smith 6109566063dSJacob Faibussowitsch while (args[i]) PetscCall(PetscFree(args[i++])); 6119566063dSJacob Faibussowitsch PetscCall(PetscFree(args)); 61239a651e2SJacob Faibussowitsch } 6133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 614793721a6SBarry Smith } 615793721a6SBarry Smith 61627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 61730befbd2SBarry Smith #include <petscconfiginfo.h> 61830befbd2SBarry Smith 619d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) 620d71ae5a4SJacob Faibussowitsch { 62127104ee2SJacob Faibussowitsch PetscFunctionBegin; 62211525c0dSBarry Smith if (!PetscGlobalRank) { 62330befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64]; 62411525c0dSBarry Smith int port; 625ffbd1cfbSBarry Smith PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE; 62611525c0dSBarry Smith size_t applinelen, introlen; 627ffbd1cfbSBarry Smith char sawsurl[256]; 62811525c0dSBarry Smith 6299566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg)); 63011525c0dSBarry Smith if (flg) { 63111525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 63211525c0dSBarry Smith 6339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL)); 63411525c0dSBarry Smith if (sawslog[0]) { 635792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog)); 63611525c0dSBarry Smith } else { 637792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL)); 63811525c0dSBarry Smith } 63911525c0dSBarry Smith } 6409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg)); 64148a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert)); 6429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL)); 643ffbd1cfbSBarry Smith if (selectport) { 644792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 645792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 646ffbd1cfbSBarry Smith } else { 6479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg)); 64848a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Port, (port)); 649ffbd1cfbSBarry Smith } 6509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg)); 65111525c0dSBarry Smith if (flg) { 652792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 6539566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(root, ".", &rootlocal)); 6549c1e0ce8SBarry Smith } else { 6559566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg)); 6569c1e0ce8SBarry Smith if (flg) { 6579566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root))); 658792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 6599c1e0ce8SBarry Smith } 66011525c0dSBarry Smith } 6619566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2)); 66211525c0dSBarry Smith if (flg2) { 66311525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 66428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option"); 6659566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root)); 6669566063dSJacob Faibussowitsch PetscCall(PetscTestDirectory(jsdir, 'r', &flg)); 66728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory"); 668792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Local_Header, ()); 66911525c0dSBarry Smith } 6709566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(programname, sizeof(programname))); 6719566063dSJacob Faibussowitsch PetscCall(PetscStrlen(help, &applinelen)); 67211525c0dSBarry Smith introlen = 4096 + applinelen; 67330a8c9c0SSurtai Han applinelen += 1024; 6749566063dSJacob Faibussowitsch PetscCall(PetscMalloc(applinelen, &appline)); 6759566063dSJacob Faibussowitsch PetscCall(PetscMalloc(introlen, &intro)); 67611525c0dSBarry Smith 67711525c0dSBarry Smith if (rootlocal) { 6789566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname)); 6799566063dSJacob Faibussowitsch PetscCall(PetscTestFile(appline, 'r', &rootlocal)); 68011525c0dSBarry Smith } 6819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetAll(NULL, &options)); 68211525c0dSBarry Smith if (rootlocal && help) { 6839566063dSJacob 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)); 68411525c0dSBarry Smith } else if (help) { 6859566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help)); 68611525c0dSBarry Smith } else { 6879566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options)); 68811525c0dSBarry Smith } 6899566063dSJacob Faibussowitsch PetscCall(PetscFree(options)); 6909566063dSJacob Faibussowitsch PetscCall(PetscGetVersion(version, sizeof(version))); 6919371c9d4SSatish Balay PetscCall(PetscSNPrintf(intro, introlen, 6929371c9d4SSatish Balay "<body>\n" 693a17b96a8SKyle 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" 694df62c222SBarry 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" 6959371c9d4SSatish Balay "%s", 6969371c9d4SSatish Balay version, petscconfigureoptions, appline)); 697792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro)); 6989566063dSJacob Faibussowitsch PetscCall(PetscFree(intro)); 6999566063dSJacob Faibussowitsch PetscCall(PetscFree(appline)); 700ffbd1cfbSBarry Smith if (selectport) { 701aa573868SBarry Smith PetscBool silent; 7027d812c46SBarry Smith 7037d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 70439a651e2SJacob Faibussowitsch while (SAWs_Initialize()) { 705792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 706792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 7077d812c46SBarry Smith } 7087d812c46SBarry Smith 7099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL)); 710aa573868SBarry Smith if (!silent) { 711792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl)); 7129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl)); 713ffbd1cfbSBarry Smith } 7147d812c46SBarry Smith } else { 715792fecdfSBarry Smith PetscCallSAWs(SAWs_Initialize, ()); 716aa573868SBarry Smith } 7179566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister("@TechReport{ saws,\n" 7180af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 7190af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 7200af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 7219371c9d4SSatish Balay " Year = 2013\n}\n", 7229371c9d4SSatish Balay NULL)); 72311525c0dSBarry Smith } 7243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72511525c0dSBarry Smith } 72611525c0dSBarry Smith #endif 72711525c0dSBarry Smith 7284dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 729d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) 730d71ae5a4SJacob Faibussowitsch { 7314dfee713SSatish Balay PetscFunctionBegin; 7324dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 7334dfee713SSatish Balay /* see MPI.py for details on this bug */ 7344dfee713SSatish Balay (void)setenv("HWLOC_COMPONENTS", "-x86", 1); 7354dfee713SSatish Balay #endif 7363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7374dfee713SSatish Balay } 7384dfee713SSatish Balay 739a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS) 740a56f64adSBarry Smith #include <adios.h> 74122580e64SBarry Smith #include <adios_read.h> 7427b56e58cSSatish Balay int64_t Petsc_adios_group; 743a56f64adSBarry Smith #endif 744a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP) 745cd1b99a6SBarry Smith #include <omp.h> 746f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 747cd1b99a6SBarry Smith #endif 748a56f64adSBarry Smith 749a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h> 750a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA) 7510e6b6b59SJacob Faibussowitsch #include <petscdevice_cuda.h> 752a4af0ceeSJacob Faibussowitsch // REMOVE ME 753a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL; 754a4af0ceeSJacob Faibussowitsch #endif 755a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP) 7560e6b6b59SJacob Faibussowitsch #include <petscdevice_hip.h> 757a4af0ceeSJacob Faibussowitsch // REMOVE ME 758a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL; 759a4af0ceeSJacob Faibussowitsch #endif 760a4af0ceeSJacob Faibussowitsch 76127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H) 762bc8a24ecSBarry Smith #include <dlfcn.h> 763bc8a24ecSBarry Smith #endif 764a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void); 765a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 7663274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscViennaCLInit(void); 767a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE; 768a4af0ceeSJacob Faibussowitsch #endif 769bc8a24ecSBarry Smith 770660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE; 771660278c0SBarry Smith 77285649d77SJunchao Zhang /* 77385649d77SJunchao Zhang PetscInitialize_Common - shared code between C and Fortran initialization 77485649d77SJunchao Zhang 77585649d77SJunchao Zhang prog: program name 77602101c96SSatish Balay file: optional PETSc database file name. Might be in Fortran string format when 'ftn' is true 77785649d77SJunchao Zhang help: program help message 778da81f932SPierre Jolivet ftn: is it called from Fortran initialization (petscinitializef_)? 779daa8fb3bSBarry Smith len: length of file string, used when Fortran is true 78085649d77SJunchao Zhang */ 781daa8fb3bSBarry Smith PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscInt len) 782d71ae5a4SJacob Faibussowitsch { 78385649d77SJunchao Zhang PetscMPIInt size; 78485649d77SJunchao Zhang PetscBool flg = PETSC_TRUE; 78585649d77SJunchao Zhang char hostname[256]; 786046ed546SBarry Smith PetscBool blas_view_flag = PETSC_FALSE; 78785649d77SJunchao Zhang 78827104ee2SJacob Faibussowitsch PetscFunctionBegin; 7893ba16761SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS); 79039a651e2SJacob Faibussowitsch /* these must be initialized in a routine, not as a constant declaration */ 79139a651e2SJacob Faibussowitsch PETSC_STDOUT = stdout; 79239a651e2SJacob Faibussowitsch PETSC_STDERR = stderr; 79339a651e2SJacob Faibussowitsch 7949566063dSJacob Faibussowitsch /* PetscCall can be used from now */ 79539a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_TRUE; 79639a651e2SJacob Faibussowitsch 79785649d77SJunchao Zhang /* 79885649d77SJunchao Zhang The checking over compatible runtime libraries is complicated by the MPI ABI initiative 79985649d77SJunchao Zhang https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 80085649d77SJunchao Zhang MPICH v3.1 (Released February 2014) 80185649d77SJunchao Zhang IBM MPI v2.1 (December 2014) 80285649d77SJunchao Zhang Intel MPI Library v5.0 (2014) 80385649d77SJunchao Zhang Cray MPT v7.0.0 (June 2014) 80485649d77SJunchao Zhang As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 80585649d77SJunchao Zhang listed above and since that time are compatible. 80685649d77SJunchao Zhang 80785649d77SJunchao Zhang Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 80885649d77SJunchao Zhang at compile time or runtime. Thus we will need to systematically track the allowed versions 80985649d77SJunchao Zhang and how they are represented in the mpi.h and MPI_Get_library_version() output in order 81085649d77SJunchao Zhang to perform the checking. 81185649d77SJunchao Zhang 81285649d77SJunchao Zhang Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 81385649d77SJunchao Zhang 81485649d77SJunchao Zhang Questions: 81585649d77SJunchao Zhang 81685649d77SJunchao Zhang Should the checks for ABI incompatibility be only on the major version number below? 81785649d77SJunchao Zhang Presumably the output to stderr will be removed before a release. 81885649d77SJunchao Zhang */ 81985649d77SJunchao Zhang 82085649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 82185649d77SJunchao Zhang { 82285649d77SJunchao Zhang char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 82385649d77SJunchao Zhang PetscMPIInt mpilibraryversionlength; 82439a651e2SJacob Faibussowitsch 8259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength)); 82685649d77SJunchao Zhang /* check for MPICH versions before MPI ABI initiative */ 82785649d77SJunchao Zhang #if defined(MPICH_VERSION) 82885649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000 82985649d77SJunchao Zhang { 83085649d77SJunchao Zhang char *ver, *lf; 83185649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 83239a651e2SJacob Faibussowitsch 8339566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver)); 83439a651e2SJacob Faibussowitsch if (ver) { 8359566063dSJacob Faibussowitsch PetscCall(PetscStrchr(ver, '\n', &lf)); 83639a651e2SJacob Faibussowitsch if (lf) { 83785649d77SJunchao Zhang *lf = 0; 8389566063dSJacob Faibussowitsch PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg)); 83985649d77SJunchao Zhang } 84085649d77SJunchao Zhang } 84185649d77SJunchao Zhang if (!flg) { 842d5b396e8SSatish Balay PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION)); 84385649d77SJunchao Zhang flg = PETSC_TRUE; 84485649d77SJunchao Zhang } 84585649d77SJunchao Zhang } 84685649d77SJunchao Zhang #endif 84785649d77SJunchao 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?) */ 848100ffedbSJunchao Zhang #elif defined(PETSC_HAVE_OPENMPI) 84985649d77SJunchao Zhang { 85085649d77SJunchao Zhang char *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf; 85185649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 85285649d77SJunchao Zhang #define PSTRSZ 2 85385649d77SJunchao Zhang char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"}; 85485649d77SJunchao Zhang char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "}; 85585649d77SJunchao Zhang int i; 85685649d77SJunchao Zhang for (i = 0; i < PSTRSZ; i++) { 8579566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver)); 85839a651e2SJacob Faibussowitsch if (ver) { 859100ffedbSJunchao Zhang PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], PETSC_PKG_OPENMPI_VERSION_MAJOR, PETSC_PKG_OPENMPI_VERSION_MINOR)); 8609566063dSJacob Faibussowitsch PetscCall(PetscStrstr(ver, bs, &bsf)); 86139a651e2SJacob Faibussowitsch if (bsf) flg = PETSC_TRUE; 86285649d77SJunchao Zhang break; 86385649d77SJunchao Zhang } 86485649d77SJunchao Zhang } 86585649d77SJunchao Zhang if (!flg) { 866100ffedbSJunchao Zhang PetscCall(PetscInfo(NULL, "PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n", mpilibraryversion, PETSC_PKG_OPENMPI_VERSION_MAJOR, PETSC_PKG_OPENMPI_VERSION_MINOR)); 86785649d77SJunchao Zhang flg = PETSC_TRUE; 86885649d77SJunchao Zhang } 86985649d77SJunchao Zhang } 87085649d77SJunchao Zhang #endif 87185649d77SJunchao Zhang } 87285649d77SJunchao Zhang #endif 87385649d77SJunchao Zhang 874e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__)) 87585649d77SJunchao 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 */ 87639a651e2SJacob 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"); 87785649d77SJunchao Zhang #endif 87885649d77SJunchao Zhang 8799566063dSJacob Faibussowitsch PetscCall(PetscOptionsCreateDefault()); 88085649d77SJunchao Zhang 88185649d77SJunchao Zhang PetscFinalizeCalled = PETSC_FALSE; 88285649d77SJunchao Zhang 8839566063dSJacob Faibussowitsch PetscCall(PetscSetProgramName(prog)); 8849566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen)); 8859566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout)); 8869566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr)); 8879566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscCommSpinLock)); 88885649d77SJunchao Zhang 88985649d77SJunchao Zhang if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 8909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN)); 89185649d77SJunchao Zhang 89285649d77SJunchao Zhang if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) { 8939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS)); 8949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE)); 89585649d77SJunchao Zhang } 89685649d77SJunchao Zhang 89785649d77SJunchao Zhang /* Done after init due to a bug in MPICH-GM? */ 8989566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 89985649d77SJunchao Zhang 9009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank)); 9019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize)); 90285649d77SJunchao Zhang 9031dc74096SMartin Diehl MPIU_BOOL = MPI_C_BOOL; 90485649d77SJunchao Zhang MPIU_ENUM = MPI_INT; 90585649d77SJunchao Zhang MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64; 90685649d77SJunchao Zhang if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 90785649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 90885649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG) 90985649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 91085649d77SJunchao Zhang #endif 91139a651e2SJacob Faibussowitsch else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t"); 91285649d77SJunchao Zhang 91385649d77SJunchao Zhang /* 91485649d77SJunchao Zhang Initialized the global complex variable; this is because with 91585649d77SJunchao Zhang shared libraries the constructors for global variables 91685649d77SJunchao Zhang are not called; at least on IRIX. 91785649d77SJunchao Zhang */ 91885649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX) 91985649d77SJunchao Zhang { 92085649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 92185649d77SJunchao Zhang PetscComplex ic(0.0, 1.0); 92285649d77SJunchao Zhang PETSC_i = ic; 92385649d77SJunchao Zhang #else 92485649d77SJunchao Zhang PETSC_i = _Complex_I; 92585649d77SJunchao Zhang #endif 92685649d77SJunchao Zhang } 92785649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */ 92885649d77SJunchao Zhang 92985649d77SJunchao Zhang /* 93085649d77SJunchao Zhang Create the PETSc MPI reduction operator that sums of the first 93185649d77SJunchao Zhang half of the entries and maxes the second half. 93285649d77SJunchao Zhang */ 9339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP)); 93485649d77SJunchao Zhang 935ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 9369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128)); 9379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128)); 93874c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 9399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128)); 9409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128)); 94185649d77SJunchao Zhang #endif 94274c01ffaSSatish Balay #endif 943ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16) 9449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16)); 9459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FP16)); 94685649d77SJunchao Zhang #endif 94785649d77SJunchao Zhang 94885649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 9499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM)); 950613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX)); 951613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN)); 952ef384fdeSPierre Jolivet #elif (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16)) 9539e517322SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128)); 95485649d77SJunchao Zhang #endif 95585649d77SJunchao Zhang 9569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR)); 95762e5d2d2SJDBetteridge PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp)); 9589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR)); 95985649d77SJunchao Zhang 96085649d77SJunchao Zhang /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */ 96185649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI) 96285649d77SJunchao Zhang { 96385649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 96493d501b3SJacob Faibussowitsch MPI_Aint blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)}; 96585649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_REAL, MPIU_INT}, tmpStruct; 96685649d77SJunchao Zhang 9679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 96893d501b3SJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT)); 9699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT)); 97185649d77SJunchao Zhang } 97285649d77SJunchao Zhang { 97385649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 97493d501b3SJacob Faibussowitsch MPI_Aint blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)}; 97585649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_SCALAR, MPIU_INT}, tmpStruct; 97685649d77SJunchao Zhang 9779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 97893d501b3SJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT)); 9799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT)); 98185649d77SJunchao Zhang } 98285649d77SJunchao Zhang #endif 98385649d77SJunchao Zhang 98485649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES) 9859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT)); 9869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2INT)); 9876497c311SBarry Smith 9886497c311SBarry Smith #if !defined(PETSC_HAVE_MPIUNI) 9896497c311SBarry Smith { 9906497c311SBarry Smith int blockSizes[] = {1, 1}; 9916497c311SBarry Smith MPI_Aint blockOffsets[] = {offsetof(struct petsc_mpiu_int_mpiint, a), offsetof(struct petsc_mpiu_int_mpiint, b)}; 9926497c311SBarry Smith MPI_Datatype blockTypes[] = {MPIU_INT, MPI_INT}, tmpStruct; 9936497c311SBarry Smith 9946497c311SBarry Smith PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 9956497c311SBarry Smith PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_int_mpiint), &MPIU_INT_MPIINT)); 9966497c311SBarry Smith PetscCallMPI(MPI_Type_free(&tmpStruct)); 9976497c311SBarry Smith PetscCallMPI(MPI_Type_commit(&MPIU_INT_MPIINT)); 9986497c311SBarry Smith } 9996497c311SBarry Smith #endif 100085649d77SJunchao Zhang #endif 10019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT)); 10029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPI_4INT)); 10039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT)); 10049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_4INT)); 100585649d77SJunchao Zhang 100685649d77SJunchao Zhang /* 100785649d77SJunchao Zhang Attributes to be set on PETSc communicators 100885649d77SJunchao Zhang */ 1009c8025a54SPierre Jolivet PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_DeleteFn, &Petsc_Counter_keyval, NULL)); 1010c8025a54SPierre Jolivet PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_DeleteFn, &Petsc_InnerComm_keyval, NULL)); 1011c8025a54SPierre Jolivet PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_DeleteFn, &Petsc_OuterComm_keyval, NULL)); 1012c8025a54SPierre Jolivet PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_DeleteFn, &Petsc_ShmComm_keyval, NULL)); 1013c8025a54SPierre Jolivet PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, NULL)); 1014c8025a54SPierre Jolivet PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, NULL)); 101585649d77SJunchao Zhang 1016fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS) 1017daa8fb3bSBarry Smith if (ftn) PetscCall(PetscInitFortran_Private(file, len)); 101885649d77SJunchao Zhang else 101985649d77SJunchao Zhang #endif 10209566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file)); 102185649d77SJunchao Zhang 1022983cdb30SBarry Smith if (PetscDefined(HAVE_MPIUNI)) { 1023983cdb30SBarry Smith const char *mpienv = getenv("PMI_SIZE"); 1024983cdb30SBarry Smith if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE"); 1025983cdb30SBarry Smith if (mpienv) { 1026983cdb30SBarry Smith PetscInt isize; 1027983cdb30SBarry Smith PetscBool mflag = PETSC_FALSE; 1028983cdb30SBarry Smith 1029983cdb30SBarry Smith PetscCall(PetscOptionsStringToInt(mpienv, &isize)); 1030983cdb30SBarry Smith PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpiuni-allow-multiprocess-launch", &mflag, NULL)); 1031983cdb30SBarry Smith PetscCheck(isize == 1 || mflag, 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. Or run with -mpiuni-allow-multiprocess-launch to allow multiple independent MPI-uni jobs."); 1032983cdb30SBarry Smith } 1033983cdb30SBarry Smith } 1034983cdb30SBarry Smith 103585649d77SJunchao Zhang /* call a second time so it can look in the options database */ 10369566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 103785649d77SJunchao Zhang 103885649d77SJunchao Zhang /* 103985649d77SJunchao Zhang Check system options and print help 104085649d77SJunchao Zhang */ 10419566063dSJacob Faibussowitsch PetscCall(PetscOptionsCheckInitial_Private(help)); 104285649d77SJunchao Zhang 1043a4af0ceeSJacob Faibussowitsch /* 10440e6b6b59SJacob Faibussowitsch Creates the logging data structures; this is enabled even if logging is not turned on 10450e6b6b59SJacob Faibussowitsch This is the last thing we do before returning to the user code to prevent having the 10460e6b6b59SJacob Faibussowitsch logging numbers contaminated by any startup time associated with MPI 10470e6b6b59SJacob Faibussowitsch */ 10480e6b6b59SJacob Faibussowitsch PetscCall(PetscLogInitialize()); 10490e6b6b59SJacob Faibussowitsch 10500e6b6b59SJacob Faibussowitsch /* 1051a4af0ceeSJacob Faibussowitsch Initialize PetscDevice and PetscDeviceContext 1052a4af0ceeSJacob Faibussowitsch 1053a4af0ceeSJacob Faibussowitsch Note to any future devs thinking of moving this, proper initialization requires: 1054a4af0ceeSJacob Faibussowitsch 1. MPI initialized 1055a4af0ceeSJacob Faibussowitsch 2. Options DB initialized 1056f0b74427SPierre Jolivet 3. PETSc error handling initialized, specifically signal handlers. This expects to set up 10570e6b6b59SJacob Faibussowitsch its own SIGSEV handler via the push/pop interface. 10580e6b6b59SJacob Faibussowitsch 4. Logging initialized 1059a4af0ceeSJacob Faibussowitsch */ 10609566063dSJacob Faibussowitsch PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD)); 1061a4af0ceeSJacob Faibussowitsch 1062a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 1063a4af0ceeSJacob Faibussowitsch flg = PETSC_FALSE; 1064d75802c7SJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg)); 10659566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL)); 1066a4af0ceeSJacob Faibussowitsch PetscViennaCLSynchronize = flg; 10679566063dSJacob Faibussowitsch PetscCall(PetscViennaCLInit()); 1068a4af0ceeSJacob Faibussowitsch #endif 1069a4af0ceeSJacob Faibussowitsch 10709566063dSJacob Faibussowitsch PetscCall(PetscCitationsInitialize()); 107185649d77SJunchao Zhang 107285649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS) 10739566063dSJacob Faibussowitsch PetscCall(PetscInitializeSAWs(ftn ? NULL : help)); 107427104ee2SJacob Faibussowitsch flg = PETSC_FALSE; 10759566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg)); 10769566063dSJacob Faibussowitsch if (flg) PetscCall(PetscStackViewSAWs()); 107785649d77SJunchao Zhang #endif 107885649d77SJunchao Zhang 107985649d77SJunchao Zhang /* 108085649d77SJunchao Zhang Load the dynamic libraries (on machines that support them), this registers all 108185649d77SJunchao Zhang the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 108285649d77SJunchao Zhang */ 10839566063dSJacob Faibussowitsch PetscCall(PetscInitialize_DynamicLibraries()); 108485649d77SJunchao Zhang 10859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 10869566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size)); 108796dcfd6cSBarry Smith PetscCall(PetscGetHostName(hostname, sizeof(hostname))); 10889566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname)); 108985649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP) 109085649d77SJunchao Zhang { 109185649d77SJunchao Zhang PetscBool omp_view_flag; 109285649d77SJunchao Zhang char *threads = getenv("OMP_NUM_THREADS"); 109385649d77SJunchao Zhang 109485649d77SJunchao Zhang if (threads) { 10959566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads)); 109685649d77SJunchao Zhang (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads); 109785649d77SJunchao Zhang } else { 10982f840973SStefano Zampini PetscNumOMPThreads = (PetscInt)omp_get_max_threads(); 10999566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads)); 110085649d77SJunchao Zhang } 1101d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys"); 11029566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg)); 11039566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag)); 1104d0609cedSBarry Smith PetscOptionsEnd(); 110585649d77SJunchao Zhang if (flg) { 11061a0ee928SPierre Jolivet PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads)); 110785649d77SJunchao Zhang omp_set_num_threads((int)PetscNumOMPThreads); 110885649d77SJunchao Zhang } 110948a46eb9SPierre Jolivet if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads)); 111085649d77SJunchao Zhang } 111185649d77SJunchao Zhang #endif 111285649d77SJunchao Zhang 1113046ed546SBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "BLAS options", "Sys"); 1114046ed546SBarry Smith PetscCall(PetscOptionsName("-blas_view", "Display number of threads to use for BLAS operations", NULL, &blas_view_flag)); 1115046ed546SBarry Smith #if defined(PETSC_HAVE_BLI_THREAD_SET_NUM_THREADS) || defined(PETSC_HAVE_MKL_SET_NUM_THREADS) || defined(PETSC_HAVE_OPENBLAS_SET_NUM_THREADS) 1116046ed546SBarry Smith { 1117046ed546SBarry Smith char *threads = NULL; 1118046ed546SBarry Smith 1119046ed546SBarry Smith /* determine any default number of threads requested in the environment; TODO: Apple libraries? */ 1120046ed546SBarry Smith #if defined(PETSC_HAVE_BLI_THREAD_SET_NUM_THREADS) 1121046ed546SBarry Smith threads = getenv("BLIS_NUM_THREADS"); 1122046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of BLIS threads %s given by BLIS_NUM_THREADS\n", threads)); 1123046ed546SBarry Smith if (!threads) { 1124046ed546SBarry Smith threads = getenv("OMP_NUM_THREADS"); 1125046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of BLIS threads %s given by OMP_NUM_THREADS\n", threads)); 1126046ed546SBarry Smith } 1127046ed546SBarry Smith #elif defined(PETSC_HAVE_MKL_SET_NUM_THREADS) 1128046ed546SBarry Smith threads = getenv("MKL_NUM_THREADS"); 1129046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of MKL threads %s given by MKL_NUM_THREADS\n", threads)); 113097873304SJunchao Zhang if (!threads) { 113197873304SJunchao Zhang threads = getenv("OMP_NUM_THREADS"); 113297873304SJunchao Zhang if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of MKL threads %s given by OMP_NUM_THREADS\n", threads)); 113397873304SJunchao Zhang } 1134046ed546SBarry Smith #elif defined(PETSC_HAVE_OPENBLAS_SET_NUM_THREADS) 1135046ed546SBarry Smith threads = getenv("OPENBLAS_NUM_THREADS"); 1136046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of OpenBLAS threads %s given by OPENBLAS_NUM_THREADS\n", threads)); 1137046ed546SBarry Smith if (!threads) { 1138046ed546SBarry Smith threads = getenv("OMP_NUM_THREADS"); 1139046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of OpenBLAS threads %s given by OMP_NUM_THREADS\n", threads)); 1140046ed546SBarry Smith } 1141046ed546SBarry Smith #endif 1142046ed546SBarry Smith if (threads) (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumBLASThreads); 1143046ed546SBarry Smith PetscCall(PetscOptionsInt("-blas_num_threads", "Number of threads to use for BLAS operations", "None", PetscNumBLASThreads, &PetscNumBLASThreads, &flg)); 114497873304SJunchao Zhang if (flg) PetscCall(PetscInfo(NULL, "BLAS: Command line number of BLAS thread %" PetscInt_FMT "given by -blas_num_threads\n", PetscNumBLASThreads)); 114597873304SJunchao Zhang if (flg || threads) { 1146046ed546SBarry Smith PetscCall(PetscBLASSetNumThreads(PetscNumBLASThreads)); 1147046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: number of threads %" PetscInt_FMT "\n", PetscNumBLASThreads)); 1148046ed546SBarry Smith } 114997873304SJunchao Zhang } 1150046ed546SBarry Smith #elif defined(PETSC_HAVE_APPLE_ACCELERATE) 1151046ed546SBarry Smith PetscCall(PetscInfo(NULL, "BLAS: Apple Accelerate library, thread support with no user control\n")); 1152046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: Apple Accelerate library, thread support with no user control\n")); 1153046ed546SBarry Smith #else 1154046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: no thread support\n")); 1155046ed546SBarry Smith #endif 1156046ed546SBarry Smith PetscOptionsEnd(); 1157046ed546SBarry Smith 115885649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 115985649d77SJunchao Zhang /* 116085649d77SJunchao Zhang Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 116185649d77SJunchao Zhang 116285649d77SJunchao Zhang Currently not used because it is not supported by MPICH. 116385649d77SJunchao Zhang */ 11649566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL)); 116585649d77SJunchao Zhang #endif 116685649d77SJunchao Zhang 116785649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS) 11689566063dSJacob Faibussowitsch PetscCall(PetscFPTCreate(10000)); 116985649d77SJunchao Zhang #endif 117085649d77SJunchao Zhang 117185649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC) 117285649d77SJunchao Zhang { 117385649d77SJunchao Zhang PetscViewer viewer; 1174648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg)); 117585649d77SJunchao Zhang if (flg) { 11769566063dSJacob Faibussowitsch PetscCall(PetscProcessPlacementView(viewer)); 1177648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 117885649d77SJunchao Zhang } 117985649d77SJunchao Zhang } 118085649d77SJunchao Zhang #endif 118185649d77SJunchao Zhang 118285649d77SJunchao Zhang flg = PETSC_TRUE; 11839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL)); 1184648c30bcSBarry Smith if (!flg) PetscCall(PetscOptionsPushCreateViewerOff(PETSC_TRUE)); 118585649d77SJunchao Zhang 118685649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS) 11873ba16761SJacob Faibussowitsch PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD); 11883ba16761SJacob Faibussowitsch PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default); 11893ba16761SJacob Faibussowitsch PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", ""); 11903ba16761SJacob Faibussowitsch PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, ""); 119185649d77SJunchao Zhang #endif 119285649d77SJunchao Zhang 119385649d77SJunchao Zhang #if defined(__VALGRIND_H) 119485649d77SJunchao Zhang PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE; 119585649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE) 1196337bb527SBarry 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")); 119785649d77SJunchao Zhang #endif 119885649d77SJunchao Zhang #endif 119985649d77SJunchao Zhang /* 120085649d77SJunchao Zhang Set flag that we are completely initialized 120185649d77SJunchao Zhang */ 120285649d77SJunchao Zhang PetscInitializeCalled = PETSC_TRUE; 120385649d77SJunchao Zhang 12049566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg)); 12059566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPythonInitialize(NULL, NULL)); 1206f1f2ae84SBarry Smith 1207f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 12089f0612e4SBarry Smith if (flg) PetscCall(PetscInfo(NULL, "Running MPI Linear Solver Server\n")); 1209f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin()); 1210f1f2ae84SBarry 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"); 12113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 121285649d77SJunchao Zhang } 121385649d77SJunchao Zhang 121410450e9eSJacob Faibussowitsch // "Unknown section 'Environmental Variables'" 121510450e9eSJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-section-header-unknown 1216e5c89e4eSSatish Balay /*@C 1217e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 1218811af0c4SBarry Smith `PetscInitialize()` calls MPI_Init() if that has yet to be called, 1219e5c89e4eSSatish Balay so this routine should always be called near the beginning of 1220e5c89e4eSSatish Balay your program -- usually the very first line! 1221e5c89e4eSSatish Balay 1222811af0c4SBarry Smith Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set 1223e5c89e4eSSatish Balay 1224e5c89e4eSSatish Balay Input Parameters: 1225e5c89e4eSSatish Balay + argc - count of number of command line arguments 1226e5c89e4eSSatish Balay . args - the command line arguments 1227be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format. 12280b4b7b1cSBarry Smith Use `NULL` or empty string to not check for code specific file. 12290b4b7b1cSBarry Smith Also checks `~/.petscrc`, `.petscrc` and `petscrc`. 12300b4b7b1cSBarry Smith Use `-skip_petscrc` in the code specific file (or command line) to skip `~/.petscrc`, `.petscrc` and `petscrc` files. 12310b4b7b1cSBarry Smith - help - [optional] Help message to print, use `NULL` for no message 1232e5c89e4eSSatish Balay 1233811af0c4SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that 12340b4b7b1cSBarry Smith communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. 1235811af0c4SBarry Smith then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even 123605827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 1237e5c89e4eSSatish Balay 1238e5c89e4eSSatish Balay Options Database Keys: 12390b4b7b1cSBarry Smith + -help [intro] - prints help method for each option; if `intro` is given the program stops after printing the introductory help message 12407ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 1241e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 12420b4b7b1cSBarry Smith . -on_error_emacs <machinename> - causes `emacsclient` to jump to error file if an error is detected 1243811af0c4SBarry Smith . -on_error_abort - calls `abort()` when error detected (no traceback) 1244811af0c4SBarry Smith . -on_error_mpiabort - calls `MPI_abort()` when error detected 12450b4b7b1cSBarry Smith . -error_output_stdout - prints PETSc error messages to `stdout` instead of the default `stderr` 12468a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 12470b4b7b1cSBarry Smith . -debugger_ranks [rank1,rank2,...] - Indicates MPI ranks to start in debugger 12480b4b7b1cSBarry Smith . -debugger_pause [sleeptime] (in seconds) - Pauses debugger, use if it takes a long time for the debugger to start up on your system 1249e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 12500b4b7b1cSBarry Smith process and wait (`-debugger_pause`) seconds for attachment 1251aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 12520b4b7b1cSBarry Smith . -malloc_test - like `-malloc_dump` `-malloc_debug`, only active for debugging build, ignored in optimized build. Often set in `PETSC_OPTIONS` environmental variable 1253811af0c4SBarry Smith . -malloc_view - show a list of all allocated memory during `PetscFinalize()` 12540b4b7b1cSBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with `-malloc_view` 12550b4b7b1cSBarry Smith . -malloc_requested_size - malloc logging will record the requested size rather than (possibly large) size after alignment 125692f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 1257e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 12580b4b7b1cSBarry Smith . -shared_tmp - indicates `/tmp` directory is known to be shared by all processors 12590b4b7b1cSBarry Smith . -not_shared_tmp - indicates each processor has own `/tmp` 12600b4b7b1cSBarry Smith . -tmp - alternative directory to use instead of `/tmp` 12615955b024SMatthew Knepley . -python <exe> - Initializes Python, and optionally takes a Python executable name 12620b4b7b1cSBarry Smith - -mpiuni-allow-multiprocess-launch - allow `mpiexec` to launch multiple independent MPI-Uni jobs, otherwise a sanity check error is invoked to prevent misuse of MPI-Uni 1263e5c89e4eSSatish Balay 1264c5b5d8d5SVaclav Hapla Options Database Keys for Option Database: 12650b4b7b1cSBarry Smith + -skip_petscrc - skip the default option files `~/.petscrc`, `.petscrc`, `petscrc` 1266c5b5d8d5SVaclav Hapla . -options_monitor - monitor all set options to standard output for the whole program run 1267811af0c4SBarry Smith - -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()` 1268c5b5d8d5SVaclav Hapla 1269c5b5d8d5SVaclav Hapla Options -options_monitor_{all,cancel} are 1270c5b5d8d5SVaclav Hapla position-independent and apply to all options set since the PETSc start. 1271c5b5d8d5SVaclav Hapla They can be used also in option files. 1272c5b5d8d5SVaclav Hapla 1273811af0c4SBarry Smith See `PetscOptionsMonitorSet()` to do monitoring programmatically. 1274c5b5d8d5SVaclav Hapla 1275e5c89e4eSSatish Balay Options Database Keys for Profiling: 1276a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 1277811af0c4SBarry Smith + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`. 1278217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 1279217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 1280495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 1281811af0c4SBarry Smith hangs without running in the debugger). See `PetscLogTraceBegin()`. 1282ad2e3d55SToby Isaac . -log_view [:filename:format][,[:filename:format]...] - Prints summary of flop and timing information to screen or file, see `PetscLogView()` (up to 4 viewers) 1283811af0c4SBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`. 1284811af0c4SBarry Smith . -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView(). 1285f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 1286b665b14eSToby Isaac . -log [filename] - Logs profiling information in a dump file, see `PetscLogDump()`. 1287b665b14eSToby Isaac . -log_all [filename] - Same as `-log`. 1288c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 12892611ad71SToby Isaac . -log_perfstubs - Starts a log handler with the perfstubs interface (which is used by TAU) 129061cc7448SToby Isaac . -log_nvtx - Starts an nvtx log handler for use with Nsight 129156a72328SZach Atkins . -log_roctx - Starts an roctx log handler for use with rocprof on AMD GPUs 1292811af0c4SBarry Smith . -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off 1293983cdb30SBarry Smith . -get_total_flops - Returns total flops done by all processors 1294983cdb30SBarry Smith . -memory_view - Print memory usage at end of run 1295c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 1296495fc317SBarry Smith 1297ffbd1cfbSBarry Smith Options Database Keys for SAWs: 1298ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 1299ffbd1cfbSBarry 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 1300ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 1301ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 1302ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 1303ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 1304ffbd1cfbSBarry Smith 1305e5c89e4eSSatish Balay Environmental Variables: 13060b4b7b1cSBarry Smith + `PETSC_TMP` - alternative directory to use instead of `/tmp` 13070b4b7b1cSBarry Smith . `PETSC_SHARED_TMP` - `/tmp` is shared by all processes 13080b4b7b1cSBarry Smith . `PETSC_NOT_SHARED_TMP` - each process has its own private `/tmp` 1309f0b74427SPierre Jolivet . `PETSC_OPTIONS` - a string containing additional options for PETSc in the form of command line "-key value" pairs 13100b4b7b1cSBarry Smith . `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml with `--download-yaml`) a string containing additional options for PETSc in the form of a YAML document 1311811af0c4SBarry Smith . `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer 1312811af0c4SBarry Smith - `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to 1313e5c89e4eSSatish Balay 1314e5c89e4eSSatish Balay Level: beginner 1315e5c89e4eSSatish Balay 1316811af0c4SBarry Smith Note: 13170b4b7b1cSBarry Smith If for some reason you must call `MPI_Init()` separately from `PetscInitialize()`, call 1318811af0c4SBarry Smith it before `PetscInitialize()`. 1319e5c89e4eSSatish Balay 1320811af0c4SBarry Smith Fortran Notes: 1321811af0c4SBarry Smith In Fortran this routine can be called with 1322811af0c4SBarry Smith .vb 1323811af0c4SBarry Smith call PetscInitialize(ierr) 1324811af0c4SBarry Smith call PetscInitialize(file,ierr) or 1325811af0c4SBarry Smith call PetscInitialize(file,help,ierr) 1326811af0c4SBarry Smith .ve 1327e5c89e4eSSatish Balay 1328811af0c4SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after 1329811af0c4SBarry Smith calling `PetscInitialize()`. 1330e5c89e4eSSatish Balay 13315dedd997SBarry Smith Options Database Key for Developers: 13325dedd997SBarry Smith . -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form: 13335dedd997SBarry Smith "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py) 13345dedd997SBarry Smith 1335db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()` 1336e5c89e4eSSatish Balay @*/ 1337d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[]) 1338d71ae5a4SJacob Faibussowitsch { 133985649d77SJunchao Zhang PetscMPIInt flag; 1340983cdb30SBarry Smith const char *prog = "Unknown Name"; 1341e5c89e4eSSatish Balay 134227104ee2SJacob Faibussowitsch PetscFunctionBegin; 13433ba16761SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS); 13449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Initialized(&flag)); 1345e5c89e4eSSatish Balay if (!flag) { 134639a651e2SJacob 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"); 13479566063dSJacob Faibussowitsch PetscCall(PetscPreMPIInit_Private()); 13485e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 13495e765c61SJed Brown { 1350835d5ba2SJunchao Zhang PetscMPIInt provided; 1351835d5ba2SJunchao Zhang PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE ? MPI_THREAD_FUNNELED : PETSC_MPI_THREAD_REQUIRED, &provided)); 1352835d5ba2SJunchao 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"); 1353835d5ba2SJunchao Zhang if (PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE) PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED; // assign it a valid value after check-up 13545e765c61SJed Brown } 13555e765c61SJed Brown #else 13569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init(argc, args)); 13575e765c61SJed Brown #endif 1358e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 1359e5c89e4eSSatish Balay } 13604dfee713SSatish Balay 136185649d77SJunchao Zhang if (argc && *argc) prog = **args; 1362e5c89e4eSSatish Balay if (argc && args) { 1363e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 1364e5c89e4eSSatish Balay PetscGlobalArgs = *args; 1365e5c89e4eSSatish Balay } 1366daa8fb3bSBarry Smith PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, 0)); 13673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1368e5c89e4eSSatish Balay } 1369e5c89e4eSSatish Balay 137095c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1371ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1372ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 137395c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 1374e5c89e4eSSatish Balay 1375008a6e76SBarry Smith /* 1376008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1377008a6e76SBarry Smith */ 1378d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void) 1379d71ae5a4SJacob Faibussowitsch { 1380008a6e76SBarry Smith PetscFunctionBegin; 1381ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128) 13829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128)); 138374c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 13849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128)); 1385008a6e76SBarry Smith #endif 138674c01ffaSSatish Balay #endif 1387ef384fdeSPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16) 13889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FP16)); 1389008a6e76SBarry Smith #endif 1390008a6e76SBarry Smith 1391de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 13929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_SUM)); 1393613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MAX)); 1394613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MIN)); 1395ef384fdeSPierre Jolivet #elif (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16)) 13969e517322SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128)); 1397008a6e76SBarry Smith #endif 1398008a6e76SBarry Smith 13999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR)); 14009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT)); 14019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT)); 140240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 14039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2INT)); 14046497c311SBarry Smith PetscCallMPI(MPI_Type_free(&MPIU_INT_MPIINT)); 1405008a6e76SBarry Smith #endif 14069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPI_4INT)); 14079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_4INT)); 14089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP)); 140962e5d2d2SJDBetteridge PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp)); 14103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1411008a6e76SBarry Smith } 1412008a6e76SBarry Smith 1413a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void); 14149f0612e4SBarry Smith PETSC_EXTERN PetscErrorCode PetscFreeAlign(void *, int, const char[], const char[]); 1415a4af0ceeSJacob Faibussowitsch 14169f0612e4SBarry Smith /*@ 14170b4b7b1cSBarry Smith PetscFinalize - Checks for options to be called at the conclusion of a PETSc program and frees any remaining PETSc objects and data structures. 14180b4b7b1cSBarry Smith of the program. Automatically calls `MPI_Finalize()` if the user had not called `MPI_Init()` before calling `PetscInitialize()`. 1419e5c89e4eSSatish Balay 1420811af0c4SBarry Smith Collective on `PETSC_COMM_WORLD` 1421e5c89e4eSSatish Balay 1422e5c89e4eSSatish Balay Options Database Keys: 1423811af0c4SBarry Smith + -options_view - Calls `PetscOptionsView()` 1424e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 14257eb1d149SBarry 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 1426e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 1427811af0c4SBarry Smith . -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed 14284bd3d7f8SBarry Smith . -memory_view - Prints total memory usage 14294bd3d7f8SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and in what functions 1430e5c89e4eSSatish Balay 1431e5c89e4eSSatish Balay Level: beginner 1432e5c89e4eSSatish Balay 1433e5c89e4eSSatish Balay Note: 1434811af0c4SBarry Smith See `PetscInitialize()` for other runtime options. 1435e5c89e4eSSatish Balay 14360b4b7b1cSBarry Smith You can call `PetscInitialize()` after `PetscFinalize()` but only with MPI-Uni or if you called `MPI_Init()` before ever calling `PetscInitialize()`. 14370b4b7b1cSBarry Smith 1438db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()` 1439e5c89e4eSSatish Balay @*/ 1440d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void) 1441d71ae5a4SJacob Faibussowitsch { 14424bb5149bSJed Brown PetscMPIInt rank; 1443a8d2bbe5SBarry Smith PetscInt nopt; 14442bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE; 1445dff31646SBarry Smith PetscBool flg; 144610463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 1447e5c89e4eSSatish Balay 14483cbbc5ffSBarry Smith PetscFunctionBegin; 144939a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()"); 14509566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PetscFinalize() called\n")); 1451b022a5c1SBarry Smith 1452f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1453f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd()); 1454f1f2ae84SBarry Smith 14559f0612e4SBarry Smith PetscCall(PetscFreeAlign(PetscGlobalArgsFortran, 0, NULL, NULL)); 14569f0612e4SBarry Smith PetscGlobalArgc = 0; 14579f0612e4SBarry Smith PetscGlobalArgs = NULL; 14589f0612e4SBarry Smith 145962e5d2d2SJDBetteridge /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */ 146062e5d2d2SJDBetteridge { 146162e5d2d2SJDBetteridge union 146262e5d2d2SJDBetteridge { 146362e5d2d2SJDBetteridge MPI_Comm comm; 146462e5d2d2SJDBetteridge void *ptr; 146562e5d2d2SJDBetteridge } ucomm; 146662e5d2d2SJDBetteridge PetscMPIInt flg; 146762e5d2d2SJDBetteridge void *tmp; 146862e5d2d2SJDBetteridge 146962e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 147062e5d2d2SJDBetteridge if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg)); 147162e5d2d2SJDBetteridge if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF)); 147262e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 147362e5d2d2SJDBetteridge if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg)); 147462e5d2d2SJDBetteridge if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD)); 147562e5d2d2SJDBetteridge } 147662e5d2d2SJDBetteridge 14779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1478a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 14793ba16761SJacob Faibussowitsch PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE); 14803ba16761SJacob Faibussowitsch PetscCallExternal(adios_finalize, rank); 1481a56f64adSBarry Smith #endif 14829566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg)); 1483dff31646SBarry Smith if (flg) { 14841f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 14851f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 14861f817a21SBarry Smith 14879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL)); 148848a46eb9SPierre Jolivet if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd)); 14899566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits)); 1490dff31646SBarry Smith cits[0] = 0; 14919566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits)); 14929566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n")); 14939566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 14949566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits)); 14959566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 14969566063dSJacob Faibussowitsch PetscCall(PetscFClose(PETSC_COMM_WORLD, fd)); 14979566063dSJacob Faibussowitsch PetscCall(PetscFree(cits)); 1498dff31646SBarry Smith } 14999566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&PetscCitationsList)); 1500dff31646SBarry Smith 15012d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 15029566063dSJacob Faibussowitsch PetscCall(PetscFPTDestroy()); 15032d53ad75SBarry Smith #endif 15042d53ad75SBarry Smith 1505681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1506681455b2SBarry Smith flg1 = PETSC_FALSE; 15079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL)); 1508681455b2SBarry Smith if (flg1) { 1509681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 1510a5af8288SJose E. Roman PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -15 Xvfb", "r", NULL)); 1511681455b2SBarry Smith } 1512681455b2SBarry Smith #endif 1513681455b2SBarry Smith 151467584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 15159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL)); 151648a46eb9SPierre Jolivet if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n")); 151767584ceeSBarry Smith #endif 1518e5c89e4eSSatish Balay 15192611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 152090d69ab7SBarry Smith flg1 = PETSC_FALSE; 15219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL)); 1522e5c89e4eSSatish Balay if (flg1) { 1523e5c89e4eSSatish Balay PetscLogDouble flops = 0; 15249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD)); 15259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops)); 1526e5c89e4eSSatish Balay } 15272611ad71SToby Isaac } 1528e5c89e4eSSatish Balay 15292611ad71SToby Isaac if (PetscDefined(USE_LOG) && PetscDefined(HAVE_MPE)) { 1530e5c89e4eSSatish Balay mname[0] = 0; 15319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1)); 15322611ad71SToby Isaac if (flg1) PetscCall(PetscLogMPEDump(mname[0] ? mname : NULL)); 1533e5c89e4eSSatish Balay } 1534a297a907SKarl Rupp 1535c9903f8fSJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 1536f0b74427SPierre Jolivet // Free PETSc/kokkos stuff before the potentially non-null PETSc default gpu stream is destroyed by PetscObjectRegisterDestroyAll 1537c9903f8fSJunchao Zhang if (PetscKokkosInitialized) { 1538c9903f8fSJunchao Zhang PetscCall(PetscKokkosFinalize_Private()); 1539c9903f8fSJunchao Zhang PetscKokkosInitialized = PETSC_FALSE; 1540c9903f8fSJunchao Zhang } 1541c9903f8fSJunchao Zhang #endif 1542c9903f8fSJunchao Zhang 15433048253cSJacob Faibussowitsch // Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 15449566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1545dd710f27SBarry Smith 15462611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 1547648c30bcSBarry Smith PetscCall(PetscOptionsPushCreateViewerOff(PETSC_FALSE)); 15489566063dSJacob Faibussowitsch PetscCall(PetscLogViewFromOptions()); 1549648c30bcSBarry Smith PetscCall(PetscOptionsPopCreateViewerOff()); 1550473903fcSJunchao Zhang // It should be turned on with PetscLogGpuTime() and never turned off except in this place 1551473903fcSJunchao Zhang PetscLogGpuTimeFlag = PETSC_FALSE; 155287c3beb6SLisandro Dalcin 15533048253cSJacob Faibussowitsch // Free any objects created by the last block of code. 15549566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1555dd710f27SBarry Smith 1556dd710f27SBarry Smith mname[0] = 0; 15579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1)); 15589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2)); 15599566063dSJacob Faibussowitsch if (flg1 || flg2) PetscCall(PetscLogDump(mname)); 15602611ad71SToby Isaac } 156110463e74SBarry Smith 156290d69ab7SBarry Smith flg1 = PETSC_FALSE; 15639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL)); 15649566063dSJacob Faibussowitsch if (!flg1) PetscCall(PetscPopSignalHandler()); 156590d69ab7SBarry Smith flg1 = PETSC_FALSE; 15669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL)); 15671baa6e33SBarry Smith if (flg1) PetscCall(PetscMPIDump(stdout)); 156890d69ab7SBarry Smith flg1 = PETSC_FALSE; 156990d69ab7SBarry Smith flg2 = PETSC_FALSE; 15708bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 15719566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1)); 15729566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 15739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL)); 1574e4c476e2SSatish Balay 15753a7d0413SPierre Jolivet if (flg2) PetscCall(PetscOptionsView(NULL, PETSC_VIEWER_STDOUT_WORLD)); 1576e5c89e4eSSatish Balay 1577e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 15789566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1)); 15799566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1)); 1580e5c89e4eSSatish Balay 158133fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 15829566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1)); 158359f199a7SJed Brown if (!flg1) flg3 = PETSC_TRUE; 1584e5c89e4eSSatish Balay if (flg3) { 15853de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 158657b1f488SBarry Smith PetscCall(PetscOptionsView(NULL, PETSC_VIEWER_STDOUT_WORLD)); 1587e5c89e4eSSatish Balay } 15889566063dSJacob Faibussowitsch PetscCall(PetscOptionsAllUsed(NULL, &nopt)); 15893de2bfdfSBarry Smith if (nopt) { 15909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n")); 15919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n")); 15923de2bfdfSBarry Smith if (nopt == 1) { 15939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n")); 1594e5c89e4eSSatish Balay } else { 15959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt)); 1596e5c89e4eSSatish Balay } 15973de2bfdfSBarry Smith } else if (flg3 && flg1) { 15989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n")); 1599df12ba86SBarry Smith } 16009566063dSJacob Faibussowitsch PetscCall(PetscOptionsLeft(NULL)); 1601e5c89e4eSSatish Balay } 1602e5c89e4eSSatish Balay 1603e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1604d45a07a7SBarry Smith if (!PetscGlobalRank) { 16059566063dSJacob Faibussowitsch PetscCall(PetscStackSAWsViewOff()); 1606792fecdfSBarry Smith PetscCallSAWs(SAWs_Finalize, ()); 1607d45a07a7SBarry Smith } 1608ec957eceSBarry Smith #endif 1609ec957eceSBarry Smith 161010463e74SBarry Smith /* 1611dbc8283eSBarry Smith List all objects the user may have forgot to free 16122eff7a51SBarry Smith */ 16132611ad71SToby Isaac if (PetscDefined(USE_LOG) && PetscObjectsLog) { 16149566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 1615a64a8e02SBarry Smith if (flg1) { 1616a64a8e02SBarry Smith MPI_Comm local_comm; 16177eb1d149SBarry Smith char string[64]; 1618a64a8e02SBarry Smith 16199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL)); 1620afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16219566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16229566063dSJacob Faibussowitsch PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE)); 16239566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 16250a1571b3SBarry Smith } 162605df10baSBarry Smith } 16274097062eSBarry Smith 1628dbc8283eSBarry Smith PetscObjectsCounts = 0; 1629dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 16309566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscObjects)); 16312eff7a51SBarry Smith 163233f85c2fSBarry Smith /* 163333f85c2fSBarry Smith Destroy any packages that registered a finalize 163433f85c2fSBarry Smith */ 16359566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalizeAll()); 163633f85c2fSBarry Smith 16379566063dSJacob Faibussowitsch PetscCall(PetscLogFinalize()); 1638101409b8SToby Isaac 163933f85c2fSBarry Smith /* 164048dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 164148dd1dffSBarry Smith */ 16422e956fe4SStefano Zampini if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll()); 164337e93019SBarry Smith 16444028d114SSatish Balay if (petsc_history) { 16459566063dSJacob Faibussowitsch PetscCall(PetscCloseHistoryFile(&petsc_history)); 164602c9f0b5SLisandro Dalcin petsc_history = NULL; 1647e5c89e4eSSatish Balay } 16489566063dSJacob Faibussowitsch PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton)); 16499566063dSJacob Faibussowitsch PetscCall(PetscInfoDestroy()); 1650e5c89e4eSSatish Balay 165167584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 165292f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1653e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 165492f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1655e5c89e4eSSatish Balay FILE *fd; 1656ed9cf6e9SBarry Smith int err; 1657e5c89e4eSSatish Balay 1658dc92acbaSJed Brown flg2 = PETSC_FALSE; 165992f119d6SBarry Smith flg3 = PETSC_FALSE; 16609566063dSJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL)); 16619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL)); 166292f119d6SBarry Smith fname[0] = 0; 16639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1)); 1664e5c89e4eSSatish Balay if (flg1 && fname[0]) { 16653ba16761SJacob Faibussowitsch PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank)); 16669371c9d4SSatish Balay fd = fopen(sname, "w"); 16679371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16689566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(fd)); 1669ed9cf6e9SBarry Smith err = fclose(fd); 167028b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 167192f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1672e5c89e4eSSatish Balay MPI_Comm local_comm; 1673e5c89e4eSSatish Balay 1674afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16759566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16769566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(stdout)); 16779566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1679e5c89e4eSSatish Balay } 1680e5c89e4eSSatish Balay fname[0] = 0; 16819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1)); 1682e5c89e4eSSatish Balay if (flg1 && fname[0]) { 16833ba16761SJacob Faibussowitsch PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank)); 16849371c9d4SSatish Balay fd = fopen(sname, "w"); 16859371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16869566063dSJacob Faibussowitsch PetscCall(PetscMallocView(fd)); 1687ed9cf6e9SBarry Smith err = fclose(fd); 168828b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 168992f119d6SBarry Smith } else if (flg1) { 169092f119d6SBarry Smith MPI_Comm local_comm; 169192f119d6SBarry Smith 1692afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16939566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16949566063dSJacob Faibussowitsch PetscCall(PetscMallocView(stdout)); 16959566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1697e5c89e4eSSatish Balay } 1698e5c89e4eSSatish Balay } 169967584ceeSBarry Smith #endif 170020e2c332SMatthew G. Knepley 17015486ca60SMatthew G. Knepley /* 17025486ca60SMatthew G. Knepley Close any open dynamic libraries 17035486ca60SMatthew G. Knepley */ 17049566063dSJacob Faibussowitsch PetscCall(PetscFinalize_DynamicLibraries()); 17055486ca60SMatthew G. Knepley 1706e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 17079566063dSJacob Faibussowitsch PetscCall(PetscOptionsDestroyDefault()); 1708e5c89e4eSSatish Balay 170971438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM) 171071438e86SJunchao Zhang if (PetscBeganNvshmem) { 17119566063dSJacob Faibussowitsch PetscCall(PetscNvshmemFinalize()); 171271438e86SJunchao Zhang PetscBeganNvshmem = PETSC_FALSE; 171371438e86SJunchao Zhang } 171471438e86SJunchao Zhang #endif 1715a0e72f99SJunchao Zhang 17169566063dSJacob Faibussowitsch PetscCall(PetscFreeMPIResources()); 1717e5c89e4eSSatish Balay 1718dbc8283eSBarry Smith /* 1719efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1720efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1721efb80d3cSBarry Smith 1722efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1723efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1724dbc8283eSBarry Smith */ 1725b770b1f6SSatish Balay { 1726dbc8283eSBarry Smith PetscCommCounter *counter; 1727dbc8283eSBarry Smith PetscMPIInt flg; 1728dbc8283eSBarry Smith MPI_Comm icomm; 17299371c9d4SSatish Balay union 17309371c9d4SSatish Balay { 17319371c9d4SSatish Balay MPI_Comm comm; 17329371c9d4SSatish Balay void *ptr; 17339371c9d4SSatish Balay } ucomm; 17349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 1735dbc8283eSBarry Smith if (flg) { 1736265f3f35SJed Brown icomm = ucomm.comm; 17379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 173828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1739dbc8283eSBarry Smith 17409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval)); 17419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 17429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1743dbc8283eSBarry Smith } 17449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 1745dbc8283eSBarry Smith if (flg) { 1746265f3f35SJed Brown icomm = ucomm.comm; 17479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 174828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1749dbc8283eSBarry Smith 17509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval)); 17519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 17529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1753dbc8283eSBarry Smith } 1754b770b1f6SSatish Balay } 1755dbc8283eSBarry Smith 17569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval)); 17579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval)); 17589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval)); 17599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval)); 176062e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval)); 176162e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval)); 1762480cf27aSJed Brown 17635ea2b939SDuncan Campbell // Free keyvals which may be silently created by some routines 17645ea2b939SDuncan Campbell if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval)); 17655ea2b939SDuncan Campbell if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval)); 17665ea2b939SDuncan Campbell 17679566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen)); 17689566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout)); 17699566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr)); 17709566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock)); 1771ef19f930SBarry Smith 1772e5c89e4eSSatish Balay if (PetscBeganMPI) { 177399b1327fSBarry Smith PetscMPIInt flag; 17749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalized(&flag)); 177539a651e2SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 177639a651e2SJacob Faibussowitsch /* wait until the very last moment to disable error handling */ 177739a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_FALSE; 17789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalize()); 177939a651e2SJacob Faibussowitsch } else PetscErrorHandlingInitialized = PETSC_FALSE; 178039a651e2SJacob Faibussowitsch 1781e5c89e4eSSatish Balay /* 1782e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1783e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1784e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1785e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1786e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 17870e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1788e5c89e4eSSatish Balay memory was not freed. 1789e5c89e4eSSatish Balay 1790e5c89e4eSSatish Balay */ 17919566063dSJacob Faibussowitsch PetscCall(PetscMallocClear()); 17929566063dSJacob Faibussowitsch PetscCall(PetscStackReset()); 1793a297a907SKarl Rupp 1794e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1795e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 17967ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE) 179756883f60SBarry Smith /* 179856883f60SBarry Smith flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the 179956883f60SBarry Smith gcov files are still being added to the directories as git tries to remove the directories. 180056883f60SBarry Smith */ 180156883f60SBarry Smith __gcov_flush(); 180256883f60SBarry Smith #endif 18031724198aSStefano Zampini /* To match PetscFunctionBegin() at the beginning of this function */ 18041724198aSStefano Zampini PetscStackClearTop; 18053ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 1806e5c89e4eSSatish Balay } 1807e5c89e4eSSatish Balay 180843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 1809d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b) 1810d71ae5a4SJacob Faibussowitsch { 181143db4dbbSBarry Smith if (*a == *b) return 1; 181243db4dbbSBarry Smith if (*a + 32 == *b) return 1; 181343db4dbbSBarry Smith if (*a - 32 == *b) return 1; 181443db4dbbSBarry Smith return 0; 181543db4dbbSBarry Smith } 1816a70650f6SBarry Smith #endif 181743db4dbbSBarry Smith 181843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 1819d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b) 1820d71ae5a4SJacob Faibussowitsch { 182143db4dbbSBarry Smith if (*a == *b) return 1; 182243db4dbbSBarry Smith if (*a + 32 == *b) return 1; 182343db4dbbSBarry Smith if (*a - 32 == *b) return 1; 182443db4dbbSBarry Smith return 0; 182543db4dbbSBarry Smith } 182643db4dbbSBarry Smith #endif 18276a210b70SBarry Smith 1828462c564dSBarry Smith static inline PetscMPIInt MPIU_Allreduce_Count(const void *inbuf, void *outbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm) 18296a210b70SBarry Smith { 1830462c564dSBarry Smith PetscMPIInt err; 18316a210b70SBarry Smith #if !defined(PETSC_HAVE_MPI_LARGE_COUNT) 18326a210b70SBarry Smith PetscMPIInt count2; 18336a210b70SBarry Smith 1834462c564dSBarry Smith PetscMPIIntCast_Internal(count, &count2); 1835462c564dSBarry Smith err = MPI_Allreduce((void *)inbuf, outbuf, count2, dtype, op, comm); 18366a210b70SBarry Smith #else 1837462c564dSBarry Smith err = MPI_Allreduce_c((void *)inbuf, outbuf, count, dtype, op, comm); 18386a210b70SBarry Smith #endif 1839462c564dSBarry Smith return err; 18406a210b70SBarry Smith } 18416a210b70SBarry Smith 1842462c564dSBarry Smith /* 1843462c564dSBarry Smith When count is 1 and dtype == MPIU_INT performs the reduction in PetscInt64 to check for integer overflow 1844462c564dSBarry Smith */ 1845462c564dSBarry Smith PetscMPIInt MPIU_Allreduce_Private(const void *inbuf, void *outbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm) 18466a210b70SBarry Smith { 1847462c564dSBarry Smith PetscMPIInt err; 18482711dd09SLisandro Dalcin if (!PetscDefined(USE_64BIT_INDICES) && count == 1 && dtype == MPIU_INT && (op == MPI_SUM || op == MPI_PROD)) { 1849462c564dSBarry Smith PetscInt64 incnt, outcnt; 18506a210b70SBarry Smith void *inbufd, *outbufd; 18516a210b70SBarry Smith 18526a210b70SBarry Smith if (inbuf != MPI_IN_PLACE) { 18536a210b70SBarry Smith incnt = *(PetscInt32 *)inbuf; 18546a210b70SBarry Smith inbufd = &incnt; 18556a210b70SBarry Smith outbufd = &outcnt; 1856462c564dSBarry Smith err = MPIU_Allreduce_Count(inbufd, outbufd, count, MPIU_INT64, op, comm); 1857835f2295SStefano Zampini } else { 1858835f2295SStefano Zampini outcnt = *(PetscInt32 *)outbuf; 1859835f2295SStefano Zampini outbufd = &outcnt; 1860835f2295SStefano Zampini err = MPIU_Allreduce_Count(MPI_IN_PLACE, outbufd, count, MPIU_INT64, op, comm); 1861835f2295SStefano Zampini } 1862462c564dSBarry Smith if (!err && outcnt > PETSC_INT_MAX) err = MPI_ERR_OTHER; 18636a210b70SBarry Smith *(PetscInt32 *)outbuf = (PetscInt32)outcnt; 18646a210b70SBarry Smith } else { 1865462c564dSBarry Smith err = MPIU_Allreduce_Count(inbuf, outbuf, count, dtype, op, comm); 18666a210b70SBarry Smith } 1867462c564dSBarry Smith return err; 18686a210b70SBarry Smith } 186949abdd8aSBarry Smith 1870*458b0db5SMartin Diehl // Check if MPIU_Allreduce() is called on the same filename:lineno and with the same data count across all processes. Error out if otherwise. 187126139085SJunchao Zhang PetscErrorCode PetscCheckAllreduceSameLineAndCount_Private(MPI_Comm comm, const char *filename, PetscMPIInt lineno, PetscMPIInt count) 187226139085SJunchao Zhang { 1873*458b0db5SMartin Diehl PetscMPIInt rbuf[4]; 187426139085SJunchao Zhang 187526139085SJunchao Zhang PetscFunctionBegin; 1876*458b0db5SMartin Diehl rbuf[0] = lineno; 1877*458b0db5SMartin Diehl rbuf[1] = -rbuf[0]; 1878*458b0db5SMartin Diehl rbuf[2] = count; 1879*458b0db5SMartin Diehl rbuf[3] = -rbuf[2]; 1880*458b0db5SMartin Diehl PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, rbuf, 4, MPI_INT, MPI_MAX, comm)); 188126139085SJunchao Zhang 188226139085SJunchao Zhang if (rbuf[0] != -rbuf[1]) { 188326139085SJunchao Zhang size_t len; 188426139085SJunchao Zhang PetscMPIInt size, rank, ilen, *recvcounts = NULL, *displs = NULL; 188526139085SJunchao Zhang char *str = NULL, *str0 = NULL; 188626139085SJunchao Zhang 188726139085SJunchao Zhang PetscCallMPI(MPI_Comm_size(comm, &size)); 188826139085SJunchao Zhang PetscCallMPI(MPI_Comm_rank(comm, &rank)); 188926139085SJunchao Zhang PetscCall(PetscStrlen(filename, &len)); 189026139085SJunchao Zhang len += 128; /* add enough space for the leading and trailing chars in PetscSNPrintf around __FILE__ */ 189126139085SJunchao Zhang PetscCall(PetscMalloc1(len, &str)); 189226139085SJunchao Zhang PetscCall(PetscSNPrintf(str, len, " On process %d, %s:%d\n", rank, filename, lineno)); 189326139085SJunchao Zhang PetscCall(PetscStrlen(str, &len)); /* string length exclusive of the NULL terminator */ 189426139085SJunchao Zhang ilen = (PetscMPIInt)len; 189526139085SJunchao Zhang if (rank == 0) PetscCall(PetscMalloc2(size, &recvcounts, size + 1, &displs)); 189626139085SJunchao Zhang PetscCallMPI(MPI_Gather(&ilen, 1, MPI_INT, recvcounts, 1, MPI_INT, 0, comm)); 189726139085SJunchao Zhang if (rank == 0) { 189826139085SJunchao Zhang displs[0] = 0; 189926139085SJunchao Zhang for (PetscMPIInt i = 0; i < size; i++) displs[i + 1] = displs[i] + recvcounts[i]; 190026139085SJunchao Zhang PetscCall(PetscMalloc1(displs[size], &str0)); 190126139085SJunchao Zhang } 190226139085SJunchao Zhang PetscCallMPI(MPI_Gatherv(str, ilen, MPI_CHAR, str0, recvcounts, displs, MPI_CHAR, 0, comm)); 190326139085SJunchao Zhang if (rank == 0) str0[displs[size] - 1] = 0; /* replace the ending \n with NULL */ 190426139085SJunchao Zhang PetscCall(PetscFree(str)); 190526139085SJunchao Zhang if (rank == 0) PetscCall(PetscFree2(recvcounts, displs)); 190626139085SJunchao Zhang SETERRQ(comm, PETSC_ERR_PLIB, "MPIU_Allreduce() called in different locations on different processes:\n%s", str0); 190726139085SJunchao Zhang } 190826139085SJunchao Zhang PetscCheck(rbuf[2] == -rbuf[3], comm, PETSC_ERR_PLIB, "MPIU_Allreduce() called with different counts %d on different processes", count); 190926139085SJunchao Zhang PetscFunctionReturn(PETSC_SUCCESS); 191026139085SJunchao Zhang } 191126139085SJunchao Zhang 191249abdd8aSBarry Smith /*@C 191349abdd8aSBarry Smith PetscCtxDestroyDefault - An implementation of a `PetscCtxDestroyFn` that uses `PetscFree()` to free the context 191449abdd8aSBarry Smith 191549abdd8aSBarry Smith Input Parameter: 191649abdd8aSBarry Smith . ctx - the context to be destroyed 191749abdd8aSBarry Smith 191849abdd8aSBarry Smith Level: intermediate 191949abdd8aSBarry Smith 192049abdd8aSBarry Smith Note: 192149abdd8aSBarry Smith This is not called directly, rather it is passed to `DMSetApplicationContextDestroy()`, `PetscContainerSetDestroy()`, 192249abdd8aSBarry Smith `PetscObjectContainterCreate()` and similar routines and then called by the destructor of the associated object. 192349abdd8aSBarry Smith 192449abdd8aSBarry Smith .seealso: `PetscObject`, `PetscCtxDestroyFn`, `PetscObjectDestroy()`, `DMSetApplicationContextDestroy()`, `PetscContainerSetDestroy()`, 192549abdd8aSBarry Smith `PetscObjectContainterCreate()` 192649abdd8aSBarry Smith @*/ 192749abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode PetscCtxDestroyDefault(void **ctx) 192849abdd8aSBarry Smith { 192949abdd8aSBarry Smith PetscFunctionBegin; 193049abdd8aSBarry Smith PetscCall(PetscFree(*ctx)); 193149abdd8aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 193249abdd8aSBarry Smith } 1933