1bc8a24ecSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay This file defines the initialization of PETSc, including PetscInitialize() 4e5c89e4eSSatish Balay */ 5af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I "petscsys.h" I*/ 6665c2dedSJed Brown #include <petscviewer.h> 7a0e72f99SJunchao Zhang 86edef35eSSatish Balay #if !defined(PETSC_HAVE_WINDOWS_COMPILERS) 96edef35eSSatish Balay #include <petsc/private/valgrind/valgrind.h> 106edef35eSSatish Balay #endif 116edef35eSSatish Balay 1285649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN) 1385649d77SJunchao Zhang #include <petsc/private/fortranimpl.h> 1485649d77SJunchao Zhang #endif 1585649d77SJunchao Zhang 1656883f60SBarry Smith #if defined(PETSC_USE_GCOV) 1756883f60SBarry Smith EXTERN_C_BEGIN 18aaf3846cSSatish Balay #if defined(PETSC_HAVE___GCOV_DUMP) 19aaf3846cSSatish Balay #define __gcov_flush(x) __gcov_dump(x) 20aaf3846cSSatish Balay #endif 2156883f60SBarry Smith void __gcov_flush(void); 2256883f60SBarry Smith EXTERN_C_END 2356883f60SBarry Smith #endif 248101f56cSMatthew Knepley 252d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 2695c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData; 272d53ad75SBarry Smith PetscFPT PetscFPTData = 0; 282d53ad75SBarry Smith #endif 292d53ad75SBarry Smith 3027104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 3116ad0300SBarry Smith #include <petscviewersaws.h> 32a6790183SBarry Smith #endif 33eb02082bSJunchao Zhang 34e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/ 35e5c89e4eSSatish Balay 3695c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history; 37e5c89e4eSSatish Balay 3895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void); 3995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void); 4095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode 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) 476de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED; 486de5d289SStefano Zampini #else 496de5d289SStefano Zampini PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0; 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; 56480cf27aSJed Brown 57e5c89e4eSSatish Balay /* 58e5c89e4eSSatish Balay Declare and set all the string names of the PETSc enums 59e5c89e4eSSatish Balay */ 6002c9f0b5SLisandro Dalcin const char *const PetscBools[] = {"FALSE", "TRUE", "PetscBool", "PETSC_", NULL}; 6102c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES", "OWN_POINTER", "USE_POINTER", "PetscCopyMode", "PETSC_", NULL}; 62e5c89e4eSSatish Balay 63ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE; 64ace3abfcSBarry Smith PetscBool PetscPreLoadingOn = PETSC_FALSE; 650f8e0872SSatish Balay 66a2f94806SJed Brown PetscInt PetscHotRegionDepth; 67a2f94806SJed Brown 686edef35eSSatish Balay PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE; 696edef35eSSatish Balay 70b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY) 71b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen; 72b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout; 73b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr; 74b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock; 75b22622e2STadeu Manoel #endif 76b22622e2STadeu Manoel 77e5c89e4eSSatish Balay /* 78945d1669SBarry Smith PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args 7972a42c3cSBarry Smith 8072a42c3cSBarry Smith Collective 8172a42c3cSBarry Smith 8272a42c3cSBarry Smith Level: advanced 8372a42c3cSBarry Smith 8495452b02SPatrick Sanan Notes: 85a56f64adSBarry Smith this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to 860f11a792SBarry Smith indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to 87a56f64adSBarry Smith be called multiple times from Julia without the problem of trying to initialize MPI more than once. 880f11a792SBarry Smith 89a56f64adSBarry Smith Developer Note: Turns off PETSc signal handling to allow Julia to manage signals 901ea5a559SBarry Smith 91db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()`, `PetscInitializeNoArguments()` 920f11a792SBarry Smith */ 939371c9d4SSatish Balay PetscErrorCode PetscInitializeNoPointers(int argc, char **args, const char *filename, const char *help) { 9472a42c3cSBarry Smith int myargc = argc; 9572a42c3cSBarry Smith char **myargs = args; 9672a42c3cSBarry Smith 9772a42c3cSBarry Smith PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&myargc, &myargs, filename, help)); 999566063dSJacob Faibussowitsch PetscCall(PetscPopSignalHandler()); 100df413903SBarry Smith PetscBeganMPI = PETSC_FALSE; 10127104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 10272a42c3cSBarry Smith } 10372a42c3cSBarry Smith 104f0865b08SBarry Smith /* 105a56f64adSBarry Smith Used by Julia interface to get communicator 106f0865b08SBarry Smith */ 1079371c9d4SSatish Balay PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm) { 108f0865b08SBarry Smith PetscFunctionBegin; 10927104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscValidPointer(comm, 1); 110f0865b08SBarry Smith *comm = PETSC_COMM_SELF; 111f0865b08SBarry Smith PetscFunctionReturn(0); 112f0865b08SBarry Smith } 113f0865b08SBarry Smith 114e5c89e4eSSatish Balay /*@C 115*811af0c4SBarry Smith PetscInitializeNoArguments - Calls `PetscInitialize()` from C/C++ without 116e5c89e4eSSatish Balay the command line arguments. 117e5c89e4eSSatish Balay 118e5c89e4eSSatish Balay Collective 119e5c89e4eSSatish Balay 120e5c89e4eSSatish Balay Level: advanced 121e5c89e4eSSatish Balay 122db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeFortran()` 123e5c89e4eSSatish Balay @*/ 1249371c9d4SSatish Balay PetscErrorCode PetscInitializeNoArguments(void) { 125e5c89e4eSSatish Balay int argc = 0; 12602c9f0b5SLisandro Dalcin char **args = NULL; 127e5c89e4eSSatish Balay 128e5c89e4eSSatish Balay PetscFunctionBegin; 1299566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, NULL)); 13039a651e2SJacob Faibussowitsch PetscFunctionReturn(0); 131e5c89e4eSSatish Balay } 132e5c89e4eSSatish Balay 133e5c89e4eSSatish Balay /*@ 134e5c89e4eSSatish Balay PetscInitialized - Determine whether PETSc is initialized. 135e5c89e4eSSatish Balay 13693b6d2d1SJed Brown Level: beginner 137e5c89e4eSSatish Balay 138db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()` 139e5c89e4eSSatish Balay @*/ 1409371c9d4SSatish Balay PetscErrorCode PetscInitialized(PetscBool *isInitialized) { 14127104ee2SJacob Faibussowitsch PetscFunctionBegin; 14227104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscValidBoolPointer(isInitialized, 1); 143e5c89e4eSSatish Balay *isInitialized = PetscInitializeCalled; 14427104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 145e5c89e4eSSatish Balay } 146e5c89e4eSSatish Balay 147e5c89e4eSSatish Balay /*@ 148*811af0c4SBarry Smith PetscFinalized - Determine whether `PetscFinalize()` has been called yet 149e5c89e4eSSatish Balay 150e5c89e4eSSatish Balay Level: developer 151e5c89e4eSSatish Balay 152db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscInitializeNoArguments()`, `PetscInitializeFortran()` 153e5c89e4eSSatish Balay @*/ 1549371c9d4SSatish Balay PetscErrorCode PetscFinalized(PetscBool *isFinalized) { 15527104ee2SJacob Faibussowitsch PetscFunctionBegin; 15627104ee2SJacob Faibussowitsch if (!PetscFinalizeCalled) PetscValidBoolPointer(isFinalized, 1); 157e5c89e4eSSatish Balay *isFinalized = PetscFinalizeCalled; 15827104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 159e5c89e4eSSatish Balay } 160e5c89e4eSSatish Balay 16157171095SVaclav Hapla PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char[]); 162e5c89e4eSSatish Balay 163e5c89e4eSSatish Balay /* 164e5c89e4eSSatish Balay This function is the MPI reduction operation used to compute the sum of the 165e5c89e4eSSatish Balay first half of the datatype and the max of the second half. 166e5c89e4eSSatish Balay */ 167367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0; 168e5c89e4eSSatish Balay 1699371c9d4SSatish Balay PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in, void *out, int *cnt, MPI_Datatype *datatype) { 170e5c89e4eSSatish Balay PetscInt *xin = (PetscInt *)in, *xout = (PetscInt *)out, i, count = *cnt; 171e5c89e4eSSatish Balay 172e5c89e4eSSatish Balay PetscFunctionBegin; 173e5c89e4eSSatish Balay if (*datatype != MPIU_2INT) { 174e5c89e4eSSatish Balay (*PetscErrorPrintf)("Can only handle MPIU_2INT data types"); 17541e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 176e5c89e4eSSatish Balay } 177e5c89e4eSSatish Balay 178e5c89e4eSSatish Balay for (i = 0; i < count; i++) { 179e5c89e4eSSatish Balay xout[2 * i] = PetscMax(xout[2 * i], xin[2 * i]); 180e5c89e4eSSatish Balay xout[2 * i + 1] += xin[2 * i + 1]; 181e5c89e4eSSatish Balay } 182812af9f3SBarry Smith PetscFunctionReturnVoid(); 183e5c89e4eSSatish Balay } 184e5c89e4eSSatish Balay 185e5c89e4eSSatish Balay /* 186e5c89e4eSSatish Balay Returns the max of the first entry owned by this processor and the 187e5c89e4eSSatish Balay sum of the second entry. 188b693b147SBarry Smith 18976f543a4SJed Brown The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero 190367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths 191b693b147SBarry Smith there would be no place to store the both needed results. 192e5c89e4eSSatish Balay */ 1939371c9d4SSatish Balay PetscErrorCode PetscMaxSum(MPI_Comm comm, const PetscInt sizes[], PetscInt *max, PetscInt *sum) { 194e5c89e4eSSatish Balay PetscFunctionBegin; 195d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 196d6e4c47cSJed Brown { 1979371c9d4SSatish Balay struct { 1989371c9d4SSatish Balay PetscInt max, sum; 1999371c9d4SSatish Balay } work; 2009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce_scatter_block((void *)sizes, &work, 1, MPIU_2INT, MPIU_MAXSUM_OP, comm)); 201d6e4c47cSJed Brown *max = work.max; 202d6e4c47cSJed Brown *sum = work.sum; 203d6e4c47cSJed Brown } 204d6e4c47cSJed Brown #else 205d6e4c47cSJed Brown { 206d6e4c47cSJed Brown PetscMPIInt size, rank; 2079371c9d4SSatish Balay struct { 2089371c9d4SSatish Balay PetscInt max, sum; 2099371c9d4SSatish Balay } * work; 2109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &work)); 2131c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((void *)sizes, work, size, MPIU_2INT, MPIU_MAXSUM_OP, comm)); 2146ac3741eSJed Brown *max = work[rank].max; 2156ac3741eSJed Brown *sum = work[rank].sum; 2169566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 217d6e4c47cSJed Brown } 218d6e4c47cSJed Brown #endif 219e5c89e4eSSatish Balay PetscFunctionReturn(0); 220e5c89e4eSSatish Balay } 221e5c89e4eSSatish Balay 222e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/ 223e5c89e4eSSatish Balay 2249e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16) 225613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 226613bf2b2SPierre Jolivet #include <quadmath.h> 227613bf2b2SPierre Jolivet #endif 2289e517322SPierre Jolivet MPI_Op MPIU_SUM___FP16___FLOAT128 = 0; 229de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 23006a205a8SBarry Smith MPI_Op MPIU_SUM = 0; 231613bf2b2SPierre Jolivet #endif 232e5c89e4eSSatish Balay 2339371c9d4SSatish Balay PETSC_EXTERN void MPIAPI PetscSum_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) { 234e5c89e4eSSatish Balay PetscInt i, count = *cnt; 235e5c89e4eSSatish Balay 236e5c89e4eSSatish Balay PetscFunctionBegin; 2377c2de775SJed Brown if (*datatype == MPIU_REAL) { 238e2e03761SBarry Smith PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 239a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] += xin[i]; 2407c2de775SJed Brown } 2417c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 242c5481aeeSJose E. Roman else if (*datatype == MPIU_COMPLEX) { 2437c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 244a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] += xin[i]; 2457c2de775SJed Brown } 2467c2de775SJed Brown #endif 247613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 248613bf2b2SPierre Jolivet else if (*datatype == MPIU___FLOAT128) { 249613bf2b2SPierre Jolivet __float128 *xin = (__float128 *)in, *xout = (__float128 *)out; 250613bf2b2SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 2519371c9d4SSatish Balay } else if (*datatype == MPIU___COMPLEX128) { 252613bf2b2SPierre Jolivet __complex128 *xin = (__complex128 *)in, *xout = (__complex128 *)out; 253613bf2b2SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 254613bf2b2SPierre Jolivet } 255613bf2b2SPierre Jolivet #endif 2569e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) 2579e517322SPierre Jolivet else if (*datatype == MPIU___FP16) { 2589e517322SPierre Jolivet __fp16 *xin = (__fp16 *)in, *xout = (__fp16 *)out; 2599e517322SPierre Jolivet for (i = 0; i < count; i++) xout[i] += xin[i]; 2609e517322SPierre Jolivet } 2619e517322SPierre Jolivet #endif 2627c2de775SJed Brown else { 2639e517322SPierre Jolivet #if !defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_HAVE_REAL___FP16) 2647c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 2659e517322SPierre Jolivet #elif !defined(PETSC_HAVE_REAL___FP16) 266613bf2b2SPierre Jolivet (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types"); 2679e517322SPierre Jolivet #elif !defined(PETSC_HAVE_REAL___FLOAT128) 2689e517322SPierre Jolivet (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, or MPIU___FP16 data types"); 2699e517322SPierre Jolivet #else 2709e517322SPierre Jolivet (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, MPIU___COMPLEX128, or MPIU___FP16 data types"); 271613bf2b2SPierre Jolivet #endif 27241e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 273e2e03761SBarry Smith } 274812af9f3SBarry Smith PetscFunctionReturnVoid(); 275e5c89e4eSSatish Balay } 276e5c89e4eSSatish Balay #endif 277e5c89e4eSSatish Balay 278570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 279d9822059SBarry Smith MPI_Op MPIU_MAX = 0; 280d9822059SBarry Smith MPI_Op MPIU_MIN = 0; 281d9822059SBarry Smith 2829371c9d4SSatish Balay PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) { 283d9822059SBarry Smith PetscInt i, count = *cnt; 284d9822059SBarry Smith 285d9822059SBarry Smith PetscFunctionBegin; 2867c2de775SJed Brown if (*datatype == MPIU_REAL) { 2878c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 288a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]); 2897c2de775SJed Brown } 2907c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2917c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2927c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 293ad540459SPierre Jolivet for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 2947c2de775SJed Brown } 2957c2de775SJed Brown #endif 2967c2de775SJed Brown else { 2977c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 29841e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 2998c764dc5SJose Roman } 300d9822059SBarry Smith PetscFunctionReturnVoid(); 301d9822059SBarry Smith } 302d9822059SBarry Smith 3039371c9d4SSatish Balay PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) { 304d9822059SBarry Smith PetscInt i, count = *cnt; 305d9822059SBarry Smith 306d9822059SBarry Smith PetscFunctionBegin; 3077c2de775SJed Brown if (*datatype == MPIU_REAL) { 3088c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 309a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]); 3107c2de775SJed Brown } 3117c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 3127c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 3137c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 314ad540459SPierre Jolivet for (i = 0; i < count; i++) xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 3157c2de775SJed Brown } 3167c2de775SJed Brown #endif 3177c2de775SJed Brown else { 3188c764dc5SJose Roman (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"); 31941e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 3208c764dc5SJose Roman } 321d9822059SBarry Smith PetscFunctionReturnVoid(); 322d9822059SBarry Smith } 323d9822059SBarry Smith #endif 324d9822059SBarry Smith 325480cf27aSJed Brown /* 326480cf27aSJed Brown Private routine to delete internal tag/name counter storage when a communicator is freed. 327480cf27aSJed Brown 328ff0e51ddSBarry 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. 329480cf27aSJed Brown 33012801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 331480cf27aSJed Brown 332480cf27aSJed Brown */ 3339371c9d4SSatish Balay PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state) { 33405342407SJunchao Zhang PetscCommCounter *counter = (PetscCommCounter *)count_val; 33557f21012SBarry Smith struct PetscCommStash *comms = counter->comms, *pcomm; 336480cf27aSJed Brown 337480cf27aSJed Brown PetscFunctionBegin; 3389566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm)); 3399566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter->iflags)); 34057f21012SBarry Smith while (comms) { 3419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&comms->comm)); 34257f21012SBarry Smith pcomm = comms; 34357f21012SBarry Smith comms = comms->next; 3449566063dSJacob Faibussowitsch PetscCall(PetscFree(pcomm)); 34557f21012SBarry Smith } 3469566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter)); 347480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 348480cf27aSJed Brown } 349480cf27aSJed Brown 350480cf27aSJed Brown /* 35147435625SJed Brown This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user 352da3039f7SJed Brown calls MPI_Comm_free(). 353da3039f7SJed Brown 354da3039f7SJed Brown This is the only entry point for breaking the links between inner and outer comms. 355480cf27aSJed Brown 356ff0e51ddSBarry Smith This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator. 357480cf27aSJed Brown 35812801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 359480cf27aSJed Brown 360480cf27aSJed Brown */ 3619371c9d4SSatish Balay PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) { 3629371c9d4SSatish Balay union 363480cf27aSJed Brown { 3649371c9d4SSatish Balay MPI_Comm comm; 3659371c9d4SSatish Balay void *ptr; 3669371c9d4SSatish Balay } icomm; 367480cf27aSJed Brown 368480cf27aSJed Brown PetscFunctionBegin; 36912801b39SBarry Smith if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval"); 370ec4fadc2SJed Brown icomm.ptr = attr_val; 37176bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 37233779a13SJunchao Zhang /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */ 37333779a13SJunchao Zhang PetscMPIInt flg; 3749371c9d4SSatish Balay union 3759371c9d4SSatish Balay { 3769371c9d4SSatish Balay MPI_Comm comm; 3779371c9d4SSatish Balay void *ptr; 3789371c9d4SSatish Balay } ocomm; 3799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg)); 38033779a13SJunchao Zhang if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute"); 38133779a13SJunchao 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"); 38233779a13SJunchao Zhang } 3839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval)); 3849566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm)); 385da3039f7SJed Brown PetscFunctionReturn(MPI_SUCCESS); 386b89831e5SBarry Smith } 387da3039f7SJed Brown 388da3039f7SJed Brown /* 38933779a13SJunchao Zhang * This is invoked on the inner comm when Petsc_InnerComm_Attr_Delete_Fn calls MPI_Comm_delete_attr(). It should not be reached any other way. 390da3039f7SJed Brown */ 3919371c9d4SSatish Balay PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) { 392da3039f7SJed Brown PetscFunctionBegin; 3939566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm)); 394480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 395480cf27aSJed Brown } 396480cf27aSJed Brown 39733779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *); 39842218b76SBarry Smith 399951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 4008cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *); 4018cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 4028cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 403e39fd77fSBarry Smith #endif 404e39fd77fSBarry Smith 4050f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE; 40612801b39SBarry Smith 407eb27c7c8SSatish Balay PETSC_INTERN int PetscGlobalArgc; 408eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs; 4096ae9a8a6SBarry Smith int PetscGlobalArgc = 0; 41002c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL; 411dff31646SBarry Smith PetscSegBuffer PetscCitationsList; 412e5c89e4eSSatish Balay 4139371c9d4SSatish Balay PetscErrorCode PetscCitationsInitialize(void) { 414051e4cf2SJed Brown PetscFunctionBegin; 4159566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList)); 41630c35bf2SSatish Balay 41730c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\ 41830c35bf2SSatish Balay Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\ 41930c35bf2SSatish Balay and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\ 4203f81df79SBarry Smith and Victor Eijkhout and Jacob Faibussowitsch and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\ 42130c35bf2SSatish Balay and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\ 42230c35bf2SSatish Balay and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\ 42330c35bf2SSatish Balay and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\ 42430c35bf2SSatish Balay and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\ 42530c35bf2SSatish Balay Title = {{PETSc/TAO} Users Manual},\n\ 42630c35bf2SSatish Balay Number = {ANL-21/39 - Revision 3.17},\n\ 42730c35bf2SSatish Balay Institution = {Argonne National Laboratory},\n\ 4289371c9d4SSatish Balay Year = {2022}\n}\n", 4299371c9d4SSatish Balay NULL)); 43030c35bf2SSatish Balay 43130c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\ 43230c35bf2SSatish Balay Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\ 43330c35bf2SSatish Balay Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\ 43430c35bf2SSatish Balay Booktitle = {Modern Software Tools in Scientific Computing},\n\ 43530c35bf2SSatish Balay Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\ 43630c35bf2SSatish Balay Pages = {163--202},\n\ 43730c35bf2SSatish Balay Publisher = {Birkh{\\\"{a}}user Press},\n\ 4389371c9d4SSatish Balay Year = {1997}\n}\n", 4399371c9d4SSatish Balay NULL)); 44030c35bf2SSatish Balay 441051e4cf2SJed Brown PetscFunctionReturn(0); 442051e4cf2SJed Brown } 443e5c89e4eSSatish Balay 4442d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */ 4452d747510SLisandro Dalcin 4469371c9d4SSatish Balay PetscErrorCode PetscSetProgramName(const char name[]) { 4472d747510SLisandro Dalcin PetscFunctionBegin; 4489566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(programname, name, sizeof(programname))); 4492d747510SLisandro Dalcin PetscFunctionReturn(0); 4502d747510SLisandro Dalcin } 4512d747510SLisandro Dalcin 4522d747510SLisandro Dalcin /*@C 4532d747510SLisandro Dalcin PetscGetProgramName - Gets the name of the running program. 4542d747510SLisandro Dalcin 4552d747510SLisandro Dalcin Not Collective 4562d747510SLisandro Dalcin 4572d747510SLisandro Dalcin Input Parameter: 4582d747510SLisandro Dalcin . len - length of the string name 4592d747510SLisandro Dalcin 4602d747510SLisandro Dalcin Output Parameter: 461*811af0c4SBarry Smith . name - the name of the running program, provide a string of length `PETSC_MAX_PATH_LEN` 4622d747510SLisandro Dalcin 4632d747510SLisandro Dalcin Level: advanced 4642d747510SLisandro Dalcin 4652d747510SLisandro Dalcin @*/ 4669371c9d4SSatish Balay PetscErrorCode PetscGetProgramName(char name[], size_t len) { 4672d747510SLisandro Dalcin PetscFunctionBegin; 4689566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(name, programname, len)); 4692d747510SLisandro Dalcin PetscFunctionReturn(0); 4702d747510SLisandro Dalcin } 4712d747510SLisandro Dalcin 472e5c89e4eSSatish Balay /*@C 473e5c89e4eSSatish Balay PetscGetArgs - Allows you to access the raw command line arguments anywhere 474*811af0c4SBarry Smith after PetscInitialize() is called but before `PetscFinalize()`. 475e5c89e4eSSatish Balay 476e5c89e4eSSatish Balay Not Collective 477e5c89e4eSSatish Balay 478e5c89e4eSSatish Balay Output Parameters: 479e5c89e4eSSatish Balay + argc - count of number of command line arguments 480e5c89e4eSSatish Balay - args - the command line arguments 481e5c89e4eSSatish Balay 482e5c89e4eSSatish Balay Level: intermediate 483e5c89e4eSSatish Balay 484e5c89e4eSSatish Balay Notes: 485e5c89e4eSSatish Balay This is usually used to pass the command line arguments into other libraries 486e5c89e4eSSatish Balay that are called internally deep in PETSc or the application. 487e5c89e4eSSatish Balay 488f177e3b1SBarry Smith The first argument contains the program name as is normal for C arguments. 489f177e3b1SBarry Smith 490db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()` 491e5c89e4eSSatish Balay @*/ 4929371c9d4SSatish Balay PetscErrorCode PetscGetArgs(int *argc, char ***args) { 493e5c89e4eSSatish Balay PetscFunctionBegin; 49439a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 495e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 496e5c89e4eSSatish Balay *args = PetscGlobalArgs; 497e5c89e4eSSatish Balay PetscFunctionReturn(0); 498e5c89e4eSSatish Balay } 499e5c89e4eSSatish Balay 500793721a6SBarry Smith /*@C 501793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 502*811af0c4SBarry Smith after `PetscInitialize()` is called but before `PetscFinalize()`. 503793721a6SBarry Smith 504793721a6SBarry Smith Not Collective 505793721a6SBarry Smith 506793721a6SBarry Smith Output Parameters: 507793721a6SBarry Smith . args - the command line arguments 508793721a6SBarry Smith 509793721a6SBarry Smith Level: intermediate 510793721a6SBarry Smith 511793721a6SBarry Smith Notes: 512793721a6SBarry Smith This does NOT start with the program name and IS null terminated (final arg is void) 513793721a6SBarry Smith 514db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()` 515793721a6SBarry Smith @*/ 5169371c9d4SSatish Balay PetscErrorCode PetscGetArguments(char ***args) { 517793721a6SBarry Smith PetscInt i, argc = PetscGlobalArgc; 518793721a6SBarry Smith 519793721a6SBarry Smith PetscFunctionBegin; 52039a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 5219371c9d4SSatish Balay if (!argc) { 5229371c9d4SSatish Balay *args = NULL; 5239371c9d4SSatish Balay PetscFunctionReturn(0); 5249371c9d4SSatish Balay } 5259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(argc, args)); 5269566063dSJacob Faibussowitsch for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i])); 5272d747510SLisandro Dalcin (*args)[argc - 1] = NULL; 528793721a6SBarry Smith PetscFunctionReturn(0); 529793721a6SBarry Smith } 530793721a6SBarry Smith 531793721a6SBarry Smith /*@C 532*811af0c4SBarry Smith PetscFreeArguments - Frees the memory obtained with `PetscGetArguments()` 533793721a6SBarry Smith 534793721a6SBarry Smith Not Collective 535793721a6SBarry Smith 536793721a6SBarry Smith Output Parameters: 537793721a6SBarry Smith . args - the command line arguments 538793721a6SBarry Smith 539793721a6SBarry Smith Level: intermediate 540793721a6SBarry Smith 541db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()` 542793721a6SBarry Smith @*/ 5439371c9d4SSatish Balay PetscErrorCode PetscFreeArguments(char **args) { 54439a651e2SJacob Faibussowitsch PetscFunctionBegin; 54539a651e2SJacob Faibussowitsch if (args) { 546793721a6SBarry Smith PetscInt i = 0; 547793721a6SBarry Smith 5489566063dSJacob Faibussowitsch while (args[i]) PetscCall(PetscFree(args[i++])); 5499566063dSJacob Faibussowitsch PetscCall(PetscFree(args)); 55039a651e2SJacob Faibussowitsch } 551793721a6SBarry Smith PetscFunctionReturn(0); 552793721a6SBarry Smith } 553793721a6SBarry Smith 55427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 55530befbd2SBarry Smith #include <petscconfiginfo.h> 55630befbd2SBarry Smith 5579371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) { 55827104ee2SJacob Faibussowitsch PetscFunctionBegin; 55911525c0dSBarry Smith if (!PetscGlobalRank) { 56030befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64]; 56111525c0dSBarry Smith int port; 562ffbd1cfbSBarry Smith PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE; 56311525c0dSBarry Smith size_t applinelen, introlen; 564ffbd1cfbSBarry Smith char sawsurl[256]; 56511525c0dSBarry Smith 5669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg)); 56711525c0dSBarry Smith if (flg) { 56811525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 56911525c0dSBarry Smith 5709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL)); 57111525c0dSBarry Smith if (sawslog[0]) { 572792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog)); 57311525c0dSBarry Smith } else { 574792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL)); 57511525c0dSBarry Smith } 57611525c0dSBarry Smith } 5779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg)); 57848a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert)); 5799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL)); 580ffbd1cfbSBarry Smith if (selectport) { 581792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 582792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 583ffbd1cfbSBarry Smith } else { 5849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg)); 58548a46eb9SPierre Jolivet if (flg) PetscCallSAWs(SAWs_Set_Port, (port)); 586ffbd1cfbSBarry Smith } 5879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg)); 58811525c0dSBarry Smith if (flg) { 589792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 5909566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(root, ".", &rootlocal)); 5919c1e0ce8SBarry Smith } else { 5929566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg)); 5939c1e0ce8SBarry Smith if (flg) { 5949566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root))); 595792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 5969c1e0ce8SBarry Smith } 59711525c0dSBarry Smith } 5989566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2)); 59911525c0dSBarry Smith if (flg2) { 60011525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 60128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option"); 6029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root)); 6039566063dSJacob Faibussowitsch PetscCall(PetscTestDirectory(jsdir, 'r', &flg)); 60428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory"); 605792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Local_Header, ()); 60611525c0dSBarry Smith } 6079566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(programname, sizeof(programname))); 6089566063dSJacob Faibussowitsch PetscCall(PetscStrlen(help, &applinelen)); 60911525c0dSBarry Smith introlen = 4096 + applinelen; 61030a8c9c0SSurtai Han applinelen += 1024; 6119566063dSJacob Faibussowitsch PetscCall(PetscMalloc(applinelen, &appline)); 6129566063dSJacob Faibussowitsch PetscCall(PetscMalloc(introlen, &intro)); 61311525c0dSBarry Smith 61411525c0dSBarry Smith if (rootlocal) { 6159566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname)); 6169566063dSJacob Faibussowitsch PetscCall(PetscTestFile(appline, 'r', &rootlocal)); 61711525c0dSBarry Smith } 6189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetAll(NULL, &options)); 61911525c0dSBarry Smith if (rootlocal && help) { 6209566063dSJacob 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)); 62111525c0dSBarry Smith } else if (help) { 6229566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help)); 62311525c0dSBarry Smith } else { 6249566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options)); 62511525c0dSBarry Smith } 6269566063dSJacob Faibussowitsch PetscCall(PetscFree(options)); 6279566063dSJacob Faibussowitsch PetscCall(PetscGetVersion(version, sizeof(version))); 6289371c9d4SSatish Balay PetscCall(PetscSNPrintf(intro, introlen, 6299371c9d4SSatish Balay "<body>\n" 630a17b96a8SKyle 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" 631df62c222SBarry 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" 6329371c9d4SSatish Balay "%s", 6339371c9d4SSatish Balay version, petscconfigureoptions, appline)); 634792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro)); 6359566063dSJacob Faibussowitsch PetscCall(PetscFree(intro)); 6369566063dSJacob Faibussowitsch PetscCall(PetscFree(appline)); 637ffbd1cfbSBarry Smith if (selectport) { 638aa573868SBarry Smith PetscBool silent; 6397d812c46SBarry Smith 6407d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 64139a651e2SJacob Faibussowitsch while (SAWs_Initialize()) { 642792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 643792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 6447d812c46SBarry Smith } 6457d812c46SBarry Smith 6469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL)); 647aa573868SBarry Smith if (!silent) { 648792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl)); 6499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl)); 650ffbd1cfbSBarry Smith } 6517d812c46SBarry Smith } else { 652792fecdfSBarry Smith PetscCallSAWs(SAWs_Initialize, ()); 653aa573868SBarry Smith } 6549566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister("@TechReport{ saws,\n" 6550af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 6560af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 6570af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 6589371c9d4SSatish Balay " Year = 2013\n}\n", 6599371c9d4SSatish Balay NULL)); 66011525c0dSBarry Smith } 66111525c0dSBarry Smith PetscFunctionReturn(0); 66211525c0dSBarry Smith } 66311525c0dSBarry Smith #endif 66411525c0dSBarry Smith 6654dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 6669371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) { 6674dfee713SSatish Balay PetscFunctionBegin; 6684dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 6694dfee713SSatish Balay /* see MPI.py for details on this bug */ 6704dfee713SSatish Balay (void)setenv("HWLOC_COMPONENTS", "-x86", 1); 6714dfee713SSatish Balay #endif 6724dfee713SSatish Balay PetscFunctionReturn(0); 6734dfee713SSatish Balay } 6744dfee713SSatish Balay 675a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS) 676a56f64adSBarry Smith #include <adios.h> 67722580e64SBarry Smith #include <adios_read.h> 6787b56e58cSSatish Balay int64_t Petsc_adios_group; 679a56f64adSBarry Smith #endif 680a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP) 681cd1b99a6SBarry Smith #include <omp.h> 682f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 683cd1b99a6SBarry Smith #endif 684a56f64adSBarry Smith 685a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_DEVICE) 686a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h> 687a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA) 688a4af0ceeSJacob Faibussowitsch // REMOVE ME 689a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL; 690a4af0ceeSJacob Faibussowitsch #endif 691a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP) 692a4af0ceeSJacob Faibussowitsch // REMOVE ME 693a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL; 694a4af0ceeSJacob Faibussowitsch #endif 695a4af0ceeSJacob Faibussowitsch #endif 696a4af0ceeSJacob Faibussowitsch 69727104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H) 698bc8a24ecSBarry Smith #include <dlfcn.h> 699bc8a24ecSBarry Smith #endif 700a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 701a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void); 702a4af0ceeSJacob Faibussowitsch #endif 703a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 704a4af0ceeSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViennaCLInit(); 705a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE; 706a4af0ceeSJacob Faibussowitsch #endif 707bc8a24ecSBarry Smith 708660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE; 709660278c0SBarry Smith 71085649d77SJunchao Zhang /* 71185649d77SJunchao Zhang PetscInitialize_Common - shared code between C and Fortran initialization 71285649d77SJunchao Zhang 71385649d77SJunchao Zhang prog: program name 71402101c96SSatish Balay file: optional PETSc database file name. Might be in Fortran string format when 'ftn' is true 71585649d77SJunchao Zhang help: program help message 71602101c96SSatish Balay ftn: is it called from Fortran initilization (petscinitializef_)? 71785649d77SJunchao Zhang readarguments,len: used when fortran is true 71885649d77SJunchao Zhang */ 7199371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len) { 72085649d77SJunchao Zhang PetscMPIInt size; 72185649d77SJunchao Zhang PetscBool flg = PETSC_TRUE; 72285649d77SJunchao Zhang char hostname[256]; 72385649d77SJunchao Zhang 72427104ee2SJacob Faibussowitsch PetscFunctionBegin; 72527104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(0); 72639a651e2SJacob Faibussowitsch /* these must be initialized in a routine, not as a constant declaration */ 72739a651e2SJacob Faibussowitsch PETSC_STDOUT = stdout; 72839a651e2SJacob Faibussowitsch PETSC_STDERR = stderr; 72939a651e2SJacob Faibussowitsch 7309566063dSJacob Faibussowitsch /* PetscCall can be used from now */ 73139a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_TRUE; 73239a651e2SJacob Faibussowitsch 73385649d77SJunchao Zhang /* 73485649d77SJunchao Zhang The checking over compatible runtime libraries is complicated by the MPI ABI initiative 73585649d77SJunchao Zhang https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 73685649d77SJunchao Zhang MPICH v3.1 (Released February 2014) 73785649d77SJunchao Zhang IBM MPI v2.1 (December 2014) 73885649d77SJunchao Zhang Intel MPI Library v5.0 (2014) 73985649d77SJunchao Zhang Cray MPT v7.0.0 (June 2014) 74085649d77SJunchao Zhang As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 74185649d77SJunchao Zhang listed above and since that time are compatible. 74285649d77SJunchao Zhang 74385649d77SJunchao Zhang Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 74485649d77SJunchao Zhang at compile time or runtime. Thus we will need to systematically track the allowed versions 74585649d77SJunchao Zhang and how they are represented in the mpi.h and MPI_Get_library_version() output in order 74685649d77SJunchao Zhang to perform the checking. 74785649d77SJunchao Zhang 74885649d77SJunchao Zhang Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 74985649d77SJunchao Zhang 75085649d77SJunchao Zhang Questions: 75185649d77SJunchao Zhang 75285649d77SJunchao Zhang Should the checks for ABI incompatibility be only on the major version number below? 75385649d77SJunchao Zhang Presumably the output to stderr will be removed before a release. 75485649d77SJunchao Zhang */ 75585649d77SJunchao Zhang 75685649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 75785649d77SJunchao Zhang { 75885649d77SJunchao Zhang char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 75985649d77SJunchao Zhang PetscMPIInt mpilibraryversionlength; 76039a651e2SJacob Faibussowitsch 7619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength)); 76285649d77SJunchao Zhang /* check for MPICH versions before MPI ABI initiative */ 76385649d77SJunchao Zhang #if defined(MPICH_VERSION) 76485649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000 76585649d77SJunchao Zhang { 76685649d77SJunchao Zhang char *ver, *lf; 76785649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 76839a651e2SJacob Faibussowitsch 7699566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver)); 77039a651e2SJacob Faibussowitsch if (ver) { 7719566063dSJacob Faibussowitsch PetscCall(PetscStrchr(ver, '\n', &lf)); 77239a651e2SJacob Faibussowitsch if (lf) { 77385649d77SJunchao Zhang *lf = 0; 7749566063dSJacob Faibussowitsch PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg)); 77585649d77SJunchao Zhang } 77685649d77SJunchao Zhang } 77785649d77SJunchao Zhang if (!flg) { 778d5b396e8SSatish Balay PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION)); 77985649d77SJunchao Zhang flg = PETSC_TRUE; 78085649d77SJunchao Zhang } 78185649d77SJunchao Zhang } 78285649d77SJunchao Zhang #endif 78385649d77SJunchao Zhang /* check for OpenMPI version, it is not part of the MPI ABI initiative (is it part of another initiative that needs to be handled?) */ 78485649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION) 78585649d77SJunchao Zhang { 78685649d77SJunchao Zhang char *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf; 78785649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 78885649d77SJunchao Zhang #define PSTRSZ 2 78985649d77SJunchao Zhang char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"}; 79085649d77SJunchao Zhang char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "}; 79185649d77SJunchao Zhang int i; 79285649d77SJunchao Zhang for (i = 0; i < PSTRSZ; i++) { 7939566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver)); 79439a651e2SJacob Faibussowitsch if (ver) { 7959566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION)); 7969566063dSJacob Faibussowitsch PetscCall(PetscStrstr(ver, bs, &bsf)); 79739a651e2SJacob Faibussowitsch if (bsf) flg = PETSC_TRUE; 79885649d77SJunchao Zhang break; 79985649d77SJunchao Zhang } 80085649d77SJunchao Zhang } 80185649d77SJunchao Zhang if (!flg) { 8027d3de750SJacob Faibussowitsch PetscInfo(NULL, "PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n", mpilibraryversion, OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION); 80385649d77SJunchao Zhang flg = PETSC_TRUE; 80485649d77SJunchao Zhang } 80585649d77SJunchao Zhang } 80685649d77SJunchao Zhang #endif 80785649d77SJunchao Zhang } 80885649d77SJunchao Zhang #endif 80985649d77SJunchao Zhang 8106f5d4113SSatish Balay #if PetscDefined(HAVE_DLSYM) && defined(__USE_GNU) 81185649d77SJunchao Zhang /* These symbols are currently in the OpenMPI and MPICH libraries; they may not always be, in that case the test will simply not detect the problem */ 81239a651e2SJacob Faibussowitsch PetscCheck(!dlsym(RTLD_DEFAULT, "ompi_mpi_init") || !dlsym(RTLD_DEFAULT, "MPID_Abort"), PETSC_COMM_SELF, PETSC_ERR_MPI_LIB_INCOMP, "Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly"); 81385649d77SJunchao Zhang #endif 81485649d77SJunchao Zhang 81585649d77SJunchao Zhang /* on Windows - set printf to default to printing 2 digit exponents */ 81685649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT) 81785649d77SJunchao Zhang _set_output_format(_TWO_DIGIT_EXPONENT); 81885649d77SJunchao Zhang #endif 81985649d77SJunchao Zhang 8209566063dSJacob Faibussowitsch PetscCall(PetscOptionsCreateDefault()); 82185649d77SJunchao Zhang 82285649d77SJunchao Zhang PetscFinalizeCalled = PETSC_FALSE; 82385649d77SJunchao Zhang 8249566063dSJacob Faibussowitsch PetscCall(PetscSetProgramName(prog)); 8259566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen)); 8269566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout)); 8279566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr)); 8289566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscCommSpinLock)); 82985649d77SJunchao Zhang 83085649d77SJunchao Zhang if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 8319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN)); 83285649d77SJunchao Zhang 83385649d77SJunchao Zhang if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) { 8349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS)); 8359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE)); 83685649d77SJunchao Zhang } 83785649d77SJunchao Zhang 83885649d77SJunchao Zhang /* Done after init due to a bug in MPICH-GM? */ 8399566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 84085649d77SJunchao Zhang 8419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank)); 8429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize)); 84385649d77SJunchao Zhang 84485649d77SJunchao Zhang MPIU_BOOL = MPI_INT; 84585649d77SJunchao Zhang MPIU_ENUM = MPI_INT; 84685649d77SJunchao Zhang MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64; 84785649d77SJunchao Zhang if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 84885649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 84985649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG) 85085649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 85185649d77SJunchao Zhang #endif 85239a651e2SJacob Faibussowitsch else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t"); 85385649d77SJunchao Zhang 85485649d77SJunchao Zhang /* 85585649d77SJunchao Zhang Initialized the global complex variable; this is because with 85685649d77SJunchao Zhang shared libraries the constructors for global variables 85785649d77SJunchao Zhang are not called; at least on IRIX. 85885649d77SJunchao Zhang */ 85985649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX) 86085649d77SJunchao Zhang { 86185649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 86285649d77SJunchao Zhang PetscComplex ic(0.0, 1.0); 86385649d77SJunchao Zhang PETSC_i = ic; 86485649d77SJunchao Zhang #else 86585649d77SJunchao Zhang PETSC_i = _Complex_I; 86685649d77SJunchao Zhang #endif 86785649d77SJunchao Zhang } 86885649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */ 86985649d77SJunchao Zhang 87085649d77SJunchao Zhang /* 87185649d77SJunchao Zhang Create the PETSc MPI reduction operator that sums of the first 87285649d77SJunchao Zhang half of the entries and maxes the second half. 87385649d77SJunchao Zhang */ 8749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP)); 87585649d77SJunchao Zhang 876613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 8779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128)); 8789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128)); 8799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128)); 8809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128)); 88185649d77SJunchao Zhang #endif 8829e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) 8839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_CHAR, &MPIU___FP16)); 8849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FP16)); 88585649d77SJunchao Zhang #endif 88685649d77SJunchao Zhang 88785649d77SJunchao Zhang #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 8889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM)); 889613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMax_Local, 1, &MPIU_MAX)); 890613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscMin_Local, 1, &MPIU_MIN)); 8919e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16) 8929e517322SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FP16___FLOAT128)); 89385649d77SJunchao Zhang #endif 89485649d77SJunchao Zhang 8959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR)); 8969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR)); 89785649d77SJunchao Zhang 89885649d77SJunchao Zhang /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */ 89985649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI) 90085649d77SJunchao Zhang { 9019371c9d4SSatish Balay struct PetscRealInt { 9029371c9d4SSatish Balay PetscReal v; 9039371c9d4SSatish Balay PetscInt i; 9049371c9d4SSatish Balay }; 90585649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 90685649d77SJunchao Zhang MPI_Aint blockOffsets[2] = {offsetof(struct PetscRealInt, v), offsetof(struct PetscRealInt, i)}; 90785649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_REAL, MPIU_INT}, tmpStruct; 90885649d77SJunchao Zhang 9099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 9109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct PetscRealInt), &MPIU_REAL_INT)); 9119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT)); 91385649d77SJunchao Zhang } 91485649d77SJunchao Zhang { 9159371c9d4SSatish Balay struct PetscScalarInt { 9169371c9d4SSatish Balay PetscScalar v; 9179371c9d4SSatish Balay PetscInt i; 9189371c9d4SSatish Balay }; 91985649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 92085649d77SJunchao Zhang MPI_Aint blockOffsets[2] = {offsetof(struct PetscScalarInt, v), offsetof(struct PetscScalarInt, i)}; 92185649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_SCALAR, MPIU_INT}, tmpStruct; 92285649d77SJunchao Zhang 9239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 9249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct PetscScalarInt), &MPIU_SCALAR_INT)); 9259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT)); 92785649d77SJunchao Zhang } 92885649d77SJunchao Zhang #endif 92985649d77SJunchao Zhang 93085649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES) 9319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT)); 9329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2INT)); 93385649d77SJunchao Zhang #endif 9349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT)); 9359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPI_4INT)); 9369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT)); 9379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_4INT)); 93885649d77SJunchao Zhang 93985649d77SJunchao Zhang /* 94085649d77SJunchao Zhang Attributes to be set on PETSc communicators 94185649d77SJunchao Zhang */ 9429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0)); 9439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0)); 9449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0)); 9459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0)); 94685649d77SJunchao Zhang 94785649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN) 9489566063dSJacob Faibussowitsch if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len)); 94985649d77SJunchao Zhang else 95085649d77SJunchao Zhang #endif 9519566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file)); 95285649d77SJunchao Zhang 95385649d77SJunchao Zhang /* call a second time so it can look in the options database */ 9549566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 95585649d77SJunchao Zhang 95685649d77SJunchao Zhang /* 95785649d77SJunchao Zhang Check system options and print help 95885649d77SJunchao Zhang */ 9599566063dSJacob Faibussowitsch PetscCall(PetscOptionsCheckInitial_Private(help)); 96085649d77SJunchao Zhang 961a4af0ceeSJacob Faibussowitsch /* 962a4af0ceeSJacob Faibussowitsch Initialize PetscDevice and PetscDeviceContext 963a4af0ceeSJacob Faibussowitsch 964a4af0ceeSJacob Faibussowitsch Note to any future devs thinking of moving this, proper initialization requires: 965a4af0ceeSJacob Faibussowitsch 1. MPI initialized 966a4af0ceeSJacob Faibussowitsch 2. Options DB initialized 967a4af0ceeSJacob Faibussowitsch 3. Petsc error handling initialized, specifically signal handlers. This expects to set up its own SIGSEV handler via 968a4af0ceeSJacob Faibussowitsch the push/pop interface. 969a4af0ceeSJacob Faibussowitsch */ 970a2158755SJunchao Zhang #if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP) || PetscDefined(HAVE_SYCL)) 9719566063dSJacob Faibussowitsch PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD)); 972a4af0ceeSJacob Faibussowitsch #endif 973a4af0ceeSJacob Faibussowitsch 974a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 975a4af0ceeSJacob Faibussowitsch flg = PETSC_FALSE; 9769566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-log_summary", &flg)); 9779566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg)); 9789566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL)); 979a4af0ceeSJacob Faibussowitsch PetscViennaCLSynchronize = flg; 9809566063dSJacob Faibussowitsch PetscCall(PetscViennaCLInit()); 981a4af0ceeSJacob Faibussowitsch #endif 982a4af0ceeSJacob Faibussowitsch 983a4af0ceeSJacob Faibussowitsch /* 984a4af0ceeSJacob Faibussowitsch Creates the logging data structures; this is enabled even if logging is not turned on 985a4af0ceeSJacob Faibussowitsch This is the last thing we do before returning to the user code to prevent having the 986a4af0ceeSJacob Faibussowitsch logging numbers contaminated by any startup time associated with MPI 987a4af0ceeSJacob Faibussowitsch */ 988a4af0ceeSJacob Faibussowitsch #if defined(PETSC_USE_LOG) 9899566063dSJacob Faibussowitsch PetscCall(PetscLogInitialize()); 990a4af0ceeSJacob Faibussowitsch #endif 991a4af0ceeSJacob Faibussowitsch 9929566063dSJacob Faibussowitsch PetscCall(PetscCitationsInitialize()); 99385649d77SJunchao Zhang 99485649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS) 9959566063dSJacob Faibussowitsch PetscCall(PetscInitializeSAWs(ftn ? NULL : help)); 99627104ee2SJacob Faibussowitsch flg = PETSC_FALSE; 9979566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg)); 9989566063dSJacob Faibussowitsch if (flg) PetscCall(PetscStackViewSAWs()); 99985649d77SJunchao Zhang #endif 100085649d77SJunchao Zhang 100185649d77SJunchao Zhang /* 100285649d77SJunchao Zhang Load the dynamic libraries (on machines that support them), this registers all 100385649d77SJunchao Zhang the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 100485649d77SJunchao Zhang */ 10059566063dSJacob Faibussowitsch PetscCall(PetscInitialize_DynamicLibraries()); 100685649d77SJunchao Zhang 10079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 10089566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size)); 10099566063dSJacob Faibussowitsch PetscCall(PetscGetHostName(hostname, 256)); 10109566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname)); 101185649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP) 101285649d77SJunchao Zhang { 101385649d77SJunchao Zhang PetscBool omp_view_flag; 101485649d77SJunchao Zhang char *threads = getenv("OMP_NUM_THREADS"); 101585649d77SJunchao Zhang 101685649d77SJunchao Zhang if (threads) { 10179566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads)); 101885649d77SJunchao Zhang (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads); 101985649d77SJunchao Zhang } else { 10202f840973SStefano Zampini PetscNumOMPThreads = (PetscInt)omp_get_max_threads(); 10219566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads)); 102285649d77SJunchao Zhang } 1023d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys"); 10249566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg)); 10259566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag)); 1026d0609cedSBarry Smith PetscOptionsEnd(); 102785649d77SJunchao Zhang if (flg) { 10289566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads)); 102985649d77SJunchao Zhang omp_set_num_threads((int)PetscNumOMPThreads); 103085649d77SJunchao Zhang } 103148a46eb9SPierre Jolivet if (omp_view_flag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads)); 103285649d77SJunchao Zhang } 103385649d77SJunchao Zhang #endif 103485649d77SJunchao Zhang 103585649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 103685649d77SJunchao Zhang /* 103785649d77SJunchao Zhang Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 103885649d77SJunchao Zhang 103985649d77SJunchao Zhang Currently not used because it is not supported by MPICH. 104085649d77SJunchao Zhang */ 10419566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL)); 104285649d77SJunchao Zhang #endif 104385649d77SJunchao Zhang 104485649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS) 10459566063dSJacob Faibussowitsch PetscCall(PetscFPTCreate(10000)); 104685649d77SJunchao Zhang #endif 104785649d77SJunchao Zhang 104885649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC) 104985649d77SJunchao Zhang { 105085649d77SJunchao Zhang PetscViewer viewer; 10519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg)); 105285649d77SJunchao Zhang if (flg) { 10539566063dSJacob Faibussowitsch PetscCall(PetscProcessPlacementView(viewer)); 10549566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 105585649d77SJunchao Zhang } 105685649d77SJunchao Zhang } 105785649d77SJunchao Zhang #endif 105885649d77SJunchao Zhang 105985649d77SJunchao Zhang flg = PETSC_TRUE; 10609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL)); 10619566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE)); 106285649d77SJunchao Zhang 106385649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS) 10649566063dSJacob Faibussowitsch PetscCall(adios_init_noxml(PETSC_COMM_WORLD)); 10659566063dSJacob Faibussowitsch PetscCall(adios_declare_group(&Petsc_adios_group, "PETSc", "", adios_stat_default)); 10669566063dSJacob Faibussowitsch PetscCall(adios_select_method(Petsc_adios_group, "MPI", "", "")); 10679566063dSJacob Faibussowitsch PetscCall(adios_read_init_method(ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, "")); 106885649d77SJunchao Zhang #endif 106985649d77SJunchao Zhang 107085649d77SJunchao Zhang #if defined(__VALGRIND_H) 107185649d77SJunchao Zhang PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE; 107285649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE) 10739566063dSJacob Faibussowitsch if (PETSC_RUNNING_ON_VALGRIND) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING: Running valgrind with the MacOS native BLAS and LAPACK can fail. If it fails suggest configuring with --download-fblaslapack or --download-f2cblaslapack")); 107485649d77SJunchao Zhang #endif 107585649d77SJunchao Zhang #endif 107685649d77SJunchao Zhang /* 107785649d77SJunchao Zhang Set flag that we are completely initialized 107885649d77SJunchao Zhang */ 107985649d77SJunchao Zhang PetscInitializeCalled = PETSC_TRUE; 108085649d77SJunchao Zhang 10819566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg)); 10829566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPythonInitialize(NULL, NULL)); 1083f1f2ae84SBarry Smith 1084f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1085f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin()); 1086f1f2ae84SBarry 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"); 108727104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 108885649d77SJunchao Zhang } 108985649d77SJunchao Zhang 1090e5c89e4eSSatish Balay /*@C 1091e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 1092*811af0c4SBarry Smith `PetscInitialize()` calls MPI_Init() if that has yet to be called, 1093e5c89e4eSSatish Balay so this routine should always be called near the beginning of 1094e5c89e4eSSatish Balay your program -- usually the very first line! 1095e5c89e4eSSatish Balay 1096*811af0c4SBarry Smith Collective on `MPI_COMM_WORLD` or `PETSC_COMM_WORLD` if it has been set 1097e5c89e4eSSatish Balay 1098e5c89e4eSSatish Balay Input Parameters: 1099e5c89e4eSSatish Balay + argc - count of number of command line arguments 1100e5c89e4eSSatish Balay . args - the command line arguments 1101be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format. 1102be10d61cSLisandro Dalcin Use NULL or empty string to not check for code specific file. 1103be10d61cSLisandro Dalcin Also checks ~/.petscrc, .petscrc and petscrc. 1104c5b5d8d5SVaclav Hapla Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files. 11050298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message 1106e5c89e4eSSatish Balay 1107*811af0c4SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of `MPI_COMM_WORLD`, create that 1108*811af0c4SBarry Smith communicator first and assign it to `PETSC_COMM_WORLD` BEFORE calling `PetscInitialize()`. Thus if you are running a 1109*811af0c4SBarry Smith four process job and two processes will run PETSc and have `PetscInitialize()` and PetscFinalize() and two process will not, 1110*811af0c4SBarry Smith then do this. If ALL processes in the job are using `PetscInitialize()` and `PetscFinalize()` then you don't need to do this, even 111105827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 1112e5c89e4eSSatish Balay 1113e5c89e4eSSatish Balay Options Database Keys: 11147ca660e7SBarry Smith + -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message 11157ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 1116e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 11178a690491SBarry Smith . -on_error_emacs <machinename> - causes emacsclient to jump to error file 1118*811af0c4SBarry Smith . -on_error_abort - calls `abort()` when error detected (no traceback) 1119*811af0c4SBarry Smith . -on_error_mpiabort - calls `MPI_abort()` when error detected 11201af3601dSBarry Smith . -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr 11218a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 1122bf4d2887SBarry Smith . -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger 1123e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds) - Pauses debugger 1124e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 1125e5c89e4eSSatish Balay process and wait (-debugger_pause) seconds for attachment 112679dccf82SBarry Smith . -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug) 112779dccf82SBarry Smith . -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no) 1128*811af0c4SBarry Smith . -malloc_debug - check for memory corruption at EVERY malloc or free, see `PetscMallocSetDebug()` 1129aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 113092f119d6SBarry 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 1131*811af0c4SBarry Smith . -malloc_view - show a list of all allocated memory during `PetscFinalize()` 113292f119d6SBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view 1133608c71bfSMatthew G. Knepley . -malloc_requested_size - malloc logging will record the requested size rather than size after alignment 113492f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 1135e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 1136e5c89e4eSSatish Balay . -shared_tmp - indicates /tmp directory is shared by all processors 1137e5c89e4eSSatish Balay . -not_shared_tmp - each processor has own /tmp 1138e5c89e4eSSatish Balay . -tmp - alternative name of /tmp directory 1139e5c89e4eSSatish Balay . -get_total_flops - returns total flops done by all processors 11400841954dSBarry Smith - -memory_view - Print memory usage at end of run 1141e5c89e4eSSatish Balay 1142c5b5d8d5SVaclav Hapla Options Database Keys for Option Database: 1143c5b5d8d5SVaclav Hapla + -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc 1144c5b5d8d5SVaclav Hapla . -options_monitor - monitor all set options to standard output for the whole program run 1145*811af0c4SBarry Smith - -options_monitor_cancel - cancel options monitoring hard-wired using `PetscOptionsMonitorSet()` 1146c5b5d8d5SVaclav Hapla 1147c5b5d8d5SVaclav Hapla Options -options_monitor_{all,cancel} are 1148c5b5d8d5SVaclav Hapla position-independent and apply to all options set since the PETSc start. 1149c5b5d8d5SVaclav Hapla They can be used also in option files. 1150c5b5d8d5SVaclav Hapla 1151*811af0c4SBarry Smith See `PetscOptionsMonitorSet()` to do monitoring programmatically. 1152c5b5d8d5SVaclav Hapla 1153e5c89e4eSSatish Balay Options Database Keys for Profiling: 1154a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 1155*811af0c4SBarry Smith + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See `PetscInfo()`. 1156217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 1157217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 1158495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 1159*811af0c4SBarry Smith hangs without running in the debugger). See `PetscLogTraceBegin()`. 1160*811af0c4SBarry Smith . -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see `PetscLogView()`. 1161*811af0c4SBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each event, see `PetscLogView()`. 1162*811af0c4SBarry Smith . -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see `PetscLogView(). 11639a9a5d4cSBarry Smith . -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the 1164495fc317SBarry Smith summary is written to the file. See PetscLogView(). 1165f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 1166*811af0c4SBarry Smith . -log_all [filename] - Logs extensive profiling information See `PetscLogDump()`. 1167*811af0c4SBarry Smith . -log [filename] - Logs basic profiline information See `PetscLogDump()`. 1168c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 1169*811af0c4SBarry Smith . -viewfromoptions on,off - Enable or disable `XXXSetFromOptions()` calls, for applications with many small solves turn this off 1170c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 1171495fc317SBarry Smith 11726a77f485SPierre Jolivet Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time 1173e5c89e4eSSatish Balay 1174ffbd1cfbSBarry Smith Options Database Keys for SAWs: 1175ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 1176ffbd1cfbSBarry 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 1177ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 1178ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 1179ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 1180ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 1181ffbd1cfbSBarry Smith 1182e5c89e4eSSatish Balay Environmental Variables: 1183*811af0c4SBarry Smith + `PETSC_TMP` - alternative tmp directory 1184*811af0c4SBarry Smith . `PETSC_SHARED_TMP` - tmp is shared by all processes 1185*811af0c4SBarry Smith . `PETSC_NOT_SHARED_TMP` - each process has its own private tmp 1186*811af0c4SBarry Smith . `PETSC_OPTIONS` - a string containing additional options for petsc in the form of command line "-key value" pairs 1187*811af0c4SBarry Smith . `PETSC_OPTIONS_YAML` - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document 1188*811af0c4SBarry Smith . `PETSC_VIEWER_SOCKET_PORT` - socket number to use for socket viewer 1189*811af0c4SBarry Smith - `PETSC_VIEWER_SOCKET_MACHINE` - machine to use for socket viewer to connect to 1190e5c89e4eSSatish Balay 1191e5c89e4eSSatish Balay Level: beginner 1192e5c89e4eSSatish Balay 1193*811af0c4SBarry Smith Note: 1194*811af0c4SBarry Smith If for some reason you must call `MPI_Init()` separately, call 1195*811af0c4SBarry Smith it before `PetscInitialize()`. 1196e5c89e4eSSatish Balay 1197*811af0c4SBarry Smith Fortran Notes: 1198*811af0c4SBarry Smith In Fortran this routine can be called with 1199*811af0c4SBarry Smith .vb 1200*811af0c4SBarry Smith call PetscInitialize(ierr) 1201*811af0c4SBarry Smith call PetscInitialize(file,ierr) or 1202*811af0c4SBarry Smith call PetscInitialize(file,help,ierr) 1203*811af0c4SBarry Smith .ve 1204e5c89e4eSSatish Balay 1205*811af0c4SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call `PetscInitializeFortran()` soon after 1206*811af0c4SBarry Smith calling `PetscInitialize()`. 1207e5c89e4eSSatish Balay 1208db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()` 1209e5c89e4eSSatish Balay @*/ 12109371c9d4SSatish Balay PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[]) { 121185649d77SJunchao Zhang PetscMPIInt flag; 121221ef0414SBarry Smith const char *prog = "Unknown Name", *mpienv; 1213e5c89e4eSSatish Balay 121427104ee2SJacob Faibussowitsch PetscFunctionBegin; 121527104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(0); 12169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Initialized(&flag)); 1217e5c89e4eSSatish Balay if (!flag) { 121839a651e2SJacob 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"); 12199566063dSJacob Faibussowitsch PetscCall(PetscPreMPIInit_Private()); 12205e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 12215e765c61SJed Brown { 122239a651e2SJacob Faibussowitsch PetscMPIInt PETSC_UNUSED provided; 12239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED, &provided)); 12245e765c61SJed Brown } 12255e765c61SJed Brown #else 12269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init(argc, args)); 12275e765c61SJed Brown #endif 122821ef0414SBarry Smith if (PetscDefined(HAVE_MPIUNI)) { 122921ef0414SBarry Smith mpienv = getenv("PMI_SIZE"); 123021ef0414SBarry Smith if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE"); 123121ef0414SBarry Smith if (mpienv) { 123221ef0414SBarry Smith PetscInt isize; 123321ef0414SBarry Smith PetscCall(PetscOptionsStringToInt(mpienv, &isize)); 123421ef0414SBarry 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"); 123521ef0414SBarry 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"); 123621ef0414SBarry Smith } 123721ef0414SBarry Smith } 1238e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 1239e5c89e4eSSatish Balay } 12404dfee713SSatish Balay 124185649d77SJunchao Zhang if (argc && *argc) prog = **args; 1242e5c89e4eSSatish Balay if (argc && args) { 1243e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 1244e5c89e4eSSatish Balay PetscGlobalArgs = *args; 1245e5c89e4eSSatish Balay } 1246*811af0c4SBarry Smith PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE, PETSC_FALSE, 0)); 124727104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 1248e5c89e4eSSatish Balay } 1249e5c89e4eSSatish Balay 1250a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 125195c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1252ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1253ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 125495c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 12554097062eSBarry Smith #endif 1256e5c89e4eSSatish Balay 1257008a6e76SBarry Smith /* 1258008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1259008a6e76SBarry Smith */ 12609371c9d4SSatish Balay PetscErrorCode PetscFreeMPIResources(void) { 1261008a6e76SBarry Smith PetscFunctionBegin; 1262613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 12639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128)); 12649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128)); 1265008a6e76SBarry Smith #endif 12669e517322SPierre Jolivet #if defined(PETSC_HAVE_REAL___FP16) 12679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FP16)); 1268008a6e76SBarry Smith #endif 1269008a6e76SBarry Smith 1270de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 12719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_SUM)); 1272613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MAX)); 1273613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MIN)); 12749e517322SPierre Jolivet #elif defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16) 12759e517322SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_SUM___FP16___FLOAT128)); 1276008a6e76SBarry Smith #endif 1277008a6e76SBarry Smith 12789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR)); 12799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT)); 12809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT)); 128140df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 12829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2INT)); 1283008a6e76SBarry Smith #endif 12849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPI_4INT)); 12859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_4INT)); 12869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP)); 1287008a6e76SBarry Smith PetscFunctionReturn(0); 1288008a6e76SBarry Smith } 1289008a6e76SBarry Smith 1290a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 1291a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void); 1292a4af0ceeSJacob Faibussowitsch #endif 1293a4af0ceeSJacob Faibussowitsch 1294e5c89e4eSSatish Balay /*@C 1295e5c89e4eSSatish Balay PetscFinalize - Checks for options to be called at the conclusion 1296*811af0c4SBarry Smith of the program. `MPI_Finalize()` is called only if the user had not 1297*811af0c4SBarry Smith called `MPI_Init()` before calling `PetscInitialize()`. 1298e5c89e4eSSatish Balay 1299*811af0c4SBarry Smith Collective on `PETSC_COMM_WORLD` 1300e5c89e4eSSatish Balay 1301e5c89e4eSSatish Balay Options Database Keys: 1302*811af0c4SBarry Smith + -options_view - Calls `PetscOptionsView()` 1303e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 13047eb1d149SBarry 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 1305e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 1306*811af0c4SBarry Smith . -malloc_dump <optional filename> - Calls `PetscMallocDump()`, displays all memory allocated that has not been freed 1307e5c89e4eSSatish Balay . -malloc_info - Prints total memory usage 130892f119d6SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and where 1309e5c89e4eSSatish Balay 1310e5c89e4eSSatish Balay Level: beginner 1311e5c89e4eSSatish Balay 1312e5c89e4eSSatish Balay Note: 1313*811af0c4SBarry Smith See `PetscInitialize()` for other runtime options. 1314e5c89e4eSSatish Balay 1315db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()` 1316e5c89e4eSSatish Balay @*/ 13179371c9d4SSatish Balay PetscErrorCode PetscFinalize(void) { 13184bb5149bSJed Brown PetscMPIInt rank; 1319a8d2bbe5SBarry Smith PetscInt nopt; 13202bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE; 1321dff31646SBarry Smith PetscBool flg; 132210463e74SBarry Smith #if defined(PETSC_USE_LOG) 132310463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 132410463e74SBarry Smith #endif 1325e5c89e4eSSatish Balay 13263cbbc5ffSBarry Smith PetscFunctionBegin; 132739a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()"); 13289566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PetscFinalize() called\n")); 1329b022a5c1SBarry Smith 1330f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1331f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd()); 1332f1f2ae84SBarry Smith 13339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1334a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 13359566063dSJacob Faibussowitsch PetscCall(adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE)); 13369566063dSJacob Faibussowitsch PetscCall(adios_finalize(rank)); 1337a56f64adSBarry Smith #endif 13389566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg)); 1339dff31646SBarry Smith if (flg) { 13401f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 13411f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 13421f817a21SBarry Smith 13439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL)); 134448a46eb9SPierre Jolivet if (filename[0]) PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd)); 13459566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits)); 1346dff31646SBarry Smith cits[0] = 0; 13479566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits)); 13489566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n")); 13499566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 13509566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits)); 13519566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 13529566063dSJacob Faibussowitsch PetscCall(PetscFClose(PETSC_COMM_WORLD, fd)); 13539566063dSJacob Faibussowitsch PetscCall(PetscFree(cits)); 1354dff31646SBarry Smith } 13559566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&PetscCitationsList)); 1356dff31646SBarry Smith 1357c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER) 135804102261SBarry Smith /* TextBelt is run for testing purposes only, please do not use this feature often */ 135904102261SBarry Smith { 136004102261SBarry Smith PetscInt nmax = 2; 136104102261SBarry Smith char **buffs; 13629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &buffs)); 13639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-textbelt", buffs, &nmax, &flg1)); 136404102261SBarry Smith if (flg1) { 136528b400f6SJacob Faibussowitsch PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-textbelt requires either the phone number or number,\"message\""); 136604102261SBarry Smith if (nmax == 1) { 13679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(128, &buffs[1])); 13689566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(buffs[1], 32)); 13699566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buffs[1], " has completed")); 137004102261SBarry Smith } 13719566063dSJacob Faibussowitsch PetscCall(PetscTextBelt(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL)); 13729566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[0])); 13739566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[1])); 137404102261SBarry Smith } 13759566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs)); 137604102261SBarry Smith } 1377763ec7b1SBarry Smith { 1378763ec7b1SBarry Smith PetscInt nmax = 2; 1379763ec7b1SBarry Smith char **buffs; 13809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &buffs)); 13819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-tellmycell", buffs, &nmax, &flg1)); 1382763ec7b1SBarry Smith if (flg1) { 138328b400f6SJacob Faibussowitsch PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-tellmycell requires either the phone number or number,\"message\""); 1384763ec7b1SBarry Smith if (nmax == 1) { 13859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(128, &buffs[1])); 13869566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(buffs[1], 32)); 13879566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buffs[1], " has completed")); 1388763ec7b1SBarry Smith } 13899566063dSJacob Faibussowitsch PetscCall(PetscTellMyCell(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL)); 13909566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[0])); 13919566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[1])); 1392763ec7b1SBarry Smith } 13939566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs)); 1394763ec7b1SBarry Smith } 139504102261SBarry Smith #endif 139604102261SBarry Smith 13972d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 13989566063dSJacob Faibussowitsch PetscCall(PetscFPTDestroy()); 13992d53ad75SBarry Smith #endif 14002d53ad75SBarry Smith 1401e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1402dff31646SBarry Smith flg = PETSC_FALSE; 14039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL)); 14041baa6e33SBarry Smith if (flg) PetscCall(PetscOptionsSAWsDestroy()); 1405d5649816SBarry Smith #endif 1406d5649816SBarry Smith 1407681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1408681455b2SBarry Smith flg1 = PETSC_FALSE; 14099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL)); 1410681455b2SBarry Smith if (flg1) { 1411681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 14129566063dSJacob Faibussowitsch PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL)); 1413681455b2SBarry Smith } 1414681455b2SBarry Smith #endif 1415681455b2SBarry Smith 141667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 14179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_info", &flg2, NULL)); 1418e5c89e4eSSatish Balay if (!flg2) { 141990d69ab7SBarry Smith flg2 = PETSC_FALSE; 14209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL)); 1421e5c89e4eSSatish Balay } 142248a46eb9SPierre Jolivet if (flg2) PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n")); 142367584ceeSBarry Smith #endif 1424e5c89e4eSSatish Balay 1425e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 142690d69ab7SBarry Smith flg1 = PETSC_FALSE; 14279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL)); 1428e5c89e4eSSatish Balay if (flg1) { 1429e5c89e4eSSatish Balay PetscLogDouble flops = 0; 14309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD)); 14319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops)); 1432e5c89e4eSSatish Balay } 1433e5c89e4eSSatish Balay #endif 1434e5c89e4eSSatish Balay 1435e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 1436e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE) 1437e5c89e4eSSatish Balay mname[0] = 0; 14389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1)); 1439e5c89e4eSSatish Balay if (flg1) { 14409566063dSJacob Faibussowitsch if (mname[0]) PetscCall(PetscLogMPEDump(mname)); 14419566063dSJacob Faibussowitsch else PetscCall(PetscLogMPEDump(0)); 1442e5c89e4eSSatish Balay } 1443e5c89e4eSSatish Balay #endif 1444356e5837SBarry Smith #endif 1445a297a907SKarl Rupp 1446dd710f27SBarry Smith /* 1447dd710f27SBarry Smith Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 1448dd710f27SBarry Smith */ 14499566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1450dd710f27SBarry Smith 1451356e5837SBarry Smith #if defined(PETSC_USE_LOG) 14529566063dSJacob Faibussowitsch PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE)); 14539566063dSJacob Faibussowitsch PetscCall(PetscLogViewFromOptions()); 14549566063dSJacob Faibussowitsch PetscCall(PetscOptionsPopGetViewerOff()); 145587c3beb6SLisandro Dalcin 1456356e5837SBarry Smith mname[0] = 0; 14579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_summary", mname, sizeof(mname), &flg1)); 1458e5c89e4eSSatish Balay if (flg1) { 145991eabc43SBarry Smith PetscViewer viewer; 14609566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD, "\n\n WARNING: -log_summary is being deprecated; switch to -log_view\n\n\n")); 146191eabc43SBarry Smith if (mname[0]) { 14629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, mname, &viewer)); 14639566063dSJacob Faibussowitsch PetscCall(PetscLogView(viewer)); 14649566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 146533f85c2fSBarry Smith } else { 146633f85c2fSBarry Smith viewer = PETSC_VIEWER_STDOUT_WORLD; 14679566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT)); 14689566063dSJacob Faibussowitsch PetscCall(PetscLogView(viewer)); 14699566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 147033f85c2fSBarry Smith } 1471e5c89e4eSSatish Balay } 1472a297a907SKarl Rupp 1473dd710f27SBarry Smith /* 1474dd710f27SBarry Smith Free any objects created by the last block of code. 1475dd710f27SBarry Smith */ 14769566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1477dd710f27SBarry Smith 1478dd710f27SBarry Smith mname[0] = 0; 14799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1)); 14809566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2)); 14819566063dSJacob Faibussowitsch if (flg1 || flg2) PetscCall(PetscLogDump(mname)); 1482e5c89e4eSSatish Balay #endif 148310463e74SBarry Smith 148490d69ab7SBarry Smith flg1 = PETSC_FALSE; 14859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL)); 14869566063dSJacob Faibussowitsch if (!flg1) PetscCall(PetscPopSignalHandler()); 148790d69ab7SBarry Smith flg1 = PETSC_FALSE; 14889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL)); 14891baa6e33SBarry Smith if (flg1) PetscCall(PetscMPIDump(stdout)); 149090d69ab7SBarry Smith flg1 = PETSC_FALSE; 149190d69ab7SBarry Smith flg2 = PETSC_FALSE; 14928bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 14939566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1)); 14949566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 14959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL)); 1496e4c476e2SSatish Balay 1497e5c89e4eSSatish Balay if (flg2) { 1498be56827dSJed Brown PetscViewer viewer; 14999566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 15009566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 15019566063dSJacob Faibussowitsch PetscCall(PetscOptionsView(NULL, viewer)); 15029566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1503e5c89e4eSSatish Balay } 1504e5c89e4eSSatish Balay 1505e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 15069566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1)); 15079566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1)); 1508e5c89e4eSSatish Balay 150933fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 15109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1)); 15119245e749SBarry Smith if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE; 1512e5c89e4eSSatish Balay if (flg3) { 15133de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 1514be56827dSJed Brown PetscViewer viewer; 15159566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 15169566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 15179566063dSJacob Faibussowitsch PetscCall(PetscOptionsView(NULL, viewer)); 15189566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1519e5c89e4eSSatish Balay } 15209566063dSJacob Faibussowitsch PetscCall(PetscOptionsAllUsed(NULL, &nopt)); 15213de2bfdfSBarry Smith if (nopt) { 15229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n")); 15239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n")); 15243de2bfdfSBarry Smith if (nopt == 1) { 15259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n")); 1526e5c89e4eSSatish Balay } else { 15279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt)); 1528e5c89e4eSSatish Balay } 15293de2bfdfSBarry Smith } else if (flg3 && flg1) { 15309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n")); 1531df12ba86SBarry Smith } 15329566063dSJacob Faibussowitsch PetscCall(PetscOptionsLeft(NULL)); 1533e5c89e4eSSatish Balay } 1534e5c89e4eSSatish Balay 1535e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1536d45a07a7SBarry Smith if (!PetscGlobalRank) { 15379566063dSJacob Faibussowitsch PetscCall(PetscStackSAWsViewOff()); 1538792fecdfSBarry Smith PetscCallSAWs(SAWs_Finalize, ()); 1539d45a07a7SBarry Smith } 1540ec957eceSBarry Smith #endif 1541ec957eceSBarry Smith 15424097062eSBarry Smith #if defined(PETSC_USE_LOG) 154310463e74SBarry Smith /* 1544dbc8283eSBarry Smith List all objects the user may have forgot to free 15452eff7a51SBarry Smith */ 154605df10baSBarry Smith if (PetscObjectsLog) { 15479566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 1548a64a8e02SBarry Smith if (flg1) { 1549a64a8e02SBarry Smith MPI_Comm local_comm; 15507eb1d149SBarry Smith char string[64]; 1551a64a8e02SBarry Smith 15529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL)); 1553afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 15549566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 15559566063dSJacob Faibussowitsch PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE)); 15569566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 15579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 15580a1571b3SBarry Smith } 155905df10baSBarry Smith } 15604097062eSBarry Smith #endif 15614097062eSBarry Smith 15624097062eSBarry Smith #if defined(PETSC_USE_LOG) 1563dbc8283eSBarry Smith PetscObjectsCounts = 0; 1564dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 15659566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscObjects)); 15664097062eSBarry Smith #endif 15672eff7a51SBarry Smith 156833f85c2fSBarry Smith /* 156933f85c2fSBarry Smith Destroy any packages that registered a finalize 157033f85c2fSBarry Smith */ 15719566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalizeAll()); 157233f85c2fSBarry Smith 1573101409b8SToby Isaac #if defined(PETSC_USE_LOG) 15749566063dSJacob Faibussowitsch PetscCall(PetscLogFinalize()); 1575101409b8SToby Isaac #endif 1576101409b8SToby Isaac 157733f85c2fSBarry Smith /* 157848dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 157948dd1dffSBarry Smith */ 15802e956fe4SStefano Zampini if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll()); 158137e93019SBarry Smith 15824028d114SSatish Balay if (petsc_history) { 15839566063dSJacob Faibussowitsch PetscCall(PetscCloseHistoryFile(&petsc_history)); 158402c9f0b5SLisandro Dalcin petsc_history = NULL; 1585e5c89e4eSSatish Balay } 15869566063dSJacob Faibussowitsch PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton)); 15879566063dSJacob Faibussowitsch PetscCall(PetscInfoDestroy()); 1588e5c89e4eSSatish Balay 158967584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 159092f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1591e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 159292f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1593e5c89e4eSSatish Balay FILE *fd; 1594ed9cf6e9SBarry Smith int err; 1595e5c89e4eSSatish Balay 1596dc92acbaSJed Brown flg2 = PETSC_FALSE; 159792f119d6SBarry Smith flg3 = PETSC_FALSE; 15989566063dSJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL)); 15999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL)); 160092f119d6SBarry Smith fname[0] = 0; 16019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1)); 1602e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1603589a23caSBarry Smith PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank); 16049371c9d4SSatish Balay fd = fopen(sname, "w"); 16059371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16069566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(fd)); 1607ed9cf6e9SBarry Smith err = fclose(fd); 160828b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 160992f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1610e5c89e4eSSatish Balay MPI_Comm local_comm; 1611e5c89e4eSSatish Balay 1612afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16139566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16149566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(stdout)); 16159566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1617e5c89e4eSSatish Balay } 1618e5c89e4eSSatish Balay fname[0] = 0; 16199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1)); 1620e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1621589a23caSBarry Smith PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank); 16229371c9d4SSatish Balay fd = fopen(sname, "w"); 16239371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16249566063dSJacob Faibussowitsch PetscCall(PetscMallocView(fd)); 1625ed9cf6e9SBarry Smith err = fclose(fd); 162628b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 162792f119d6SBarry Smith } else if (flg1) { 162892f119d6SBarry Smith MPI_Comm local_comm; 162992f119d6SBarry Smith 1630afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16319566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16329566063dSJacob Faibussowitsch PetscCall(PetscMallocView(stdout)); 16339566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1635e5c89e4eSSatish Balay } 1636e5c89e4eSSatish Balay } 163767584ceeSBarry Smith #endif 163820e2c332SMatthew G. Knepley 16395486ca60SMatthew G. Knepley /* 16405486ca60SMatthew G. Knepley Close any open dynamic libraries 16415486ca60SMatthew G. Knepley */ 16429566063dSJacob Faibussowitsch PetscCall(PetscFinalize_DynamicLibraries()); 16435486ca60SMatthew G. Knepley 1644e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 16459566063dSJacob Faibussowitsch PetscCall(PetscOptionsDestroyDefault()); 1646e5c89e4eSSatish Balay 1647e5c89e4eSSatish Balay PetscGlobalArgc = 0; 164802c9f0b5SLisandro Dalcin PetscGlobalArgs = NULL; 1649e5c89e4eSSatish Balay 1650c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 1651c2b86a48SJunchao Zhang if (PetscBeganKokkos) { 16529566063dSJacob Faibussowitsch PetscCall(PetscKokkosFinalize_Private()); 1653c2b86a48SJunchao Zhang PetscBeganKokkos = PETSC_FALSE; 165445639126SStefano Zampini PetscKokkosInitialized = PETSC_FALSE; 1655c2b86a48SJunchao Zhang } 1656c2b86a48SJunchao Zhang #endif 1657c2b86a48SJunchao Zhang 165871438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM) 165971438e86SJunchao Zhang if (PetscBeganNvshmem) { 16609566063dSJacob Faibussowitsch PetscCall(PetscNvshmemFinalize()); 166171438e86SJunchao Zhang PetscBeganNvshmem = PETSC_FALSE; 166271438e86SJunchao Zhang } 166371438e86SJunchao Zhang #endif 1664a0e72f99SJunchao Zhang 16659566063dSJacob Faibussowitsch PetscCall(PetscFreeMPIResources()); 1666e5c89e4eSSatish Balay 1667dbc8283eSBarry Smith /* 1668efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1669efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1670efb80d3cSBarry Smith 1671efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1672efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1673dbc8283eSBarry Smith */ 1674b770b1f6SSatish Balay { 1675dbc8283eSBarry Smith PetscCommCounter *counter; 1676dbc8283eSBarry Smith PetscMPIInt flg; 1677dbc8283eSBarry Smith MPI_Comm icomm; 16789371c9d4SSatish Balay union 16799371c9d4SSatish Balay { 16809371c9d4SSatish Balay MPI_Comm comm; 16819371c9d4SSatish Balay void *ptr; 16829371c9d4SSatish Balay } ucomm; 16839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 1684dbc8283eSBarry Smith if (flg) { 1685265f3f35SJed Brown icomm = ucomm.comm; 16869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 168728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1688dbc8283eSBarry Smith 16899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval)); 16909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 16919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1692dbc8283eSBarry Smith } 16939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 1694dbc8283eSBarry Smith if (flg) { 1695265f3f35SJed Brown icomm = ucomm.comm; 16969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 169728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1698dbc8283eSBarry Smith 16999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval)); 17009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 17019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1702dbc8283eSBarry Smith } 1703b770b1f6SSatish Balay } 1704dbc8283eSBarry Smith 17059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval)); 17069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval)); 17079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval)); 17089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval)); 1709480cf27aSJed Brown 17109566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen)); 17119566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout)); 17129566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr)); 17139566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock)); 1714ef19f930SBarry Smith 1715e5c89e4eSSatish Balay if (PetscBeganMPI) { 171699b1327fSBarry Smith PetscMPIInt flag; 17179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalized(&flag)); 171839a651e2SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 171939a651e2SJacob Faibussowitsch /* wait until the very last moment to disable error handling */ 172039a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_FALSE; 17219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalize()); 172239a651e2SJacob Faibussowitsch } else PetscErrorHandlingInitialized = PETSC_FALSE; 172339a651e2SJacob Faibussowitsch 1724e5c89e4eSSatish Balay /* 1725e5c89e4eSSatish Balay 1726e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1727e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1728e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1729e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1730e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 17310e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1732e5c89e4eSSatish Balay memory was not freed. 1733e5c89e4eSSatish Balay 1734e5c89e4eSSatish Balay */ 17359566063dSJacob Faibussowitsch PetscCall(PetscMallocClear()); 17369566063dSJacob Faibussowitsch PetscCall(PetscStackReset()); 1737a297a907SKarl Rupp 1738e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1739e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 174056883f60SBarry Smith #if defined(PETSC_USE_GCOV) 174156883f60SBarry Smith /* 174256883f60SBarry Smith flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the 174356883f60SBarry Smith gcov files are still being added to the directories as git tries to remove the directories. 174456883f60SBarry Smith */ 174556883f60SBarry Smith __gcov_flush(); 174656883f60SBarry Smith #endif 17471724198aSStefano Zampini /* To match PetscFunctionBegin() at the beginning of this function */ 17481724198aSStefano Zampini PetscStackClearTop; 174927104ee2SJacob Faibussowitsch return 0; 1750e5c89e4eSSatish Balay } 1751e5c89e4eSSatish Balay 175243db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 17539371c9d4SSatish Balay PETSC_EXTERN int lsame_(char *a, char *b) { 175443db4dbbSBarry Smith if (*a == *b) return 1; 175543db4dbbSBarry Smith if (*a + 32 == *b) return 1; 175643db4dbbSBarry Smith if (*a - 32 == *b) return 1; 175743db4dbbSBarry Smith return 0; 175843db4dbbSBarry Smith } 1759a70650f6SBarry Smith #endif 176043db4dbbSBarry Smith 176143db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 17629371c9d4SSatish Balay PETSC_EXTERN int lsame(char *a, char *b) { 176343db4dbbSBarry Smith if (*a == *b) return 1; 176443db4dbbSBarry Smith if (*a + 32 == *b) return 1; 176543db4dbbSBarry Smith if (*a - 32 == *b) return 1; 176643db4dbbSBarry Smith return 0; 176743db4dbbSBarry Smith } 176843db4dbbSBarry Smith #endif 1769