1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay This file defines the initialization of PETSc, including PetscInitialize() 4e5c89e4eSSatish Balay */ 5af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I "petscsys.h" I*/ 6473903fcSJunchao Zhang #include <petsc/private/logimpl.h> 7665c2dedSJed Brown #include <petscviewer.h> 862e5d2d2SJDBetteridge #include <petsc/private/garbagecollector.h> 9a0e72f99SJunchao Zhang 106edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS) 116edef35eSSatish Balay #include <petsc/private/valgrind/valgrind.h> 126edef35eSSatish Balay #endif 136edef35eSSatish Balay 14fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS) 1585649d77SJunchao Zhang #include <petsc/private/fortranimpl.h> 1685649d77SJunchao Zhang #endif 1785649d77SJunchao Zhang 187ce81a4bSJacob Faibussowitsch #if PetscDefined(USE_COVERAGE) 1956883f60SBarry Smith EXTERN_C_BEGIN 20aaf3846cSSatish Balay #if defined(PETSC_HAVE___GCOV_DUMP) 21aaf3846cSSatish Balay #define __gcov_flush(x) __gcov_dump(x) 22aaf3846cSSatish Balay #endif 2356883f60SBarry Smith void __gcov_flush(void); 2456883f60SBarry Smith EXTERN_C_END 2556883f60SBarry Smith #endif 268101f56cSMatthew Knepley 272d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 2895c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData; 292d53ad75SBarry Smith PetscFPT PetscFPTData = 0; 302d53ad75SBarry Smith #endif 312d53ad75SBarry Smith 3227104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 3316ad0300SBarry Smith #include <petscviewersaws.h> 34a6790183SBarry Smith #endif 35eb02082bSJunchao Zhang 3695c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history; 37e5c89e4eSSatish Balay 3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void); 3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void); 4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm, int); 4195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm, int); 4295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE **); 430069ddf5SShri Abhyankar 446de5d289SStefano Zampini /* user may set these BEFORE calling PetscInitialize() */ 45e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL; 4627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_MPI_INIT_THREAD) 47835d5ba2SJunchao Zhang PetscMPIInt PETSC_MPI_THREAD_REQUIRED = PETSC_DECIDE; 486de5d289SStefano Zampini #else 49835d5ba2SJunchao Zhang PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_SINGLE; 506de5d289SStefano Zampini #endif 51e5c89e4eSSatish Balay 52480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval = MPI_KEYVAL_INVALID; 53480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID; 54480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID; 555f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval = MPI_KEYVAL_INVALID; 5662e5d2d2SJDBetteridge PetscMPIInt Petsc_CreationIdx_keyval = MPI_KEYVAL_INVALID; 5762e5d2d2SJDBetteridge PetscMPIInt Petsc_Garbage_HMap_keyval = MPI_KEYVAL_INVALID; 58480cf27aSJed Brown 595ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedWD_keyval = MPI_KEYVAL_INVALID; 605ea2b939SDuncan Campbell PetscMPIInt Petsc_SharedTmp_keyval = MPI_KEYVAL_INVALID; 615ea2b939SDuncan Campbell 62e5c89e4eSSatish Balay /* 63e5c89e4eSSatish Balay Declare and set all the string names of the PETSc enums 64e5c89e4eSSatish Balay */ 6502c9f0b5SLisandro Dalcin const char *const PetscBools[] = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL}; 6602c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL}; 67e5c89e4eSSatish Balay 68ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE; 69ace3abfcSBarry Smith PetscBool PetscPreLoadingOn = PETSC_FALSE; 700f8e0872SSatish Balay 71a2f94806SJed Brown PetscInt PetscHotRegionDepth; 72a2f94806SJed Brown 736edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE; 746edef35eSSatish Balay 75b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY) 76b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen; 77b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout; 78b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr; 79b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock; 80b22622e2STadeu Manoel #endif 81b22622e2STadeu Manoel 82046ed546SBarry Smith extern PetscInt PetscNumBLASThreads; 83046ed546SBarry Smith 84aec76313SJacob Faibussowitsch /*@C 85945d1669SBarry Smith PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args 8672a42c3cSBarry Smith 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 @*/ 129d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitializeNoArguments(void) 130d71ae5a4SJacob Faibussowitsch { 131e5c89e4eSSatish Balay int argc = 0; 13202c9f0b5SLisandro Dalcin char **args = NULL; 133e5c89e4eSSatish Balay 134e5c89e4eSSatish Balay PetscFunctionBegin; 1359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, NULL)); 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137e5c89e4eSSatish Balay } 138e5c89e4eSSatish Balay 139e5c89e4eSSatish Balay /*@ 140e5c89e4eSSatish Balay PetscInitialized - Determine whether PETSc is initialized. 141e5c89e4eSSatish Balay 14210450e9eSJacob Faibussowitsch Output Parameter: 14310450e9eSJacob Faibussowitsch . isInitialized - `PETSC_TRUE` if PETSc is initialized, `PETSC_FALSE` otherwise 14410450e9eSJacob Faibussowitsch 14593b6d2d1SJed Brown Level: beginner 146e5c89e4eSSatish Balay 147db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()` 148e5c89e4eSSatish Balay @*/ 149d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialized(PetscBool *isInitialized) 150d71ae5a4SJacob Faibussowitsch { 15127104ee2SJacob Faibussowitsch PetscFunctionBegin; 1524f572ea9SToby Isaac if (PetscInitializeCalled) PetscAssertPointer(isInitialized, 1); 153e5c89e4eSSatish Balay *isInitialized = PetscInitializeCalled; 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155e5c89e4eSSatish Balay } 156e5c89e4eSSatish Balay 157e5c89e4eSSatish Balay /*@ 158811af0c4SBarry Smith PetscFinalized - Determine whether `PetscFinalize()` has been called yet 159e5c89e4eSSatish Balay 16010450e9eSJacob Faibussowitsch Output Parameter: 16110450e9eSJacob Faibussowitsch . isFinalized - `PETSC_TRUE` if PETSc is finalized, `PETSC_FALSE` otherwise 16210450e9eSJacob Faibussowitsch 163e5c89e4eSSatish Balay Level: developer 164e5c89e4eSSatish Balay 165db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()` 166e5c89e4eSSatish Balay @*/ 167d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalized(PetscBool *isFinalized) 168d71ae5a4SJacob Faibussowitsch { 16927104ee2SJacob Faibussowitsch PetscFunctionBegin; 1704f572ea9SToby Isaac if (!PetscFinalizeCalled) PetscAssertPointer(isFinalized, 1); 171e5c89e4eSSatish Balay *isFinalized = PetscFinalizeCalled; 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173e5c89e4eSSatish Balay } 174e5c89e4eSSatish Balay 17557171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]); 176e5c89e4eSSatish Balay 177e5c89e4eSSatish Balay /* 178e5c89e4eSSatish Balay This function is the MPI reduction operation used to compute the sum of the 179e5c89e4eSSatish Balay first half of the datatype and the max of the second half. 180e5c89e4eSSatish Balay */ 181367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0; 18262e5d2d2SJDBetteridge MPI_Op Petsc_Garbage_SetIntersectOp = 0; 183e5c89e4eSSatish Balay 184*6497c311SBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) 185d71ae5a4SJacob Faibussowitsch { 186e5c89e4eSSatish Balay PetscFunctionBegin; 187*6497c311SBarry Smith if (*datatype == MPIU_INT_MPIINT && PetscDefined(USE_64BIT_INDICES)) { 188*6497c311SBarry Smith #if defined(PETSC_USE_64BIT_INDICES) 189*6497c311SBarry Smith struct petsc_mpiu_int_mpiint *xin = (struct petsc_mpiu_int_mpiint *)in, *xout = (struct petsc_mpiu_int_mpiint *)out; 190*6497c311SBarry Smith PetscMPIInt count = *cnt; 191e5c89e4eSSatish Balay 192*6497c311SBarry Smith for (PetscMPIInt i = 0; i < count; i++) { 193*6497c311SBarry Smith xout[i].a = PetscMax(xout[i].a, xin[i].a); 194*6497c311SBarry Smith xout[i].b += xin[i].b; 195*6497c311SBarry Smith } 196*6497c311SBarry Smith #endif 197*6497c311SBarry Smith } else if (*datatype == MPIU_2INT || *datatype == MPIU_INT_MPIINT) { 198*6497c311SBarry Smith PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out; 199*6497c311SBarry Smith PetscMPIInt count = *cnt; 200*6497c311SBarry Smith 201*6497c311SBarry 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 } 205*6497c311SBarry Smith } else { 206*6497c311SBarry Smith PetscErrorCode ierr = (*PetscErrorPrintf)("Can only handle MPIU_2INT and MPIU_INT_MPIINT data types"); 207*6497c311SBarry Smith (void)ierr; 208*6497c311SBarry Smith PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 209*6497c311SBarry Smith } 210812af9f3SBarry Smith PetscFunctionReturnVoid(); 211e5c89e4eSSatish Balay } 212e5c89e4eSSatish Balay 213*6497c311SBarry Smith /*@ 214*6497c311SBarry Smith PetscMaxSum - Returns the max of the first entry over all MPI processes and the sum of the second entry. 215b693b147SBarry Smith 216*6497c311SBarry Smith Collective 217*6497c311SBarry Smith 218*6497c311SBarry Smith Input Parameters: 219*6497c311SBarry Smith + comm - the communicator 220*6497c311SBarry Smith - array - an arry of length 2 times `size`, the number of MPI processes 221*6497c311SBarry Smith 222*6497c311SBarry Smith Output Parameters: 223*6497c311SBarry Smith + max - the maximum of `array[2*rank]` over all MPI processes 224*6497c311SBarry Smith - sum - the sum of the `array[2*rank + 1]` over all MPI processes 225*6497c311SBarry Smith 226*6497c311SBarry Smith Level: developer 227*6497c311SBarry Smith 228*6497c311SBarry Smith .seealso: `PetscInitialize()` 229*6497c311SBarry Smith @*/ 230*6497c311SBarry 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; 238*6497c311SBarry 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)); 251*6497c311SBarry Smith PetscCall(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 2609e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16) 261613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_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 { 271*6497c311SBarry 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 284613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_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 2959e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) 2969e517322SPierre Jolivet else if (*datatype == MPIU___FP16) { 2979e517322SPierre Jolivet __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out; 298*6497c311SBarry Smith for (i = 0; i < count; i++) xout[i] = (__fp16)(xin[i] + xout[i]); 2999e517322SPierre Jolivet } 3009e517322SPierre Jolivet #endif 3017c2de775SJed Brown else { 3029e517322SPierre Jolivet #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16) 3033ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SElF, (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types")); 3049e517322SPierre Jolivet #elif !defined(PETSC_HAVE_REAL___FP16) 3053ba16761SJacob Faibussowitsch PetscCallAbort(MPI_COMM_SELF, (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types")); 3069e517322SPierre Jolivet #elif !defined(PETSC_HAVE_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; 3809566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm)); 3819566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter->iflags)); 38257f21012SBarry Smith while (comms) { 3839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&comms->comm)); 38457f21012SBarry Smith pcomm = comms; 38557f21012SBarry Smith comms = comms->next; 3869566063dSJacob Faibussowitsch PetscCall(PetscFree(pcomm)); 38757f21012SBarry Smith } 3889566063dSJacob Faibussowitsch PetscCallMPI(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; 41212801b39SBarry Smith if (keyval != Petsc_InnerComm_keyval) SETERRMPI(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; 4229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg)); 42333779a13SJunchao Zhang if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute"); 42433779a13SJunchao Zhang if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm's OuterComm attribute does not point to outer PETSc comm"); 42533779a13SJunchao Zhang } 4269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval)); 4279566063dSJacob Faibussowitsch PetscCallMPI(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; 4379566063dSJacob Faibussowitsch PetscCallMPI(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\ 46930c35bf2SSatish Balay and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\ 47030c35bf2SSatish Balay and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\ 47130c35bf2SSatish Balay Title = {{PETSc/TAO} Users Manual},\n\ 47230694c2dSSatish Balay Number = {ANL-21/39 - Revision 3.21},\n\ 473f5a9e7c7SSatish Balay Doi = {10.2172/2205494},\n\ 47430c35bf2SSatish Balay Institution = {Argonne National Laboratory},\n\ 47530694c2dSSatish Balay Year = {2024}\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 523811af0c4SBarry Smith after PetscInitialize() is called but before `PetscFinalize()`. 524e5c89e4eSSatish Balay 525cc4c1da9SBarry Smith Not Collective, No Fortran Support 526e5c89e4eSSatish Balay 527e5c89e4eSSatish Balay Output Parameters: 528e5c89e4eSSatish Balay + argc - count of 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 53921532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()`, `PetscInitialize()` 540e5c89e4eSSatish Balay @*/ 541d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArgs(int *argc, char ***args) 542d71ae5a4SJacob Faibussowitsch { 543e5c89e4eSSatish Balay PetscFunctionBegin; 54439a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 545e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 546e5c89e4eSSatish Balay *args = PetscGlobalArgs; 5473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 548e5c89e4eSSatish Balay } 549e5c89e4eSSatish Balay 550793721a6SBarry Smith /*@C 551793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 552811af0c4SBarry Smith after `PetscInitialize()` is called but before `PetscFinalize()`. 553793721a6SBarry Smith 554cc4c1da9SBarry Smith Not Collective, No Fortran Support 555793721a6SBarry Smith 5562fe279fdSBarry Smith Output Parameter: 557793721a6SBarry Smith . args - the command line arguments 558793721a6SBarry Smith 559793721a6SBarry Smith Level: intermediate 560793721a6SBarry Smith 56121532e8aSBarry Smith Note: 56221532e8aSBarry Smith This does NOT start with the program name and IS `NULL` terminated (final arg is void) 563793721a6SBarry Smith 56421532e8aSBarry Smith .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()`, `PetscInitialize()` 565793721a6SBarry Smith @*/ 566d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetArguments(char ***args) 567d71ae5a4SJacob Faibussowitsch { 568793721a6SBarry Smith PetscInt i, argc = PetscGlobalArgc; 569793721a6SBarry Smith 570793721a6SBarry Smith PetscFunctionBegin; 57139a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 5729371c9d4SSatish Balay if (!argc) { 5739371c9d4SSatish Balay *args = NULL; 5743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5759371c9d4SSatish Balay } 5769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(argc, args)); 5779566063dSJacob Faibussowitsch for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i])); 5782d747510SLisandro Dalcin (*args)[argc - 1] = NULL; 5793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 580793721a6SBarry Smith } 581793721a6SBarry Smith 582793721a6SBarry Smith /*@C 583811af0c4SBarry Smith PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()` 584793721a6SBarry Smith 585cc4c1da9SBarry Smith Not Collective, No Fortran Support 586793721a6SBarry Smith 5872fe279fdSBarry Smith Output Parameter: 588793721a6SBarry Smith . args - the command line arguments 589793721a6SBarry Smith 590793721a6SBarry Smith Level: intermediate 591793721a6SBarry Smith 592db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()` 593793721a6SBarry Smith @*/ 594d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeArguments(char **args) 595d71ae5a4SJacob Faibussowitsch { 59639a651e2SJacob Faibussowitsch PetscFunctionBegin; 59739a651e2SJacob Faibussowitsch if (args) { 598793721a6SBarry Smith PetscInt i = 0; 599793721a6SBarry Smith 6009566063dSJacob Faibussowitsch while (args[i]) PetscCall(PetscFree(args[i++])); 6019566063dSJacob Faibussowitsch PetscCall(PetscFree(args)); 60239a651e2SJacob Faibussowitsch } 6033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 604793721a6SBarry Smith } 605793721a6SBarry Smith 60627104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 60730befbd2SBarry Smith #include <petscconfiginfo.h> 60830befbd2SBarry Smith 609d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) 610d71ae5a4SJacob Faibussowitsch { 61127104ee2SJacob Faibussowitsch PetscFunctionBegin; 61211525c0dSBarry Smith if (!PetscGlobalRank) { 61330befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64]; 61411525c0dSBarry Smith int port; 615ffbd1cfbSBarry Smith PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE; 61611525c0dSBarry Smith size_t applinelen, introlen; 617ffbd1cfbSBarry Smith char sawsurl[256]; 61811525c0dSBarry Smith 6199566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg)); 62011525c0dSBarry Smith if (flg) { 62111525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 62211525c0dSBarry Smith 6239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL)); 62411525c0dSBarry Smith if (sawslog[0]) { 625792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog)); 62611525c0dSBarry Smith } else { 627792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL)); 62811525c0dSBarry Smith } 62911525c0dSBarry Smith } 6309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg)); 63148a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert)); 6329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL)); 633ffbd1cfbSBarry Smith if (selectport) { 634792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 635792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 636ffbd1cfbSBarry Smith } else { 6379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg)); 63848a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Port, (port)); 639ffbd1cfbSBarry Smith } 6409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg)); 64111525c0dSBarry Smith if (flg) { 642792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 6439566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(root, ".", &rootlocal)); 6449c1e0ce8SBarry Smith } else { 6459566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg)); 6469c1e0ce8SBarry Smith if (flg) { 6479566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root))); 648792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 6499c1e0ce8SBarry Smith } 65011525c0dSBarry Smith } 6519566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2)); 65211525c0dSBarry Smith if (flg2) { 65311525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 65428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option"); 6559566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root)); 6569566063dSJacob Faibussowitsch PetscCall(PetscTestDirectory(jsdir, 'r', &flg)); 65728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory"); 658792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Local_Header, ()); 65911525c0dSBarry Smith } 6609566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(programname, sizeof(programname))); 6619566063dSJacob Faibussowitsch PetscCall(PetscStrlen(help, &applinelen)); 66211525c0dSBarry Smith introlen = 4096 + applinelen; 66330a8c9c0SSurtai Han applinelen += 1024; 6649566063dSJacob Faibussowitsch PetscCall(PetscMalloc(applinelen, &appline)); 6659566063dSJacob Faibussowitsch PetscCall(PetscMalloc(introlen, &intro)); 66611525c0dSBarry Smith 66711525c0dSBarry Smith if (rootlocal) { 6689566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname)); 6699566063dSJacob Faibussowitsch PetscCall(PetscTestFile(appline, 'r', &rootlocal)); 67011525c0dSBarry Smith } 6719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetAll(NULL, &options)); 67211525c0dSBarry Smith if (rootlocal && help) { 6739566063dSJacob 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)); 67411525c0dSBarry Smith } else if (help) { 6759566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help)); 67611525c0dSBarry Smith } else { 6779566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options)); 67811525c0dSBarry Smith } 6799566063dSJacob Faibussowitsch PetscCall(PetscFree(options)); 6809566063dSJacob Faibussowitsch PetscCall(PetscGetVersion(version, sizeof(version))); 6819371c9d4SSatish Balay PetscCall(PetscSNPrintf(intro, introlen, 6829371c9d4SSatish Balay "<body>\n" 683a17b96a8SKyle 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" 684df62c222SBarry 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" 6859371c9d4SSatish Balay "%s", 6869371c9d4SSatish Balay version, petscconfigureoptions, appline)); 687792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro)); 6889566063dSJacob Faibussowitsch PetscCall(PetscFree(intro)); 6899566063dSJacob Faibussowitsch PetscCall(PetscFree(appline)); 690ffbd1cfbSBarry Smith if (selectport) { 691aa573868SBarry Smith PetscBool silent; 6927d812c46SBarry Smith 6937d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 69439a651e2SJacob Faibussowitsch while (SAWs_Initialize()) { 695792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 696792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 6977d812c46SBarry Smith } 6987d812c46SBarry Smith 6999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL)); 700aa573868SBarry Smith if (!silent) { 701792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl)); 7029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl)); 703ffbd1cfbSBarry Smith } 7047d812c46SBarry Smith } else { 705792fecdfSBarry Smith PetscCallSAWs(SAWs_Initialize, ()); 706aa573868SBarry Smith } 7079566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister("@TechReport{ saws,\n" 7080af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 7090af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 7100af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 7119371c9d4SSatish Balay " Year = 2013\n}\n", 7129371c9d4SSatish Balay NULL)); 71311525c0dSBarry Smith } 7143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71511525c0dSBarry Smith } 71611525c0dSBarry Smith #endif 71711525c0dSBarry Smith 7184dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 719d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) 720d71ae5a4SJacob Faibussowitsch { 7214dfee713SSatish Balay PetscFunctionBegin; 7224dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 7234dfee713SSatish Balay /* see MPI.py for details on this bug */ 7244dfee713SSatish Balay (void)setenv("HWLOC_COMPONENTS", "-x86", 1); 7254dfee713SSatish Balay #endif 7263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7274dfee713SSatish Balay } 7284dfee713SSatish Balay 729a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS) 730a56f64adSBarry Smith #include <adios.h> 73122580e64SBarry Smith #include <adios_read.h> 7327b56e58cSSatish Balay int64_t Petsc_adios_group; 733a56f64adSBarry Smith #endif 734a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP) 735cd1b99a6SBarry Smith #include <omp.h> 736f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 737cd1b99a6SBarry Smith #endif 738a56f64adSBarry Smith 739a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h> 740a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA) 7410e6b6b59SJacob Faibussowitsch #include <petscdevice_cuda.h> 742a4af0ceeSJacob Faibussowitsch // REMOVE ME 743a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL; 744a4af0ceeSJacob Faibussowitsch #endif 745a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP) 7460e6b6b59SJacob Faibussowitsch #include <petscdevice_hip.h> 747a4af0ceeSJacob Faibussowitsch // REMOVE ME 748a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL; 749a4af0ceeSJacob Faibussowitsch #endif 750a4af0ceeSJacob Faibussowitsch 75127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H) 752bc8a24ecSBarry Smith #include <dlfcn.h> 753bc8a24ecSBarry Smith #endif 754a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void); 755a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 7563274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscViennaCLInit(void); 757a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE; 758a4af0ceeSJacob Faibussowitsch #endif 759bc8a24ecSBarry Smith 760660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE; 761660278c0SBarry Smith 76285649d77SJunchao Zhang /* 76385649d77SJunchao Zhang PetscInitialize_Common - shared code between C and Fortran initialization 76485649d77SJunchao Zhang 76585649d77SJunchao Zhang prog: program name 76602101c96SSatish Balay file: optional PETSc database file name. Might be in Fortran string format when 'ftn' is true 76785649d77SJunchao Zhang help: program help message 768da81f932SPierre Jolivet ftn: is it called from Fortran initialization (petscinitializef_)? 76985649d77SJunchao Zhang readarguments,len: used when fortran is true 77085649d77SJunchao Zhang */ 771d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len) 772d71ae5a4SJacob Faibussowitsch { 77385649d77SJunchao Zhang PetscMPIInt size; 77485649d77SJunchao Zhang PetscBool flg = PETSC_TRUE; 77585649d77SJunchao Zhang char hostname[256]; 776046ed546SBarry Smith PetscBool blas_view_flag = PETSC_FALSE; 77785649d77SJunchao Zhang 77827104ee2SJacob Faibussowitsch PetscFunctionBegin; 7793ba16761SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS); 78039a651e2SJacob Faibussowitsch /* these must be initialized in a routine, not as a constant declaration */ 78139a651e2SJacob Faibussowitsch PETSC_STDOUT = stdout; 78239a651e2SJacob Faibussowitsch PETSC_STDERR = stderr; 78339a651e2SJacob Faibussowitsch 7849566063dSJacob Faibussowitsch /* PetscCall can be used from now */ 78539a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_TRUE; 78639a651e2SJacob Faibussowitsch 78785649d77SJunchao Zhang /* 78885649d77SJunchao Zhang The checking over compatible runtime libraries is complicated by the MPI ABI initiative 78985649d77SJunchao Zhang https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 79085649d77SJunchao Zhang MPICH v3.1 (Released February 2014) 79185649d77SJunchao Zhang IBM MPI v2.1 (December 2014) 79285649d77SJunchao Zhang Intel MPI Library v5.0 (2014) 79385649d77SJunchao Zhang Cray MPT v7.0.0 (June 2014) 79485649d77SJunchao Zhang As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 79585649d77SJunchao Zhang listed above and since that time are compatible. 79685649d77SJunchao Zhang 79785649d77SJunchao Zhang Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 79885649d77SJunchao Zhang at compile time or runtime. Thus we will need to systematically track the allowed versions 79985649d77SJunchao Zhang and how they are represented in the mpi.h and MPI_Get_library_version() output in order 80085649d77SJunchao Zhang to perform the checking. 80185649d77SJunchao Zhang 80285649d77SJunchao Zhang Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 80385649d77SJunchao Zhang 80485649d77SJunchao Zhang Questions: 80585649d77SJunchao Zhang 80685649d77SJunchao Zhang Should the checks for ABI incompatibility be only on the major version number below? 80785649d77SJunchao Zhang Presumably the output to stderr will be removed before a release. 80885649d77SJunchao Zhang */ 80985649d77SJunchao Zhang 81085649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 81185649d77SJunchao Zhang { 81285649d77SJunchao Zhang char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 81385649d77SJunchao Zhang PetscMPIInt mpilibraryversionlength; 81439a651e2SJacob Faibussowitsch 8159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength)); 81685649d77SJunchao Zhang /* check for MPICH versions before MPI ABI initiative */ 81785649d77SJunchao Zhang #if defined(MPICH_VERSION) 81885649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000 81985649d77SJunchao Zhang { 82085649d77SJunchao Zhang char *ver, *lf; 82185649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 82239a651e2SJacob Faibussowitsch 8239566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver)); 82439a651e2SJacob Faibussowitsch if (ver) { 8259566063dSJacob Faibussowitsch PetscCall(PetscStrchr(ver, '\n', &lf)); 82639a651e2SJacob Faibussowitsch if (lf) { 82785649d77SJunchao Zhang *lf = 0; 8289566063dSJacob Faibussowitsch PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg)); 82985649d77SJunchao Zhang } 83085649d77SJunchao Zhang } 83185649d77SJunchao Zhang if (!flg) { 832d5b396e8SSatish Balay PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION)); 83385649d77SJunchao Zhang flg = PETSC_TRUE; 83485649d77SJunchao Zhang } 83585649d77SJunchao Zhang } 83685649d77SJunchao Zhang #endif 83785649d77SJunchao 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?) */ 838100ffedbSJunchao Zhang #elif defined(PETSC_HAVE_OPENMPI) 83985649d77SJunchao Zhang { 84085649d77SJunchao Zhang char *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf; 84185649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 84285649d77SJunchao Zhang #define PSTRSZ 2 84385649d77SJunchao Zhang char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"}; 84485649d77SJunchao Zhang char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "}; 84585649d77SJunchao Zhang int i; 84685649d77SJunchao Zhang for (i = 0; i < PSTRSZ; i++) { 8479566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver)); 84839a651e2SJacob Faibussowitsch if (ver) { 849100ffedbSJunchao Zhang PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], PETSC_PKG_OPENMPI_VERSION_MAJOR, PETSC_PKG_OPENMPI_VERSION_MINOR)); 8509566063dSJacob Faibussowitsch PetscCall(PetscStrstr(ver, bs, &bsf)); 85139a651e2SJacob Faibussowitsch if (bsf) flg = PETSC_TRUE; 85285649d77SJunchao Zhang break; 85385649d77SJunchao Zhang } 85485649d77SJunchao Zhang } 85585649d77SJunchao Zhang if (!flg) { 856100ffedbSJunchao 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)); 85785649d77SJunchao Zhang flg = PETSC_TRUE; 85885649d77SJunchao Zhang } 85985649d77SJunchao Zhang } 86085649d77SJunchao Zhang #endif 86185649d77SJunchao Zhang } 86285649d77SJunchao Zhang #endif 86385649d77SJunchao Zhang 864e8953f01SSatish Balay #if defined(PETSC_HAVE_DLADDR) && !(defined(__cray__) && defined(__clang__)) 86585649d77SJunchao 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 */ 86639a651e2SJacob 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"); 86785649d77SJunchao Zhang #endif 86885649d77SJunchao Zhang 86985649d77SJunchao Zhang /* on Windows - set printf to default to printing 2 digit exponents */ 87085649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT) 87185649d77SJunchao Zhang _set_output_format(_TWO_DIGIT_EXPONENT); 87285649d77SJunchao Zhang #endif 87385649d77SJunchao Zhang 8749566063dSJacob Faibussowitsch PetscCall(PetscOptionsCreateDefault()); 87585649d77SJunchao Zhang 87685649d77SJunchao Zhang PetscFinalizeCalled = PETSC_FALSE; 87785649d77SJunchao Zhang 8789566063dSJacob Faibussowitsch PetscCall(PetscSetProgramName(prog)); 8799566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen)); 8809566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout)); 8819566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr)); 8829566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscCommSpinLock)); 88385649d77SJunchao Zhang 88485649d77SJunchao Zhang if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 8859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN)); 88685649d77SJunchao Zhang 88785649d77SJunchao Zhang if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) { 8889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS)); 8899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE)); 89085649d77SJunchao Zhang } 89185649d77SJunchao Zhang 89285649d77SJunchao Zhang /* Done after init due to a bug in MPICH-GM? */ 8939566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 89485649d77SJunchao Zhang 8959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank)); 8969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize)); 89785649d77SJunchao Zhang 89885649d77SJunchao Zhang MPIU_BOOL = MPI_INT; 89985649d77SJunchao Zhang MPIU_ENUM = MPI_INT; 90085649d77SJunchao Zhang MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64; 90185649d77SJunchao Zhang if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 90285649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 90385649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG) 90485649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 90585649d77SJunchao Zhang #endif 90639a651e2SJacob Faibussowitsch else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t"); 90785649d77SJunchao Zhang 90885649d77SJunchao Zhang /* 90985649d77SJunchao Zhang Initialized the global complex variable; this is because with 91085649d77SJunchao Zhang shared libraries the constructors for global variables 91185649d77SJunchao Zhang are not called; at least on IRIX. 91285649d77SJunchao Zhang */ 91385649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX) 91485649d77SJunchao Zhang { 91585649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 91685649d77SJunchao Zhang PetscComplex ic(0.0, 1.0); 91785649d77SJunchao Zhang PETSC_i = ic; 91885649d77SJunchao Zhang #else 91985649d77SJunchao Zhang PETSC_i = _Complex_I; 92085649d77SJunchao Zhang #endif 92185649d77SJunchao Zhang } 92285649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */ 92385649d77SJunchao Zhang 92485649d77SJunchao Zhang /* 92585649d77SJunchao Zhang Create the PETSc MPI reduction operator that sums of the first 92685649d77SJunchao Zhang half of the entries and maxes the second half. 92785649d77SJunchao Zhang */ 9289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP)); 92985649d77SJunchao Zhang 930613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 9319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128)); 9329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128)); 93374c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 9349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128)); 9359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128)); 93685649d77SJunchao Zhang #endif 93774c01ffaSSatish Balay #endif 9389e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) 9399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16)); 9409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FP16)); 94185649d77SJunchao Zhang #endif 94285649d77SJunchao Zhang 94385649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 9449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM)); 945613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX)); 946613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN)); 9479e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16) 9489e517322SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128)); 94985649d77SJunchao Zhang #endif 95085649d77SJunchao Zhang 9519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR)); 95262e5d2d2SJDBetteridge PetscCallMPI(MPI_Op_create(PetscGarbageKeySortedIntersect, 1, &Petsc_Garbage_SetIntersectOp)); 9539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR)); 95485649d77SJunchao Zhang 95585649d77SJunchao Zhang /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */ 95685649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI) 95785649d77SJunchao Zhang { 95885649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 95993d501b3SJacob Faibussowitsch MPI_Aint blockOffsets[2] = {offsetof(struct petsc_mpiu_real_int, v), offsetof(struct petsc_mpiu_real_int, i)}; 96085649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_REAL, MPIU_INT}, tmpStruct; 96185649d77SJunchao Zhang 9629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 96393d501b3SJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_real_int), &MPIU_REAL_INT)); 9649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT)); 96685649d77SJunchao Zhang } 96785649d77SJunchao Zhang { 96885649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 96993d501b3SJacob Faibussowitsch MPI_Aint blockOffsets[2] = {offsetof(struct petsc_mpiu_scalar_int, v), offsetof(struct petsc_mpiu_scalar_int, i)}; 97085649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_SCALAR, MPIU_INT}, tmpStruct; 97185649d77SJunchao Zhang 9729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 97393d501b3SJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_scalar_int), &MPIU_SCALAR_INT)); 9749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT)); 97685649d77SJunchao Zhang } 97785649d77SJunchao Zhang #endif 97885649d77SJunchao Zhang 97985649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES) 9809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT)); 9819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2INT)); 982*6497c311SBarry Smith 983*6497c311SBarry Smith #if !defined(PETSC_HAVE_MPIUNI) 984*6497c311SBarry Smith { 985*6497c311SBarry Smith int blockSizes[] = {1, 1}; 986*6497c311SBarry Smith MPI_Aint blockOffsets[] = {offsetof(struct petsc_mpiu_int_mpiint, a), offsetof(struct petsc_mpiu_int_mpiint, b)}; 987*6497c311SBarry Smith MPI_Datatype blockTypes[] = {MPIU_INT, MPI_INT}, tmpStruct; 988*6497c311SBarry Smith 989*6497c311SBarry Smith PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 990*6497c311SBarry Smith PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct petsc_mpiu_int_mpiint), &MPIU_INT_MPIINT)); 991*6497c311SBarry Smith PetscCallMPI(MPI_Type_free(&tmpStruct)); 992*6497c311SBarry Smith PetscCallMPI(MPI_Type_commit(&MPIU_INT_MPIINT)); 993*6497c311SBarry Smith } 994*6497c311SBarry Smith #endif 99585649d77SJunchao Zhang #endif 9969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT)); 9979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPI_4INT)); 9989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT)); 9999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_4INT)); 100085649d77SJunchao Zhang 100185649d77SJunchao Zhang /* 100285649d77SJunchao Zhang Attributes to be set on PETSc communicators 100385649d77SJunchao Zhang */ 10048434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_DeleteFn, &Petsc_Counter_keyval, (void *)0)); 10058434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_DeleteFn, &Petsc_InnerComm_keyval, (void *)0)); 10068434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_DeleteFn, &Petsc_OuterComm_keyval, (void *)0)); 10078434afd1SBarry Smith PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_DeleteFn, &Petsc_ShmComm_keyval, (void *)0)); 100862e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_CreationIdx_keyval, (void *)0)); 100962e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Garbage_HMap_keyval, (void *)0)); 101085649d77SJunchao Zhang 1011fbf9dbe5SBarry Smith #if defined(PETSC_USE_FORTRAN_BINDINGS) 10129566063dSJacob Faibussowitsch if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len)); 101385649d77SJunchao Zhang else 101485649d77SJunchao Zhang #endif 10159566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file)); 101685649d77SJunchao Zhang 101785649d77SJunchao Zhang /* call a second time so it can look in the options database */ 10189566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 101985649d77SJunchao Zhang 102085649d77SJunchao Zhang /* 102185649d77SJunchao Zhang Check system options and print help 102285649d77SJunchao Zhang */ 10239566063dSJacob Faibussowitsch PetscCall(PetscOptionsCheckInitial_Private(help)); 102485649d77SJunchao Zhang 1025a4af0ceeSJacob Faibussowitsch /* 10260e6b6b59SJacob Faibussowitsch Creates the logging data structures; this is enabled even if logging is not turned on 10270e6b6b59SJacob Faibussowitsch This is the last thing we do before returning to the user code to prevent having the 10280e6b6b59SJacob Faibussowitsch logging numbers contaminated by any startup time associated with MPI 10290e6b6b59SJacob Faibussowitsch */ 10300e6b6b59SJacob Faibussowitsch PetscCall(PetscLogInitialize()); 10310e6b6b59SJacob Faibussowitsch 10320e6b6b59SJacob Faibussowitsch /* 1033a4af0ceeSJacob Faibussowitsch Initialize PetscDevice and PetscDeviceContext 1034a4af0ceeSJacob Faibussowitsch 1035a4af0ceeSJacob Faibussowitsch Note to any future devs thinking of moving this, proper initialization requires: 1036a4af0ceeSJacob Faibussowitsch 1. MPI initialized 1037a4af0ceeSJacob Faibussowitsch 2. Options DB initialized 10380e6b6b59SJacob Faibussowitsch 3. Petsc error handling initialized, specifically signal handlers. This expects to set up 10390e6b6b59SJacob Faibussowitsch its own SIGSEV handler via the push/pop interface. 10400e6b6b59SJacob Faibussowitsch 4. Logging initialized 1041a4af0ceeSJacob Faibussowitsch */ 10429566063dSJacob Faibussowitsch PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD)); 1043a4af0ceeSJacob Faibussowitsch 1044a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 1045a4af0ceeSJacob Faibussowitsch flg = PETSC_FALSE; 1046d75802c7SJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg)); 10479566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL)); 1048a4af0ceeSJacob Faibussowitsch PetscViennaCLSynchronize = flg; 10499566063dSJacob Faibussowitsch PetscCall(PetscViennaCLInit()); 1050a4af0ceeSJacob Faibussowitsch #endif 1051a4af0ceeSJacob Faibussowitsch 10529566063dSJacob Faibussowitsch PetscCall(PetscCitationsInitialize()); 105385649d77SJunchao Zhang 105485649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS) 10559566063dSJacob Faibussowitsch PetscCall(PetscInitializeSAWs(ftn ? NULL : help)); 105627104ee2SJacob Faibussowitsch flg = PETSC_FALSE; 10579566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg)); 10589566063dSJacob Faibussowitsch if (flg) PetscCall(PetscStackViewSAWs()); 105985649d77SJunchao Zhang #endif 106085649d77SJunchao Zhang 106185649d77SJunchao Zhang /* 106285649d77SJunchao Zhang Load the dynamic libraries (on machines that support them), this registers all 106385649d77SJunchao Zhang the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 106485649d77SJunchao Zhang */ 10659566063dSJacob Faibussowitsch PetscCall(PetscInitialize_DynamicLibraries()); 106685649d77SJunchao Zhang 10679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 10689566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size)); 106996dcfd6cSBarry Smith PetscCall(PetscGetHostName(hostname, sizeof(hostname))); 10709566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname)); 107185649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP) 107285649d77SJunchao Zhang { 107385649d77SJunchao Zhang PetscBool omp_view_flag; 107485649d77SJunchao Zhang char *threads = getenv("OMP_NUM_THREADS"); 107585649d77SJunchao Zhang 107685649d77SJunchao Zhang if (threads) { 10779566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads)); 107885649d77SJunchao Zhang (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads); 107985649d77SJunchao Zhang } else { 10802f840973SStefano Zampini PetscNumOMPThreads = (PetscInt)omp_get_max_threads(); 10819566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads)); 108285649d77SJunchao Zhang } 1083d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys"); 10849566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg)); 10859566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag)); 1086d0609cedSBarry Smith PetscOptionsEnd(); 108785649d77SJunchao Zhang if (flg) { 10881a0ee928SPierre Jolivet PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads)); 108985649d77SJunchao Zhang omp_set_num_threads((int)PetscNumOMPThreads); 109085649d77SJunchao Zhang } 109148a46eb9SPierre Jolivet if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads)); 109285649d77SJunchao Zhang } 109385649d77SJunchao Zhang #endif 109485649d77SJunchao Zhang 1095046ed546SBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "BLAS options", "Sys"); 1096046ed546SBarry Smith PetscCall(PetscOptionsName("-blas_view", "Display number of threads to use for BLAS operations", NULL, &blas_view_flag)); 1097046ed546SBarry Smith #if defined(PETSC_HAVE_BLI_THREAD_SET_NUM_THREADS) || defined(PETSC_HAVE_MKL_SET_NUM_THREADS) || defined(PETSC_HAVE_OPENBLAS_SET_NUM_THREADS) 1098046ed546SBarry Smith { 1099046ed546SBarry Smith char *threads = NULL; 1100046ed546SBarry Smith 1101046ed546SBarry Smith /* determine any default number of threads requested in the environment; TODO: Apple libraries? */ 1102046ed546SBarry Smith #if defined(PETSC_HAVE_BLI_THREAD_SET_NUM_THREADS) 1103046ed546SBarry Smith threads = getenv("BLIS_NUM_THREADS"); 1104046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of BLIS threads %s given by BLIS_NUM_THREADS\n", threads)); 1105046ed546SBarry Smith if (!threads) { 1106046ed546SBarry Smith threads = getenv("OMP_NUM_THREADS"); 1107046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of BLIS threads %s given by OMP_NUM_THREADS\n", threads)); 1108046ed546SBarry Smith } 1109046ed546SBarry Smith #elif defined(PETSC_HAVE_MKL_SET_NUM_THREADS) 1110046ed546SBarry Smith threads = getenv("MKL_NUM_THREADS"); 1111046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of MKL threads %s given by MKL_NUM_THREADS\n", threads)); 111297873304SJunchao Zhang if (!threads) { 111397873304SJunchao Zhang threads = getenv("OMP_NUM_THREADS"); 111497873304SJunchao Zhang if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of MKL threads %s given by OMP_NUM_THREADS\n", threads)); 111597873304SJunchao Zhang } 1116046ed546SBarry Smith #elif defined(PETSC_HAVE_OPENBLAS_SET_NUM_THREADS) 1117046ed546SBarry Smith threads = getenv("OPENBLAS_NUM_THREADS"); 1118046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of OpenBLAS threads %s given by OPENBLAS_NUM_THREADS\n", threads)); 1119046ed546SBarry Smith if (!threads) { 1120046ed546SBarry Smith threads = getenv("OMP_NUM_THREADS"); 1121046ed546SBarry Smith if (threads) PetscCall(PetscInfo(NULL, "BLAS: Environment number of OpenBLAS threads %s given by OMP_NUM_THREADS\n", threads)); 1122046ed546SBarry Smith } 1123046ed546SBarry Smith #endif 1124046ed546SBarry Smith if (threads) (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumBLASThreads); 1125046ed546SBarry Smith PetscCall(PetscOptionsInt("-blas_num_threads", "Number of threads to use for BLAS operations", "None", PetscNumBLASThreads, &PetscNumBLASThreads, &flg)); 112697873304SJunchao Zhang if (flg) PetscCall(PetscInfo(NULL, "BLAS: Command line number of BLAS thread %" PetscInt_FMT "given by -blas_num_threads\n", PetscNumBLASThreads)); 112797873304SJunchao Zhang if (flg || threads) { 1128046ed546SBarry Smith PetscCall(PetscBLASSetNumThreads(PetscNumBLASThreads)); 1129046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: number of threads %" PetscInt_FMT "\n", PetscNumBLASThreads)); 1130046ed546SBarry Smith } 113197873304SJunchao Zhang } 1132046ed546SBarry Smith #elif defined(PETSC_HAVE_APPLE_ACCELERATE) 1133046ed546SBarry Smith PetscCall(PetscInfo(NULL, "BLAS: Apple Accelerate library, thread support with no user control\n")); 1134046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: Apple Accelerate library, thread support with no user control\n")); 1135046ed546SBarry Smith #else 1136046ed546SBarry Smith if (blas_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BLAS: no thread support\n")); 1137046ed546SBarry Smith #endif 1138046ed546SBarry Smith PetscOptionsEnd(); 1139046ed546SBarry Smith 114085649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 114185649d77SJunchao Zhang /* 114285649d77SJunchao Zhang Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 114385649d77SJunchao Zhang 114485649d77SJunchao Zhang Currently not used because it is not supported by MPICH. 114585649d77SJunchao Zhang */ 11469566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL)); 114785649d77SJunchao Zhang #endif 114885649d77SJunchao Zhang 114985649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS) 11509566063dSJacob Faibussowitsch PetscCall(PetscFPTCreate(10000)); 115185649d77SJunchao Zhang #endif 115285649d77SJunchao Zhang 115385649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC) 115485649d77SJunchao Zhang { 115585649d77SJunchao Zhang PetscViewer viewer; 1156648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg)); 115785649d77SJunchao Zhang if (flg) { 11589566063dSJacob Faibussowitsch PetscCall(PetscProcessPlacementView(viewer)); 1159648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 116085649d77SJunchao Zhang } 116185649d77SJunchao Zhang } 116285649d77SJunchao Zhang #endif 116385649d77SJunchao Zhang 116485649d77SJunchao Zhang flg = PETSC_TRUE; 11659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL)); 1166648c30bcSBarry Smith if (!flg) PetscCall(PetscOptionsPushCreateViewerOff(PETSC_TRUE)); 116785649d77SJunchao Zhang 116885649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS) 11693ba16761SJacob Faibussowitsch PetscCallExternal(adios_init_noxml, PETSC_COMM_WORLD); 11703ba16761SJacob Faibussowitsch PetscCallExternal(adios_declare_group, &Petsc_adios_group, "PETSc", "", adios_stat_default); 11713ba16761SJacob Faibussowitsch PetscCallExternal(adios_select_method, Petsc_adios_group, "MPI", "", ""); 11723ba16761SJacob Faibussowitsch PetscCallExternal(adios_read_init_method, ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, ""); 117385649d77SJunchao Zhang #endif 117485649d77SJunchao Zhang 117585649d77SJunchao Zhang #if defined(__VALGRIND_H) 117685649d77SJunchao Zhang PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE; 117785649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE) 1178337bb527SBarry 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")); 117985649d77SJunchao Zhang #endif 118085649d77SJunchao Zhang #endif 118185649d77SJunchao Zhang /* 118285649d77SJunchao Zhang Set flag that we are completely initialized 118385649d77SJunchao Zhang */ 118485649d77SJunchao Zhang PetscInitializeCalled = PETSC_TRUE; 118585649d77SJunchao Zhang 11869566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg)); 11879566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPythonInitialize(NULL, NULL)); 1188f1f2ae84SBarry Smith 1189f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 11909f0612e4SBarry Smith if (flg) PetscCall(PetscInfo(NULL, "Running MPI Linear Solver Server\n")); 1191f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin()); 1192f1f2ae84SBarry 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"); 11933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119485649d77SJunchao Zhang } 119585649d77SJunchao Zhang 119610450e9eSJacob Faibussowitsch // "Unknown section 'Environmental Variables'" 119710450e9eSJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-section-header-unknown 1198e5c89e4eSSatish Balay /*@C 1199e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 1200811af0c4SBarry Smith `PetscInitialize()` calls MPI_Init() if that has yet to be called, 1201e5c89e4eSSatish Balay so this routine should always be called near the beginning of 1202e5c89e4eSSatish Balay your program -- usually the very first line! 1203e5c89e4eSSatish Balay 1204811af0c4SBarry Smith Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set 1205e5c89e4eSSatish Balay 1206e5c89e4eSSatish Balay Input Parameters: 1207e5c89e4eSSatish Balay + argc - count of number of command line arguments 1208e5c89e4eSSatish Balay . args - the command line arguments 1209be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format. 1210be10d61cSLisandro Dalcin Use NULL or empty string to not check for code specific file. 1211be10d61cSLisandro Dalcin Also checks ~/.petscrc, .petscrc and petscrc. 1212c5b5d8d5SVaclav Hapla Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files. 12130298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message 1214e5c89e4eSSatish Balay 1215811af0c4SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that 1216811af0c4SBarry Smith communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a 1217811af0c4SBarry Smith four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not, 1218811af0c4SBarry Smith then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even 121905827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 1220e5c89e4eSSatish Balay 1221e5c89e4eSSatish Balay Options Database Keys: 12227ca660e7SBarry Smith + -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message 12237ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 1224e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 12258a690491SBarry Smith . -on_error_emacs <machinename> - causes emacsclient to jump to error file 1226811af0c4SBarry Smith . -on_error_abort - calls `abort()` when error detected (no traceback) 1227811af0c4SBarry Smith . -on_error_mpiabort - calls `MPI_abort()` when error detected 12281af3601dSBarry Smith . -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr 12298a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 1230bf4d2887SBarry Smith . -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger 1231e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds) - Pauses debugger 1232e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 1233e5c89e4eSSatish Balay process and wait (-debugger_pause) seconds for attachment 1234aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 123592f119d6SBarry Smith . -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds, ignored in optimized build. May want to set in PETSC_OPTIONS environmental variable 1236811af0c4SBarry Smith . -malloc_view - show a list of all allocated memory during `PetscFinalize()` 123792f119d6SBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view 1238608c71bfSMatthew G. Knepley . -malloc_requested_size - malloc logging will record the requested size rather than size after alignment 123992f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 1240e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 1241e5c89e4eSSatish Balay . -shared_tmp - indicates /tmp directory is shared by all processors 1242e5c89e4eSSatish Balay . -not_shared_tmp - each processor has own /tmp 1243e5c89e4eSSatish Balay . -tmp - alternative name of /tmp directory 1244e5c89e4eSSatish Balay . -get_total_flops - returns total flops done by all processors 12450841954dSBarry Smith - -memory_view - Print memory usage at end of run 1246e5c89e4eSSatish Balay 1247c5b5d8d5SVaclav Hapla Options Database Keys for Option Database: 1248c5b5d8d5SVaclav Hapla + -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc 1249c5b5d8d5SVaclav Hapla . -options_monitor - monitor all set options to standard output for the whole program run 1250811af0c4SBarry Smith - -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()` 1251c5b5d8d5SVaclav Hapla 1252c5b5d8d5SVaclav Hapla Options -options_monitor_{all,cancel} are 1253c5b5d8d5SVaclav Hapla position-independent and apply to all options set since the PETSc start. 1254c5b5d8d5SVaclav Hapla They can be used also in option files. 1255c5b5d8d5SVaclav Hapla 1256811af0c4SBarry Smith See `PetscOptionsMonitorSet()` to do monitoring programmatically. 1257c5b5d8d5SVaclav Hapla 1258e5c89e4eSSatish Balay Options Database Keys for Profiling: 1259a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 1260811af0c4SBarry Smith + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`. 1261217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 1262217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 1263495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 1264811af0c4SBarry Smith hangs without running in the debugger). See `PetscLogTraceBegin()`. 1265ad2e3d55SToby Isaac . -log_view [:filename:format][,[:filename:format]...] - Prints summary of flop and timing information to screen or file, see `PetscLogView()` (up to 4 viewers) 1266811af0c4SBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`. 1267811af0c4SBarry Smith . -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView(). 1268f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 1269b665b14eSToby Isaac . -log [filename] - Logs profiling information in a dump file, see `PetscLogDump()`. 1270b665b14eSToby Isaac . -log_all [filename] - Same as `-log`. 1271c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 12722611ad71SToby Isaac . -log_perfstubs - Starts a log handler with the perfstubs interface (which is used by TAU) 127361cc7448SToby Isaac . -log_nvtx - Starts an nvtx log handler for use with Nsight 1274811af0c4SBarry Smith . -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off 1275c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 1276495fc317SBarry Smith 1277ffbd1cfbSBarry Smith Options Database Keys for SAWs: 1278ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 1279ffbd1cfbSBarry 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 1280ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 1281ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 1282ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 1283ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 1284ffbd1cfbSBarry Smith 1285e5c89e4eSSatish Balay Environmental Variables: 1286811af0c4SBarry Smith + `PETSC_TMP` - alternative tmp directory 1287811af0c4SBarry Smith . `PETSC_SHARED_TMP` - tmp is shared by all processes 1288811af0c4SBarry Smith . `PETSC_NOT_SHARED_TMP` - each process has its own private tmp 1289811af0c4SBarry Smith . `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs 1290811af0c4SBarry Smith . `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document 1291811af0c4SBarry Smith . `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer 1292811af0c4SBarry Smith - `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to 1293e5c89e4eSSatish Balay 1294e5c89e4eSSatish Balay Level: beginner 1295e5c89e4eSSatish Balay 1296811af0c4SBarry Smith Note: 1297811af0c4SBarry Smith If for some reason you must call `MPI_Init()` separately, call 1298811af0c4SBarry Smith it before `PetscInitialize()`. 1299e5c89e4eSSatish Balay 1300811af0c4SBarry Smith Fortran Notes: 1301811af0c4SBarry Smith In Fortran this routine can be called with 1302811af0c4SBarry Smith .vb 1303811af0c4SBarry Smith call PetscInitialize(ierr) 1304811af0c4SBarry Smith call PetscInitialize(file,ierr) or 1305811af0c4SBarry Smith call PetscInitialize(file,help,ierr) 1306811af0c4SBarry Smith .ve 1307e5c89e4eSSatish Balay 1308811af0c4SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after 1309811af0c4SBarry Smith calling `PetscInitialize()`. 1310e5c89e4eSSatish Balay 13115dedd997SBarry Smith Options Database Key for Developers: 13125dedd997SBarry Smith . -checkfunctionlist - automatically checks that function lists associated with objects are correctly cleaned up. Produces messages of the form: 13135dedd997SBarry Smith "function name: MatInodeGetInodeSizes_C" if they are not cleaned up. This flag is always set for the test harness (in framework.py) 13145dedd997SBarry Smith 1315db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()` 1316e5c89e4eSSatish Balay @*/ 1317d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[]) 1318d71ae5a4SJacob Faibussowitsch { 131985649d77SJunchao Zhang PetscMPIInt flag; 132021ef0414SBarry Smith const char *prog = "Unknown Name", *mpienv; 1321e5c89e4eSSatish Balay 132227104ee2SJacob Faibussowitsch PetscFunctionBegin; 13233ba16761SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS); 13249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Initialized(&flag)); 1325e5c89e4eSSatish Balay if (!flag) { 132639a651e2SJacob 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"); 13279566063dSJacob Faibussowitsch PetscCall(PetscPreMPIInit_Private()); 13285e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 13295e765c61SJed Brown { 1330835d5ba2SJunchao Zhang PetscMPIInt provided; 1331835d5ba2SJunchao Zhang PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE ? MPI_THREAD_FUNNELED : PETSC_MPI_THREAD_REQUIRED, &provided)); 1332835d5ba2SJunchao 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"); 1333835d5ba2SJunchao Zhang if (PETSC_MPI_THREAD_REQUIRED == PETSC_DECIDE) PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED; // assign it a valid value after check-up 13345e765c61SJed Brown } 13355e765c61SJed Brown #else 13369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init(argc, args)); 13375e765c61SJed Brown #endif 133821ef0414SBarry Smith if (PetscDefined(HAVE_MPIUNI)) { 133921ef0414SBarry Smith mpienv = getenv("PMI_SIZE"); 134021ef0414SBarry Smith if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE"); 134121ef0414SBarry Smith if (mpienv) { 134221ef0414SBarry Smith PetscInt isize; 134321ef0414SBarry Smith PetscCall(PetscOptionsStringToInt(mpienv, &isize)); 134421ef0414SBarry Smith if (isize != 1) printf("You are using an MPI-uni (sequential) install of PETSc but trying to launch parallel jobs; you need full MPI version of PETSc\n"); 134521ef0414SBarry Smith PetscCheck(isize == 1, MPI_COMM_SELF, PETSC_ERR_MPI, "You are using an MPI-uni (sequential) install of PETSc but trying to launch parallel jobs; you need full MPI version of PETSc"); 134621ef0414SBarry Smith } 134721ef0414SBarry Smith } 1348e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 1349e5c89e4eSSatish Balay } 13504dfee713SSatish Balay 135185649d77SJunchao Zhang if (argc && *argc) prog = **args; 1352e5c89e4eSSatish Balay if (argc && args) { 1353e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 1354e5c89e4eSSatish Balay PetscGlobalArgs = *args; 1355e5c89e4eSSatish Balay } 1356811af0c4SBarry Smith PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0)); 13573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1358e5c89e4eSSatish Balay } 1359e5c89e4eSSatish Balay 136095c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1361ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1362ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 136395c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 1364e5c89e4eSSatish Balay 1365008a6e76SBarry Smith /* 1366008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1367008a6e76SBarry Smith */ 1368d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeMPIResources(void) 1369d71ae5a4SJacob Faibussowitsch { 1370008a6e76SBarry Smith PetscFunctionBegin; 1371613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 13729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128)); 137374c01ffaSSatish Balay #if defined(PETSC_HAVE_COMPLEX) 13749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128)); 1375008a6e76SBarry Smith #endif 137674c01ffaSSatish Balay #endif 13779e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) 13789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FP16)); 1379008a6e76SBarry Smith #endif 1380008a6e76SBarry Smith 1381de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 13829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_SUM)); 1383613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MAX)); 1384613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MIN)); 13859e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16) 13869e517322SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128)); 1387008a6e76SBarry Smith #endif 1388008a6e76SBarry Smith 13899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR)); 13909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT)); 13919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT)); 139240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 13939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2INT)); 1394*6497c311SBarry Smith PetscCallMPI(MPI_Type_free(&MPIU_INT_MPIINT)); 1395008a6e76SBarry Smith #endif 13969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPI_4INT)); 13979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_4INT)); 13989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP)); 139962e5d2d2SJDBetteridge PetscCallMPI(MPI_Op_free(&Petsc_Garbage_SetIntersectOp)); 14003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1401008a6e76SBarry Smith } 1402008a6e76SBarry Smith 1403a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void); 14049f0612e4SBarry Smith PETSC_EXTERN PetscErrorCode PetscFreeAlign(void *, int, const char[], const char[]); 1405a4af0ceeSJacob Faibussowitsch 14069f0612e4SBarry Smith /*@ 1407e5c89e4eSSatish Balay PetscFinalize - Checks for options to be called at the conclusion 1408811af0c4SBarry Smith of the program. `MPI_Finalize()` is called only if the user had not 1409811af0c4SBarry Smith called `MPI_Init()` before calling `PetscInitialize()`. 1410e5c89e4eSSatish Balay 1411811af0c4SBarry Smith Collective on `PETSC_COMM_WORLD` 1412e5c89e4eSSatish Balay 1413e5c89e4eSSatish Balay Options Database Keys: 1414811af0c4SBarry Smith + -options_view - Calls `PetscOptionsView()` 1415e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 14167eb1d149SBarry 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 1417e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 1418811af0c4SBarry Smith . -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed 14194bd3d7f8SBarry Smith . -memory_view - Prints total memory usage 14204bd3d7f8SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and in what functions 1421e5c89e4eSSatish Balay 1422e5c89e4eSSatish Balay Level: beginner 1423e5c89e4eSSatish Balay 1424e5c89e4eSSatish Balay Note: 1425811af0c4SBarry Smith See `PetscInitialize()` for other runtime options. 1426e5c89e4eSSatish Balay 1427db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()` 1428e5c89e4eSSatish Balay @*/ 1429d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFinalize(void) 1430d71ae5a4SJacob Faibussowitsch { 14314bb5149bSJed Brown PetscMPIInt rank; 1432a8d2bbe5SBarry Smith PetscInt nopt; 14332bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE; 1434dff31646SBarry Smith PetscBool flg; 143510463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 1436e5c89e4eSSatish Balay 14373cbbc5ffSBarry Smith PetscFunctionBegin; 143839a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()"); 14399566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PetscFinalize() called\n")); 1440b022a5c1SBarry Smith 1441f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1442f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd()); 1443f1f2ae84SBarry Smith 14449f0612e4SBarry Smith PetscCall(PetscFreeAlign(PetscGlobalArgsFortran, 0, NULL, NULL)); 14459f0612e4SBarry Smith PetscGlobalArgc = 0; 14469f0612e4SBarry Smith PetscGlobalArgs = NULL; 14479f0612e4SBarry Smith 144862e5d2d2SJDBetteridge /* Clean up Garbage automatically on COMM_SELF and COMM_WORLD at finalize */ 144962e5d2d2SJDBetteridge { 145062e5d2d2SJDBetteridge union 145162e5d2d2SJDBetteridge { 145262e5d2d2SJDBetteridge MPI_Comm comm; 145362e5d2d2SJDBetteridge void *ptr; 145462e5d2d2SJDBetteridge } ucomm; 145562e5d2d2SJDBetteridge PetscMPIInt flg; 145662e5d2d2SJDBetteridge void *tmp; 145762e5d2d2SJDBetteridge 145862e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 145962e5d2d2SJDBetteridge if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg)); 146062e5d2d2SJDBetteridge if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_SELF)); 146162e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 146262e5d2d2SJDBetteridge if (flg) PetscCallMPI(MPI_Comm_get_attr(ucomm.comm, Petsc_Garbage_HMap_keyval, &tmp, &flg)); 146362e5d2d2SJDBetteridge if (flg) PetscCall(PetscGarbageCleanup(PETSC_COMM_WORLD)); 146462e5d2d2SJDBetteridge } 146562e5d2d2SJDBetteridge 14669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1467a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 14683ba16761SJacob Faibussowitsch PetscCallExternal(adios_read_finalize_method, ADIOS_READ_METHOD_BP_AGGREGATE); 14693ba16761SJacob Faibussowitsch PetscCallExternal(adios_finalize, rank); 1470a56f64adSBarry Smith #endif 14719566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg)); 1472dff31646SBarry Smith if (flg) { 14731f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 14741f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 14751f817a21SBarry Smith 14769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL)); 147748a46eb9SPierre Jolivet if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd)); 14789566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits)); 1479dff31646SBarry Smith cits[0] = 0; 14809566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits)); 14819566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n")); 14829566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 14839566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits)); 14849566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 14859566063dSJacob Faibussowitsch PetscCall(PetscFClose(PETSC_COMM_WORLD, fd)); 14869566063dSJacob Faibussowitsch PetscCall(PetscFree(cits)); 1487dff31646SBarry Smith } 14889566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&PetscCitationsList)); 1489dff31646SBarry Smith 14902d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 14919566063dSJacob Faibussowitsch PetscCall(PetscFPTDestroy()); 14922d53ad75SBarry Smith #endif 14932d53ad75SBarry Smith 1494e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1495dff31646SBarry Smith flg = PETSC_FALSE; 14969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL)); 14971baa6e33SBarry Smith if (flg) PetscCall(PetscOptionsSAWsDestroy()); 1498d5649816SBarry Smith #endif 1499d5649816SBarry Smith 1500681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1501681455b2SBarry Smith flg1 = PETSC_FALSE; 15029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL)); 1503681455b2SBarry Smith if (flg1) { 1504681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 1505a5af8288SJose E. Roman PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -15 Xvfb", "r", NULL)); 1506681455b2SBarry Smith } 1507681455b2SBarry Smith #endif 1508681455b2SBarry Smith 150967584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 15109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL)); 151148a46eb9SPierre Jolivet if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n")); 151267584ceeSBarry Smith #endif 1513e5c89e4eSSatish Balay 15142611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 151590d69ab7SBarry Smith flg1 = PETSC_FALSE; 15169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL)); 1517e5c89e4eSSatish Balay if (flg1) { 1518e5c89e4eSSatish Balay PetscLogDouble flops = 0; 15199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD)); 15209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops)); 1521e5c89e4eSSatish Balay } 15222611ad71SToby Isaac } 1523e5c89e4eSSatish Balay 15242611ad71SToby Isaac if (PetscDefined(USE_LOG) && PetscDefined(HAVE_MPE)) { 1525e5c89e4eSSatish Balay mname[0] = 0; 15269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1)); 15272611ad71SToby Isaac if (flg1) PetscCall(PetscLogMPEDump(mname[0] ? mname : NULL)); 1528e5c89e4eSSatish Balay } 1529a297a907SKarl Rupp 1530c9903f8fSJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 1531c9903f8fSJunchao Zhang // Free petsc/kokkos stuff before the potentially non-null petsc default gpu stream is destroyed by PetscObjectRegisterDestroyAll 1532c9903f8fSJunchao Zhang if (PetscKokkosInitialized) { 1533c9903f8fSJunchao Zhang PetscCall(PetscKokkosFinalize_Private()); 1534c9903f8fSJunchao Zhang PetscKokkosInitialized = PETSC_FALSE; 1535c9903f8fSJunchao Zhang } 1536c9903f8fSJunchao Zhang #endif 1537c9903f8fSJunchao Zhang 15383048253cSJacob Faibussowitsch // Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 15399566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1540dd710f27SBarry Smith 15412611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 1542648c30bcSBarry Smith PetscCall(PetscOptionsPushCreateViewerOff(PETSC_FALSE)); 15439566063dSJacob Faibussowitsch PetscCall(PetscLogViewFromOptions()); 1544648c30bcSBarry Smith PetscCall(PetscOptionsPopCreateViewerOff()); 1545473903fcSJunchao Zhang // It should be turned on with PetscLogGpuTime() and never turned off except in this place 1546473903fcSJunchao Zhang PetscLogGpuTimeFlag = PETSC_FALSE; 154787c3beb6SLisandro Dalcin 15483048253cSJacob Faibussowitsch // Free any objects created by the last block of code. 15499566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1550dd710f27SBarry Smith 1551dd710f27SBarry Smith mname[0] = 0; 15529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1)); 15539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2)); 15549566063dSJacob Faibussowitsch if (flg1 || flg2) PetscCall(PetscLogDump(mname)); 15552611ad71SToby Isaac } 155610463e74SBarry Smith 155790d69ab7SBarry Smith flg1 = PETSC_FALSE; 15589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL)); 15599566063dSJacob Faibussowitsch if (!flg1) PetscCall(PetscPopSignalHandler()); 156090d69ab7SBarry Smith flg1 = PETSC_FALSE; 15619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL)); 15621baa6e33SBarry Smith if (flg1) PetscCall(PetscMPIDump(stdout)); 156390d69ab7SBarry Smith flg1 = PETSC_FALSE; 156490d69ab7SBarry Smith flg2 = PETSC_FALSE; 15658bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 15669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1)); 15679566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 15689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL)); 1569e4c476e2SSatish Balay 157057b1f488SBarry Smith if (flg2) { PetscCall(PetscOptionsView(NULL, PETSC_VIEWER_STDOUT_WORLD)); } 1571e5c89e4eSSatish Balay 1572e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 15739566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1)); 15749566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1)); 1575e5c89e4eSSatish Balay 157633fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 15779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1)); 157859f199a7SJed Brown if (!flg1) flg3 = PETSC_TRUE; 1579e5c89e4eSSatish Balay if (flg3) { 15803de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 158157b1f488SBarry Smith PetscCall(PetscOptionsView(NULL, PETSC_VIEWER_STDOUT_WORLD)); 1582e5c89e4eSSatish Balay } 15839566063dSJacob Faibussowitsch PetscCall(PetscOptionsAllUsed(NULL, &nopt)); 15843de2bfdfSBarry Smith if (nopt) { 15859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n")); 15869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n")); 15873de2bfdfSBarry Smith if (nopt == 1) { 15889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n")); 1589e5c89e4eSSatish Balay } else { 15909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt)); 1591e5c89e4eSSatish Balay } 15923de2bfdfSBarry Smith } else if (flg3 && flg1) { 15939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n")); 1594df12ba86SBarry Smith } 15959566063dSJacob Faibussowitsch PetscCall(PetscOptionsLeft(NULL)); 1596e5c89e4eSSatish Balay } 1597e5c89e4eSSatish Balay 1598e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1599d45a07a7SBarry Smith if (!PetscGlobalRank) { 16009566063dSJacob Faibussowitsch PetscCall(PetscStackSAWsViewOff()); 1601792fecdfSBarry Smith PetscCallSAWs(SAWs_Finalize, ()); 1602d45a07a7SBarry Smith } 1603ec957eceSBarry Smith #endif 1604ec957eceSBarry Smith 160510463e74SBarry Smith /* 1606dbc8283eSBarry Smith List all objects the user may have forgot to free 16072eff7a51SBarry Smith */ 16082611ad71SToby Isaac if (PetscDefined(USE_LOG) && PetscObjectsLog) { 16099566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 1610a64a8e02SBarry Smith if (flg1) { 1611a64a8e02SBarry Smith MPI_Comm local_comm; 16127eb1d149SBarry Smith char string[64]; 1613a64a8e02SBarry Smith 16149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL)); 1615afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16169566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16179566063dSJacob Faibussowitsch PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE)); 16189566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 16200a1571b3SBarry Smith } 162105df10baSBarry Smith } 16224097062eSBarry Smith 1623dbc8283eSBarry Smith PetscObjectsCounts = 0; 1624dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 16259566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscObjects)); 16262eff7a51SBarry Smith 162733f85c2fSBarry Smith /* 162833f85c2fSBarry Smith Destroy any packages that registered a finalize 162933f85c2fSBarry Smith */ 16309566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalizeAll()); 163133f85c2fSBarry Smith 16329566063dSJacob Faibussowitsch PetscCall(PetscLogFinalize()); 1633101409b8SToby Isaac 163433f85c2fSBarry Smith /* 163548dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 163648dd1dffSBarry Smith */ 16372e956fe4SStefano Zampini if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll()); 163837e93019SBarry Smith 16394028d114SSatish Balay if (petsc_history) { 16409566063dSJacob Faibussowitsch PetscCall(PetscCloseHistoryFile(&petsc_history)); 164102c9f0b5SLisandro Dalcin petsc_history = NULL; 1642e5c89e4eSSatish Balay } 16439566063dSJacob Faibussowitsch PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton)); 16449566063dSJacob Faibussowitsch PetscCall(PetscInfoDestroy()); 1645e5c89e4eSSatish Balay 164667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 164792f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1648e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 164992f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1650e5c89e4eSSatish Balay FILE *fd; 1651ed9cf6e9SBarry Smith int err; 1652e5c89e4eSSatish Balay 1653dc92acbaSJed Brown flg2 = PETSC_FALSE; 165492f119d6SBarry Smith flg3 = PETSC_FALSE; 16559566063dSJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL)); 16569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL)); 165792f119d6SBarry Smith fname[0] = 0; 16589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1)); 1659e5c89e4eSSatish Balay if (flg1 && fname[0]) { 16603ba16761SJacob Faibussowitsch PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank)); 16619371c9d4SSatish Balay fd = fopen(sname, "w"); 16629371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16639566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(fd)); 1664ed9cf6e9SBarry Smith err = fclose(fd); 166528b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 166692f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1667e5c89e4eSSatish Balay MPI_Comm local_comm; 1668e5c89e4eSSatish Balay 1669afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16709566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16719566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(stdout)); 16729566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1674e5c89e4eSSatish Balay } 1675e5c89e4eSSatish Balay fname[0] = 0; 16769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1)); 1677e5c89e4eSSatish Balay if (flg1 && fname[0]) { 16783ba16761SJacob Faibussowitsch PetscCall(PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank)); 16799371c9d4SSatish Balay fd = fopen(sname, "w"); 16809371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16819566063dSJacob Faibussowitsch PetscCall(PetscMallocView(fd)); 1682ed9cf6e9SBarry Smith err = fclose(fd); 168328b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 168492f119d6SBarry Smith } else if (flg1) { 168592f119d6SBarry Smith MPI_Comm local_comm; 168692f119d6SBarry Smith 1687afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16889566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16899566063dSJacob Faibussowitsch PetscCall(PetscMallocView(stdout)); 16909566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1692e5c89e4eSSatish Balay } 1693e5c89e4eSSatish Balay } 169467584ceeSBarry Smith #endif 169520e2c332SMatthew G. Knepley 16965486ca60SMatthew G. Knepley /* 16975486ca60SMatthew G. Knepley Close any open dynamic libraries 16985486ca60SMatthew G. Knepley */ 16999566063dSJacob Faibussowitsch PetscCall(PetscFinalize_DynamicLibraries()); 17005486ca60SMatthew G. Knepley 1701e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 17029566063dSJacob Faibussowitsch PetscCall(PetscOptionsDestroyDefault()); 1703e5c89e4eSSatish Balay 170471438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM) 170571438e86SJunchao Zhang if (PetscBeganNvshmem) { 17069566063dSJacob Faibussowitsch PetscCall(PetscNvshmemFinalize()); 170771438e86SJunchao Zhang PetscBeganNvshmem = PETSC_FALSE; 170871438e86SJunchao Zhang } 170971438e86SJunchao Zhang #endif 1710a0e72f99SJunchao Zhang 17119566063dSJacob Faibussowitsch PetscCall(PetscFreeMPIResources()); 1712e5c89e4eSSatish Balay 1713dbc8283eSBarry Smith /* 1714efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1715efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1716efb80d3cSBarry Smith 1717efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1718efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1719dbc8283eSBarry Smith */ 1720b770b1f6SSatish Balay { 1721dbc8283eSBarry Smith PetscCommCounter *counter; 1722dbc8283eSBarry Smith PetscMPIInt flg; 1723dbc8283eSBarry Smith MPI_Comm icomm; 17249371c9d4SSatish Balay union 17259371c9d4SSatish Balay { 17269371c9d4SSatish Balay MPI_Comm comm; 17279371c9d4SSatish Balay void *ptr; 17289371c9d4SSatish Balay } ucomm; 17299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 1730dbc8283eSBarry Smith if (flg) { 1731265f3f35SJed Brown icomm = ucomm.comm; 17329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 173328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1734dbc8283eSBarry Smith 17359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval)); 17369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 17379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1738dbc8283eSBarry Smith } 17399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 1740dbc8283eSBarry Smith if (flg) { 1741265f3f35SJed Brown icomm = ucomm.comm; 17429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 174328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1744dbc8283eSBarry Smith 17459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval)); 17469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 17479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1748dbc8283eSBarry Smith } 1749b770b1f6SSatish Balay } 1750dbc8283eSBarry Smith 17519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval)); 17529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval)); 17539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval)); 17549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval)); 175562e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_free_keyval(&Petsc_CreationIdx_keyval)); 175662e5d2d2SJDBetteridge PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Garbage_HMap_keyval)); 1757480cf27aSJed Brown 17585ea2b939SDuncan Campbell // Free keyvals which may be silently created by some routines 17595ea2b939SDuncan Campbell if (Petsc_SharedWD_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedWD_keyval)); 17605ea2b939SDuncan Campbell if (Petsc_SharedTmp_keyval != MPI_KEYVAL_INVALID) PetscCallMPI(MPI_Comm_free_keyval(&Petsc_SharedTmp_keyval)); 17615ea2b939SDuncan Campbell 17629566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen)); 17639566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout)); 17649566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr)); 17659566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock)); 1766ef19f930SBarry Smith 1767e5c89e4eSSatish Balay if (PetscBeganMPI) { 176899b1327fSBarry Smith PetscMPIInt flag; 17699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalized(&flag)); 177039a651e2SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 177139a651e2SJacob Faibussowitsch /* wait until the very last moment to disable error handling */ 177239a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_FALSE; 17739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalize()); 177439a651e2SJacob Faibussowitsch } else PetscErrorHandlingInitialized = PETSC_FALSE; 177539a651e2SJacob Faibussowitsch 1776e5c89e4eSSatish Balay /* 1777e5c89e4eSSatish Balay 1778e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1779e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1780e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1781e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1782e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 17830e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1784e5c89e4eSSatish Balay memory was not freed. 1785e5c89e4eSSatish Balay 1786e5c89e4eSSatish Balay */ 17879566063dSJacob Faibussowitsch PetscCall(PetscMallocClear()); 17889566063dSJacob Faibussowitsch PetscCall(PetscStackReset()); 1789a297a907SKarl Rupp 1790e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1791e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 17927ce81a4bSJacob Faibussowitsch #if defined(PETSC_USE_COVERAGE) 179356883f60SBarry Smith /* 179456883f60SBarry Smith flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the 179556883f60SBarry Smith gcov files are still being added to the directories as git tries to remove the directories. 179656883f60SBarry Smith */ 179756883f60SBarry Smith __gcov_flush(); 179856883f60SBarry Smith #endif 17991724198aSStefano Zampini /* To match PetscFunctionBegin() at the beginning of this function */ 18001724198aSStefano Zampini PetscStackClearTop; 18013ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 1802e5c89e4eSSatish Balay } 1803e5c89e4eSSatish Balay 180443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 1805d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame_(char *a, char *b) 1806d71ae5a4SJacob Faibussowitsch { 180743db4dbbSBarry Smith if (*a == *b) return 1; 180843db4dbbSBarry Smith if (*a + 32 == *b) return 1; 180943db4dbbSBarry Smith if (*a - 32 == *b) return 1; 181043db4dbbSBarry Smith return 0; 181143db4dbbSBarry Smith } 1812a70650f6SBarry Smith #endif 181343db4dbbSBarry Smith 181443db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 1815d71ae5a4SJacob Faibussowitsch PETSC_EXTERN int lsame(char *a, char *b) 1816d71ae5a4SJacob Faibussowitsch { 181743db4dbbSBarry Smith if (*a == *b) return 1; 181843db4dbbSBarry Smith if (*a + 32 == *b) return 1; 181943db4dbbSBarry Smith if (*a - 32 == *b) return 1; 182043db4dbbSBarry Smith return 0; 182143db4dbbSBarry Smith } 182243db4dbbSBarry Smith #endif 1823