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 115e5c89e4eSSatish Balay 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 /*@ 148e5c89e4eSSatish Balay 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 224613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 225613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 226613bf2b2SPierre Jolivet #include <quadmath.h> 227613bf2b2SPierre Jolivet MPI_Op MPIU_SUM___FLOAT128 = 0; 228613bf2b2SPierre Jolivet #endif 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 2567c2de775SJed Brown else { 257613bf2b2SPierre Jolivet #if !defined(PETSC_HAVE_REAL___FLOAT128) 2587c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 259613bf2b2SPierre Jolivet #else 260613bf2b2SPierre Jolivet (*PetscErrorPrintf)("Can only handle MPIU_REAL, MPIU_COMPLEX, MPIU___FLOAT128, or MPIU___COMPLEX128 data types"); 261613bf2b2SPierre Jolivet #endif 26241e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 263e2e03761SBarry Smith } 264812af9f3SBarry Smith PetscFunctionReturnVoid(); 265e5c89e4eSSatish Balay } 266e5c89e4eSSatish Balay #endif 267e5c89e4eSSatish Balay 268570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 269d9822059SBarry Smith MPI_Op MPIU_MAX = 0; 270d9822059SBarry Smith MPI_Op MPIU_MIN = 0; 271d9822059SBarry Smith 2729371c9d4SSatish Balay PETSC_EXTERN void MPIAPI PetscMax_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) { 273d9822059SBarry Smith PetscInt i, count = *cnt; 274d9822059SBarry Smith 275d9822059SBarry Smith PetscFunctionBegin; 2767c2de775SJed Brown if (*datatype == MPIU_REAL) { 2778c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 278a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMax(xout[i], xin[i]); 2797c2de775SJed Brown } 2807c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2817c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2827c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 2839371c9d4SSatish Balay for (i = 0; i < count; i++) { xout[i] = PetscRealPartComplex(xout[i]) < PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; } 2847c2de775SJed Brown } 2857c2de775SJed Brown #endif 2867c2de775SJed Brown else { 2877c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 28841e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 2898c764dc5SJose Roman } 290d9822059SBarry Smith PetscFunctionReturnVoid(); 291d9822059SBarry Smith } 292d9822059SBarry Smith 2939371c9d4SSatish Balay PETSC_EXTERN void MPIAPI PetscMin_Local(void *in, void *out, PetscMPIInt *cnt, MPI_Datatype *datatype) { 294d9822059SBarry Smith PetscInt i, count = *cnt; 295d9822059SBarry Smith 296d9822059SBarry Smith PetscFunctionBegin; 2977c2de775SJed Brown if (*datatype == MPIU_REAL) { 2988c764dc5SJose Roman PetscReal *xin = (PetscReal *)in, *xout = (PetscReal *)out; 299a297a907SKarl Rupp for (i = 0; i < count; i++) xout[i] = PetscMin(xout[i], xin[i]); 3007c2de775SJed Brown } 3017c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 3027c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 3037c2de775SJed Brown PetscComplex *xin = (PetscComplex *)in, *xout = (PetscComplex *)out; 3049371c9d4SSatish Balay for (i = 0; i < count; i++) { xout[i] = PetscRealPartComplex(xout[i]) > PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; } 3057c2de775SJed Brown } 3067c2de775SJed Brown #endif 3077c2de775SJed Brown else { 3088c764dc5SJose Roman (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"); 30941e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF, PETSC_ERR_ARG_WRONG); 3108c764dc5SJose Roman } 311d9822059SBarry Smith PetscFunctionReturnVoid(); 312d9822059SBarry Smith } 313d9822059SBarry Smith #endif 314d9822059SBarry Smith 315480cf27aSJed Brown /* 316480cf27aSJed Brown Private routine to delete internal tag/name counter storage when a communicator is freed. 317480cf27aSJed Brown 318ff0e51ddSBarry 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. 319480cf27aSJed Brown 32012801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 321480cf27aSJed Brown 322480cf27aSJed Brown */ 3239371c9d4SSatish Balay PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *count_val, void *extra_state) { 32405342407SJunchao Zhang PetscCommCounter *counter = (PetscCommCounter *)count_val; 32557f21012SBarry Smith struct PetscCommStash *comms = counter->comms, *pcomm; 326480cf27aSJed Brown 327480cf27aSJed Brown PetscFunctionBegin; 3289566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "Deleting counter data in an MPI_Comm %ld\n", (long)comm)); 3299566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter->iflags)); 33057f21012SBarry Smith while (comms) { 3319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&comms->comm)); 33257f21012SBarry Smith pcomm = comms; 33357f21012SBarry Smith comms = comms->next; 3349566063dSJacob Faibussowitsch PetscCall(PetscFree(pcomm)); 33557f21012SBarry Smith } 3369566063dSJacob Faibussowitsch PetscCallMPI(PetscFree(counter)); 337480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 338480cf27aSJed Brown } 339480cf27aSJed Brown 340480cf27aSJed Brown /* 34147435625SJed Brown This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user 342da3039f7SJed Brown calls MPI_Comm_free(). 343da3039f7SJed Brown 344da3039f7SJed Brown This is the only entry point for breaking the links between inner and outer comms. 345480cf27aSJed Brown 346ff0e51ddSBarry Smith This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator. 347480cf27aSJed Brown 34812801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 349480cf27aSJed Brown 350480cf27aSJed Brown */ 3519371c9d4SSatish Balay PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) { 3529371c9d4SSatish Balay union 353480cf27aSJed Brown { 3549371c9d4SSatish Balay MPI_Comm comm; 3559371c9d4SSatish Balay void *ptr; 3569371c9d4SSatish Balay } icomm; 357480cf27aSJed Brown 358480cf27aSJed Brown PetscFunctionBegin; 35912801b39SBarry Smith if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Unexpected keyval"); 360ec4fadc2SJed Brown icomm.ptr = attr_val; 36176bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 36233779a13SJunchao Zhang /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */ 36333779a13SJunchao Zhang PetscMPIInt flg; 3649371c9d4SSatish Balay union 3659371c9d4SSatish Balay { 3669371c9d4SSatish Balay MPI_Comm comm; 3679371c9d4SSatish Balay void *ptr; 3689371c9d4SSatish Balay } ocomm; 3699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm.comm, Petsc_OuterComm_keyval, &ocomm, &flg)); 37033779a13SJunchao Zhang if (!flg) SETERRMPI(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner comm does not have OuterComm attribute"); 37133779a13SJunchao 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"); 37233779a13SJunchao Zhang } 3739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm.comm, Petsc_OuterComm_keyval)); 3749566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n", (long)comm, (long)icomm.comm)); 375da3039f7SJed Brown PetscFunctionReturn(MPI_SUCCESS); 376b89831e5SBarry Smith } 377da3039f7SJed Brown 378da3039f7SJed Brown /* 37933779a13SJunchao 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. 380da3039f7SJed Brown */ 3819371c9d4SSatish Balay PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm, PetscMPIInt keyval, void *attr_val, void *extra_state) { 382da3039f7SJed Brown PetscFunctionBegin; 3839566063dSJacob Faibussowitsch PetscCallMPI(PetscInfo(NULL, "Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n", (long)comm)); 384480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 385480cf27aSJed Brown } 386480cf27aSJed Brown 38733779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm, PetscMPIInt, void *, void *); 38842218b76SBarry Smith 389951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 3908cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype, MPI_Aint *, void *); 3918cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 3928cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void *, MPI_Datatype, PetscMPIInt, void *, MPI_Offset, void *); 393e39fd77fSBarry Smith #endif 394e39fd77fSBarry Smith 3950f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS = MPI_ERR_LASTCODE, PETSC_MPI_ERROR_CODE; 39612801b39SBarry Smith 397eb27c7c8SSatish Balay PETSC_INTERN int PetscGlobalArgc; 398eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs; 3996ae9a8a6SBarry Smith int PetscGlobalArgc = 0; 40002c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL; 401dff31646SBarry Smith PetscSegBuffer PetscCitationsList; 402e5c89e4eSSatish Balay 4039371c9d4SSatish Balay PetscErrorCode PetscCitationsInitialize(void) { 404051e4cf2SJed Brown PetscFunctionBegin; 4059566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(1, 10000, &PetscCitationsList)); 40630c35bf2SSatish Balay 40730c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@TechReport{petsc-user-ref,\n\ 40830c35bf2SSatish Balay Author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Steven Benson and Jed Brown\n\ 40930c35bf2SSatish Balay and Peter Brune and Kris Buschelman and Emil Constantinescu and Lisandro Dalcin and Alp Dener\n\ 41030c35bf2SSatish Balay and Victor Eijkhout and William~D. Gropp and V\'{a}clav Hapla and Tobin Isaac and Pierre Jolivet\n\ 41130c35bf2SSatish Balay and Dmitry Karpeev and Dinesh Kaushik and Matthew~G. Knepley and Fande Kong and Scott Kruger\n\ 41230c35bf2SSatish Balay and Dave~A. May and Lois Curfman McInnes and Richard Tran Mills and Lawrence Mitchell and Todd Munson\n\ 41330c35bf2SSatish Balay and Jose~E. Roman and Karl Rupp and Patrick Sanan and Jason Sarich and Barry~F. Smith\n\ 41430c35bf2SSatish Balay and Stefano Zampini and Hong Zhang and Hong Zhang and Junchao Zhang},\n\ 41530c35bf2SSatish Balay Title = {{PETSc/TAO} Users Manual},\n\ 41630c35bf2SSatish Balay Number = {ANL-21/39 - Revision 3.17},\n\ 41730c35bf2SSatish Balay Institution = {Argonne National Laboratory},\n\ 4189371c9d4SSatish Balay Year = {2022}\n}\n", 4199371c9d4SSatish Balay NULL)); 42030c35bf2SSatish Balay 42130c35bf2SSatish Balay PetscCall(PetscCitationsRegister("@InProceedings{petsc-efficient,\n\ 42230c35bf2SSatish Balay Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n\ 42330c35bf2SSatish Balay Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n\ 42430c35bf2SSatish Balay Booktitle = {Modern Software Tools in Scientific Computing},\n\ 42530c35bf2SSatish Balay Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n\ 42630c35bf2SSatish Balay Pages = {163--202},\n\ 42730c35bf2SSatish Balay Publisher = {Birkh{\\\"{a}}user Press},\n\ 4289371c9d4SSatish Balay Year = {1997}\n}\n", 4299371c9d4SSatish Balay NULL)); 43030c35bf2SSatish Balay 431051e4cf2SJed Brown PetscFunctionReturn(0); 432051e4cf2SJed Brown } 433e5c89e4eSSatish Balay 4342d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */ 4352d747510SLisandro Dalcin 4369371c9d4SSatish Balay PetscErrorCode PetscSetProgramName(const char name[]) { 4372d747510SLisandro Dalcin PetscFunctionBegin; 4389566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(programname, name, sizeof(programname))); 4392d747510SLisandro Dalcin PetscFunctionReturn(0); 4402d747510SLisandro Dalcin } 4412d747510SLisandro Dalcin 4422d747510SLisandro Dalcin /*@C 4432d747510SLisandro Dalcin PetscGetProgramName - Gets the name of the running program. 4442d747510SLisandro Dalcin 4452d747510SLisandro Dalcin Not Collective 4462d747510SLisandro Dalcin 4472d747510SLisandro Dalcin Input Parameter: 4482d747510SLisandro Dalcin . len - length of the string name 4492d747510SLisandro Dalcin 4502d747510SLisandro Dalcin Output Parameter: 4512d747510SLisandro Dalcin . name - the name of the running program 4522d747510SLisandro Dalcin 4532d747510SLisandro Dalcin Level: advanced 4542d747510SLisandro Dalcin 4552d747510SLisandro Dalcin Notes: 4562d747510SLisandro Dalcin The name of the program is copied into the user-provided character 4572d747510SLisandro Dalcin array of length len. On some machines the program name includes 4582d747510SLisandro Dalcin its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN. 4592d747510SLisandro Dalcin @*/ 4609371c9d4SSatish Balay PetscErrorCode PetscGetProgramName(char name[], size_t len) { 4612d747510SLisandro Dalcin PetscFunctionBegin; 4629566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(name, programname, len)); 4632d747510SLisandro Dalcin PetscFunctionReturn(0); 4642d747510SLisandro Dalcin } 4652d747510SLisandro Dalcin 466e5c89e4eSSatish Balay /*@C 467e5c89e4eSSatish Balay PetscGetArgs - Allows you to access the raw command line arguments anywhere 468e5c89e4eSSatish Balay after PetscInitialize() is called but before PetscFinalize(). 469e5c89e4eSSatish Balay 470e5c89e4eSSatish Balay Not Collective 471e5c89e4eSSatish Balay 472e5c89e4eSSatish Balay Output Parameters: 473e5c89e4eSSatish Balay + argc - count of number of command line arguments 474e5c89e4eSSatish Balay - args - the command line arguments 475e5c89e4eSSatish Balay 476e5c89e4eSSatish Balay Level: intermediate 477e5c89e4eSSatish Balay 478e5c89e4eSSatish Balay Notes: 479e5c89e4eSSatish Balay This is usually used to pass the command line arguments into other libraries 480e5c89e4eSSatish Balay that are called internally deep in PETSc or the application. 481e5c89e4eSSatish Balay 482f177e3b1SBarry Smith The first argument contains the program name as is normal for C arguments. 483f177e3b1SBarry Smith 484db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArguments()` 485e5c89e4eSSatish Balay 486e5c89e4eSSatish Balay @*/ 4879371c9d4SSatish Balay PetscErrorCode PetscGetArgs(int *argc, char ***args) { 488e5c89e4eSSatish Balay PetscFunctionBegin; 48939a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 490e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 491e5c89e4eSSatish Balay *args = PetscGlobalArgs; 492e5c89e4eSSatish Balay PetscFunctionReturn(0); 493e5c89e4eSSatish Balay } 494e5c89e4eSSatish Balay 495793721a6SBarry Smith /*@C 496793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 497793721a6SBarry Smith after PetscInitialize() is called but before PetscFinalize(). 498793721a6SBarry Smith 499793721a6SBarry Smith Not Collective 500793721a6SBarry Smith 501793721a6SBarry Smith Output Parameters: 502793721a6SBarry Smith . args - the command line arguments 503793721a6SBarry Smith 504793721a6SBarry Smith Level: intermediate 505793721a6SBarry Smith 506793721a6SBarry Smith Notes: 507793721a6SBarry Smith This does NOT start with the program name and IS null terminated (final arg is void) 508793721a6SBarry Smith 509db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscFreeArguments()` 510793721a6SBarry Smith 511793721a6SBarry Smith @*/ 5129371c9d4SSatish Balay PetscErrorCode PetscGetArguments(char ***args) { 513793721a6SBarry Smith PetscInt i, argc = PetscGlobalArgc; 514793721a6SBarry Smith 515793721a6SBarry Smith PetscFunctionBegin; 51639a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled || !PetscFinalizeCalled, PETSC_COMM_SELF, PETSC_ERR_ORDER, "You must call after PetscInitialize() but before PetscFinalize()"); 5179371c9d4SSatish Balay if (!argc) { 5189371c9d4SSatish Balay *args = NULL; 5199371c9d4SSatish Balay PetscFunctionReturn(0); 5209371c9d4SSatish Balay } 5219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(argc, args)); 5229566063dSJacob Faibussowitsch for (i = 0; i < argc - 1; i++) PetscCall(PetscStrallocpy(PetscGlobalArgs[i + 1], &(*args)[i])); 5232d747510SLisandro Dalcin (*args)[argc - 1] = NULL; 524793721a6SBarry Smith PetscFunctionReturn(0); 525793721a6SBarry Smith } 526793721a6SBarry Smith 527793721a6SBarry Smith /*@C 528793721a6SBarry Smith PetscFreeArguments - Frees the memory obtained with PetscGetArguments() 529793721a6SBarry Smith 530793721a6SBarry Smith Not Collective 531793721a6SBarry Smith 532793721a6SBarry Smith Output Parameters: 533793721a6SBarry Smith . args - the command line arguments 534793721a6SBarry Smith 535793721a6SBarry Smith Level: intermediate 536793721a6SBarry Smith 537db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscGetArguments()` 538793721a6SBarry Smith 539793721a6SBarry Smith @*/ 5409371c9d4SSatish Balay PetscErrorCode PetscFreeArguments(char **args) { 54139a651e2SJacob Faibussowitsch PetscFunctionBegin; 54239a651e2SJacob Faibussowitsch if (args) { 543793721a6SBarry Smith PetscInt i = 0; 544793721a6SBarry Smith 5459566063dSJacob Faibussowitsch while (args[i]) PetscCall(PetscFree(args[i++])); 5469566063dSJacob Faibussowitsch PetscCall(PetscFree(args)); 54739a651e2SJacob Faibussowitsch } 548793721a6SBarry Smith PetscFunctionReturn(0); 549793721a6SBarry Smith } 550793721a6SBarry Smith 55127104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_SAWS) 55230befbd2SBarry Smith #include <petscconfiginfo.h> 55330befbd2SBarry Smith 5549371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) { 55527104ee2SJacob Faibussowitsch PetscFunctionBegin; 55611525c0dSBarry Smith if (!PetscGlobalRank) { 55730befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN], root[PETSC_MAX_PATH_LEN], *intro, programname[64], *appline, *options, version[64]; 55811525c0dSBarry Smith int port; 559ffbd1cfbSBarry Smith PetscBool flg, rootlocal = PETSC_FALSE, flg2, selectport = PETSC_FALSE; 56011525c0dSBarry Smith size_t applinelen, introlen; 561ffbd1cfbSBarry Smith char sawsurl[256]; 56211525c0dSBarry Smith 5639566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_log", &flg)); 56411525c0dSBarry Smith if (flg) { 56511525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 56611525c0dSBarry Smith 5679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_log", sawslog, sizeof(sawslog), NULL)); 56811525c0dSBarry Smith if (sawslog[0]) { 569792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (sawslog)); 57011525c0dSBarry Smith } else { 571792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Use_Logfile, (NULL)); 57211525c0dSBarry Smith } 57311525c0dSBarry Smith } 5749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_https", cert, sizeof(cert), &flg)); 5759371c9d4SSatish Balay if (flg) { PetscCallSAWs(SAWs_Set_Use_HTTPS, (cert)); } 5769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select", &selectport, NULL)); 577ffbd1cfbSBarry Smith if (selectport) { 578792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 579792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 580ffbd1cfbSBarry Smith } else { 5819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-saws_port", &port, &flg)); 5829371c9d4SSatish Balay if (flg) { PetscCallSAWs(SAWs_Set_Port, (port)); } 583ffbd1cfbSBarry Smith } 5849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-saws_root", root, sizeof(root), &flg)); 58511525c0dSBarry Smith if (flg) { 586792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 5879566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(root, ".", &rootlocal)); 5889c1e0ce8SBarry Smith } else { 5899566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_options", &flg)); 5909c1e0ce8SBarry Smith if (flg) { 5919566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/saws", root, sizeof(root))); 592792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Document_Root, (root)); 5939c1e0ce8SBarry Smith } 59411525c0dSBarry Smith } 5959566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-saws_local", &flg2)); 59611525c0dSBarry Smith if (flg2) { 59711525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 59828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "-saws_local option requires -saws_root option"); 5999566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(jsdir, sizeof(jsdir), "%s/js", root)); 6009566063dSJacob Faibussowitsch PetscCall(PetscTestDirectory(jsdir, 'r', &flg)); 60128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "-saws_local option requires js directory in root directory"); 602792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Local_Header, ()); 60311525c0dSBarry Smith } 6049566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(programname, sizeof(programname))); 6059566063dSJacob Faibussowitsch PetscCall(PetscStrlen(help, &applinelen)); 60611525c0dSBarry Smith introlen = 4096 + applinelen; 60730a8c9c0SSurtai Han applinelen += 1024; 6089566063dSJacob Faibussowitsch PetscCall(PetscMalloc(applinelen, &appline)); 6099566063dSJacob Faibussowitsch PetscCall(PetscMalloc(introlen, &intro)); 61011525c0dSBarry Smith 61111525c0dSBarry Smith if (rootlocal) { 6129566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "%s.c.html", programname)); 6139566063dSJacob Faibussowitsch PetscCall(PetscTestFile(appline, 'r', &rootlocal)); 61411525c0dSBarry Smith } 6159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetAll(NULL, &options)); 61611525c0dSBarry Smith if (rootlocal && help) { 6179566063dSJacob 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)); 61811525c0dSBarry Smith } else if (help) { 6199566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center>Running %s %s</center><br><center><pre>%s</pre></center><br>", programname, options, help)); 62011525c0dSBarry Smith } else { 6219566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(appline, applinelen, "<center> Running %s %s</center><br>\n", programname, options)); 62211525c0dSBarry Smith } 6239566063dSJacob Faibussowitsch PetscCall(PetscFree(options)); 6249566063dSJacob Faibussowitsch PetscCall(PetscGetVersion(version, sizeof(version))); 6259371c9d4SSatish Balay PetscCall(PetscSNPrintf(intro, introlen, 6269371c9d4SSatish Balay "<body>\n" 627a17b96a8SKyle 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" 628df62c222SBarry 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" 6299371c9d4SSatish Balay "%s", 6309371c9d4SSatish Balay version, petscconfigureoptions, appline)); 631792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Body, ("index.html", 0, intro)); 6329566063dSJacob Faibussowitsch PetscCall(PetscFree(intro)); 6339566063dSJacob Faibussowitsch PetscCall(PetscFree(appline)); 634ffbd1cfbSBarry Smith if (selectport) { 635aa573868SBarry Smith PetscBool silent; 6367d812c46SBarry Smith 6377d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 63839a651e2SJacob Faibussowitsch while (SAWs_Initialize()) { 639792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_Available_Port, (&port)); 640792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Port, (port)); 6417d812c46SBarry Smith } 6427d812c46SBarry Smith 6439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saws_port_auto_select_silent", &silent, NULL)); 644aa573868SBarry Smith if (!silent) { 645792fecdfSBarry Smith PetscCallSAWs(SAWs_Get_FullURL, (sizeof(sawsurl), sawsurl)); 6469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Point your browser to %s for SAWs\n", sawsurl)); 647ffbd1cfbSBarry Smith } 6487d812c46SBarry Smith } else { 649792fecdfSBarry Smith PetscCallSAWs(SAWs_Initialize, ()); 650aa573868SBarry Smith } 6519566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister("@TechReport{ saws,\n" 6520af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 6530af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 6540af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 6559371c9d4SSatish Balay " Year = 2013\n}\n", 6569371c9d4SSatish Balay NULL)); 65711525c0dSBarry Smith } 65811525c0dSBarry Smith PetscFunctionReturn(0); 65911525c0dSBarry Smith } 66011525c0dSBarry Smith #endif 66111525c0dSBarry Smith 6624dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 6639371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) { 6644dfee713SSatish Balay PetscFunctionBegin; 6654dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 6664dfee713SSatish Balay /* see MPI.py for details on this bug */ 6674dfee713SSatish Balay (void)setenv("HWLOC_COMPONENTS", "-x86", 1); 6684dfee713SSatish Balay #endif 6694dfee713SSatish Balay PetscFunctionReturn(0); 6704dfee713SSatish Balay } 6714dfee713SSatish Balay 672a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_ADIOS) 673a56f64adSBarry Smith #include <adios.h> 67422580e64SBarry Smith #include <adios_read.h> 6757b56e58cSSatish Balay int64_t Petsc_adios_group; 676a56f64adSBarry Smith #endif 677a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_OPENMP) 678cd1b99a6SBarry Smith #include <omp.h> 679f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 680cd1b99a6SBarry Smith #endif 681a56f64adSBarry Smith 682a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_DEVICE) 683a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h> 684a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA) 685a4af0ceeSJacob Faibussowitsch // REMOVE ME 686a4af0ceeSJacob Faibussowitsch cudaStream_t PetscDefaultCudaStream = NULL; 687a4af0ceeSJacob Faibussowitsch #endif 688a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_HIP) 689a4af0ceeSJacob Faibussowitsch // REMOVE ME 690a4af0ceeSJacob Faibussowitsch hipStream_t PetscDefaultHipStream = NULL; 691a4af0ceeSJacob Faibussowitsch #endif 692a4af0ceeSJacob Faibussowitsch #endif 693a4af0ceeSJacob Faibussowitsch 69427104ee2SJacob Faibussowitsch #if PetscDefined(HAVE_DLFCN_H) 695bc8a24ecSBarry Smith #include <dlfcn.h> 696bc8a24ecSBarry Smith #endif 697a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 698a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogInitialize(void); 699a4af0ceeSJacob Faibussowitsch #endif 700a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 701a4af0ceeSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViennaCLInit(); 702a4af0ceeSJacob Faibussowitsch PetscBool PetscViennaCLSynchronize = PETSC_FALSE; 703a4af0ceeSJacob Faibussowitsch #endif 704bc8a24ecSBarry Smith 705660278c0SBarry Smith PetscBool PetscCIEnabled = PETSC_FALSE, PetscCIEnabledPortableErrorOutput = PETSC_FALSE; 706660278c0SBarry Smith 70785649d77SJunchao Zhang /* 70885649d77SJunchao Zhang PetscInitialize_Common - shared code between C and Fortran initialization 70985649d77SJunchao Zhang 71085649d77SJunchao Zhang prog: program name 71102101c96SSatish Balay file: optional PETSc database file name. Might be in Fortran string format when 'ftn' is true 71285649d77SJunchao Zhang help: program help message 71302101c96SSatish Balay ftn: is it called from Fortran initilization (petscinitializef_)? 71485649d77SJunchao Zhang readarguments,len: used when fortran is true 71585649d77SJunchao Zhang */ 7169371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char *prog, const char *file, const char *help, PetscBool ftn, PetscBool readarguments, PetscInt len) { 71785649d77SJunchao Zhang PetscMPIInt size; 71885649d77SJunchao Zhang PetscBool flg = PETSC_TRUE; 71985649d77SJunchao Zhang char hostname[256]; 72085649d77SJunchao Zhang 72127104ee2SJacob Faibussowitsch PetscFunctionBegin; 72227104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(0); 72339a651e2SJacob Faibussowitsch /* these must be initialized in a routine, not as a constant declaration */ 72439a651e2SJacob Faibussowitsch PETSC_STDOUT = stdout; 72539a651e2SJacob Faibussowitsch PETSC_STDERR = stderr; 72639a651e2SJacob Faibussowitsch 7279566063dSJacob Faibussowitsch /* PetscCall can be used from now */ 72839a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_TRUE; 72939a651e2SJacob Faibussowitsch 73085649d77SJunchao Zhang /* 73185649d77SJunchao Zhang The checking over compatible runtime libraries is complicated by the MPI ABI initiative 73285649d77SJunchao Zhang https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 73385649d77SJunchao Zhang MPICH v3.1 (Released February 2014) 73485649d77SJunchao Zhang IBM MPI v2.1 (December 2014) 73585649d77SJunchao Zhang Intel MPI Library v5.0 (2014) 73685649d77SJunchao Zhang Cray MPT v7.0.0 (June 2014) 73785649d77SJunchao Zhang As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 73885649d77SJunchao Zhang listed above and since that time are compatible. 73985649d77SJunchao Zhang 74085649d77SJunchao Zhang Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 74185649d77SJunchao Zhang at compile time or runtime. Thus we will need to systematically track the allowed versions 74285649d77SJunchao Zhang and how they are represented in the mpi.h and MPI_Get_library_version() output in order 74385649d77SJunchao Zhang to perform the checking. 74485649d77SJunchao Zhang 74585649d77SJunchao Zhang Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 74685649d77SJunchao Zhang 74785649d77SJunchao Zhang Questions: 74885649d77SJunchao Zhang 74985649d77SJunchao Zhang Should the checks for ABI incompatibility be only on the major version number below? 75085649d77SJunchao Zhang Presumably the output to stderr will be removed before a release. 75185649d77SJunchao Zhang */ 75285649d77SJunchao Zhang 75385649d77SJunchao Zhang #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 75485649d77SJunchao Zhang { 75585649d77SJunchao Zhang char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 75685649d77SJunchao Zhang PetscMPIInt mpilibraryversionlength; 75739a651e2SJacob Faibussowitsch 7589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_library_version(mpilibraryversion, &mpilibraryversionlength)); 75985649d77SJunchao Zhang /* check for MPICH versions before MPI ABI initiative */ 76085649d77SJunchao Zhang #if defined(MPICH_VERSION) 76185649d77SJunchao Zhang #if MPICH_NUMVERSION < 30100000 76285649d77SJunchao Zhang { 76385649d77SJunchao Zhang char *ver, *lf; 76485649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 76539a651e2SJacob Faibussowitsch 7669566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, "MPICH Version:", &ver)); 76739a651e2SJacob Faibussowitsch if (ver) { 7689566063dSJacob Faibussowitsch PetscCall(PetscStrchr(ver, '\n', &lf)); 76939a651e2SJacob Faibussowitsch if (lf) { 77085649d77SJunchao Zhang *lf = 0; 7719566063dSJacob Faibussowitsch PetscCall(PetscStrendswith(ver, MPICH_VERSION, &flg)); 77285649d77SJunchao Zhang } 77385649d77SJunchao Zhang } 77485649d77SJunchao Zhang if (!flg) { 775d5b396e8SSatish Balay PetscCall(PetscInfo(NULL, "PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n", mpilibraryversion, MPICH_VERSION)); 77685649d77SJunchao Zhang flg = PETSC_TRUE; 77785649d77SJunchao Zhang } 77885649d77SJunchao Zhang } 77985649d77SJunchao Zhang #endif 78085649d77SJunchao 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?) */ 78185649d77SJunchao Zhang #elif defined(OMPI_MAJOR_VERSION) 78285649d77SJunchao Zhang { 78385649d77SJunchao Zhang char *ver, bs[MPI_MAX_LIBRARY_VERSION_STRING], *bsf; 78485649d77SJunchao Zhang PetscBool flg = PETSC_FALSE; 78585649d77SJunchao Zhang #define PSTRSZ 2 78685649d77SJunchao Zhang char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI", "FUJITSU MPI"}; 78785649d77SJunchao Zhang char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v", "Library "}; 78885649d77SJunchao Zhang int i; 78985649d77SJunchao Zhang for (i = 0; i < PSTRSZ; i++) { 7909566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mpilibraryversion, ompistr1[i], &ver)); 79139a651e2SJacob Faibussowitsch if (ver) { 7929566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(bs, MPI_MAX_LIBRARY_VERSION_STRING, "%s%d.%d", ompistr2[i], OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION)); 7939566063dSJacob Faibussowitsch PetscCall(PetscStrstr(ver, bs, &bsf)); 79439a651e2SJacob Faibussowitsch if (bsf) flg = PETSC_TRUE; 79585649d77SJunchao Zhang break; 79685649d77SJunchao Zhang } 79785649d77SJunchao Zhang } 79885649d77SJunchao Zhang if (!flg) { 7997d3de750SJacob 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); 80085649d77SJunchao Zhang flg = PETSC_TRUE; 80185649d77SJunchao Zhang } 80285649d77SJunchao Zhang } 80385649d77SJunchao Zhang #endif 80485649d77SJunchao Zhang } 80585649d77SJunchao Zhang #endif 80685649d77SJunchao Zhang 8076f5d4113SSatish Balay #if PetscDefined(HAVE_DLSYM) && defined(__USE_GNU) 80885649d77SJunchao 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 */ 80939a651e2SJacob 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"); 81085649d77SJunchao Zhang #endif 81185649d77SJunchao Zhang 81285649d77SJunchao Zhang /* on Windows - set printf to default to printing 2 digit exponents */ 81385649d77SJunchao Zhang #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT) 81485649d77SJunchao Zhang _set_output_format(_TWO_DIGIT_EXPONENT); 81585649d77SJunchao Zhang #endif 81685649d77SJunchao Zhang 8179566063dSJacob Faibussowitsch PetscCall(PetscOptionsCreateDefault()); 81885649d77SJunchao Zhang 81985649d77SJunchao Zhang PetscFinalizeCalled = PETSC_FALSE; 82085649d77SJunchao Zhang 8219566063dSJacob Faibussowitsch PetscCall(PetscSetProgramName(prog)); 8229566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen)); 8239566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout)); 8249566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr)); 8259566063dSJacob Faibussowitsch PetscCall(PetscSpinlockCreate(&PetscCommSpinLock)); 82685649d77SJunchao Zhang 82785649d77SJunchao Zhang if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 8289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_set_errhandler(PETSC_COMM_WORLD, MPI_ERRORS_RETURN)); 82985649d77SJunchao Zhang 83085649d77SJunchao Zhang if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) { 8319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS)); 8329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Add_error_code(PETSC_MPI_ERROR_CLASS, &PETSC_MPI_ERROR_CODE)); 83385649d77SJunchao Zhang } 83485649d77SJunchao Zhang 83585649d77SJunchao Zhang /* Done after init due to a bug in MPICH-GM? */ 8369566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 83785649d77SJunchao Zhang 8389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PetscGlobalRank)); 8399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PetscGlobalSize)); 84085649d77SJunchao Zhang 84185649d77SJunchao Zhang MPIU_BOOL = MPI_INT; 84285649d77SJunchao Zhang MPIU_ENUM = MPI_INT; 84385649d77SJunchao Zhang MPIU_FORTRANADDR = (sizeof(void *) == sizeof(int)) ? MPI_INT : MPIU_INT64; 84485649d77SJunchao Zhang if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 84585649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 84685649d77SJunchao Zhang #if defined(PETSC_SIZEOF_LONG_LONG) 84785649d77SJunchao Zhang else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 84885649d77SJunchao Zhang #endif 84939a651e2SJacob Faibussowitsch else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "Could not find MPI type for size_t"); 85085649d77SJunchao Zhang 85185649d77SJunchao Zhang /* 85285649d77SJunchao Zhang Initialized the global complex variable; this is because with 85385649d77SJunchao Zhang shared libraries the constructors for global variables 85485649d77SJunchao Zhang are not called; at least on IRIX. 85585649d77SJunchao Zhang */ 85685649d77SJunchao Zhang #if defined(PETSC_HAVE_COMPLEX) 85785649d77SJunchao Zhang { 85885649d77SJunchao Zhang #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 85985649d77SJunchao Zhang PetscComplex ic(0.0, 1.0); 86085649d77SJunchao Zhang PETSC_i = ic; 86185649d77SJunchao Zhang #else 86285649d77SJunchao Zhang PETSC_i = _Complex_I; 86385649d77SJunchao Zhang #endif 86485649d77SJunchao Zhang } 86585649d77SJunchao Zhang #endif /* PETSC_HAVE_COMPLEX */ 86685649d77SJunchao Zhang 86785649d77SJunchao Zhang /* 86885649d77SJunchao Zhang Create the PETSc MPI reduction operator that sums of the first 86985649d77SJunchao Zhang half of the entries and maxes the second half. 87085649d77SJunchao Zhang */ 8719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_create(MPIU_MaxSum_Local, 1, &MPIU_MAXSUM_OP)); 87285649d77SJunchao Zhang 873613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 8749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPI_DOUBLE, &MPIU___FLOAT128)); 8759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___FLOAT128)); 8769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_DOUBLE, &MPIU___COMPLEX128)); 8779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU___COMPLEX128)); 878613bf2b2SPierre Jolivet #if !defined(PETSC_USE_REAL___FLOAT128) 879613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_create(PetscSum_Local, 1, &MPIU_SUM___FLOAT128)); 88085649d77SJunchao Zhang #endif 881613bf2b2SPierre Jolivet #endif 882613bf2b2SPierre Jolivet #if defined(PETSC_USE_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)); 89185649d77SJunchao Zhang #endif 89285649d77SJunchao Zhang 8939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_SCALAR, &MPIU_2SCALAR)); 8949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2SCALAR)); 89585649d77SJunchao Zhang 89685649d77SJunchao Zhang /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */ 89785649d77SJunchao Zhang #if !defined(PETSC_HAVE_MPIUNI) 89885649d77SJunchao Zhang { 8999371c9d4SSatish Balay struct PetscRealInt { 9009371c9d4SSatish Balay PetscReal v; 9019371c9d4SSatish Balay PetscInt i; 9029371c9d4SSatish Balay }; 90385649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 90485649d77SJunchao Zhang MPI_Aint blockOffsets[2] = {offsetof(struct PetscRealInt, v), offsetof(struct PetscRealInt, i)}; 90585649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_REAL, MPIU_INT}, tmpStruct; 90685649d77SJunchao Zhang 9079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 9089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct PetscRealInt), &MPIU_REAL_INT)); 9099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_REAL_INT)); 91185649d77SJunchao Zhang } 91285649d77SJunchao Zhang { 9139371c9d4SSatish Balay struct PetscScalarInt { 9149371c9d4SSatish Balay PetscScalar v; 9159371c9d4SSatish Balay PetscInt i; 9169371c9d4SSatish Balay }; 91785649d77SJunchao Zhang PetscMPIInt blockSizes[2] = {1, 1}; 91885649d77SJunchao Zhang MPI_Aint blockOffsets[2] = {offsetof(struct PetscScalarInt, v), offsetof(struct PetscScalarInt, i)}; 91985649d77SJunchao Zhang MPI_Datatype blockTypes[2] = {MPIU_SCALAR, MPIU_INT}, tmpStruct; 92085649d77SJunchao Zhang 9219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(2, blockSizes, blockOffsets, blockTypes, &tmpStruct)); 9229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmpStruct, 0, sizeof(struct PetscScalarInt), &MPIU_SCALAR_INT)); 9239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmpStruct)); 9249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_SCALAR_INT)); 92585649d77SJunchao Zhang } 92685649d77SJunchao Zhang #endif 92785649d77SJunchao Zhang 92885649d77SJunchao Zhang #if defined(PETSC_USE_64BIT_INDICES) 9299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(2, MPIU_INT, &MPIU_2INT)); 9309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_2INT)); 93185649d77SJunchao Zhang #endif 9329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPI_INT, &MPI_4INT)); 9339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPI_4INT)); 9349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_contiguous(4, MPIU_INT, &MPIU_4INT)); 9359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&MPIU_4INT)); 93685649d77SJunchao Zhang 93785649d77SJunchao Zhang /* 93885649d77SJunchao Zhang Attributes to be set on PETSc communicators 93985649d77SJunchao Zhang */ 9409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_Counter_Attr_Delete_Fn, &Petsc_Counter_keyval, (void *)0)); 9419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_InnerComm_Attr_Delete_Fn, &Petsc_InnerComm_keyval, (void *)0)); 9429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_OuterComm_Attr_Delete_Fn, &Petsc_OuterComm_keyval, (void *)0)); 9439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, Petsc_ShmComm_Attr_Delete_Fn, &Petsc_ShmComm_keyval, (void *)0)); 94485649d77SJunchao Zhang 94585649d77SJunchao Zhang #if defined(PETSC_HAVE_FORTRAN) 9469566063dSJacob Faibussowitsch if (ftn) PetscCall(PetscInitFortran_Private(readarguments, file, len)); 94785649d77SJunchao Zhang else 94885649d77SJunchao Zhang #endif 9499566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL, &PetscGlobalArgc, &PetscGlobalArgs, file)); 95085649d77SJunchao Zhang 95185649d77SJunchao Zhang /* call a second time so it can look in the options database */ 9529566063dSJacob Faibussowitsch PetscCall(PetscErrorPrintfInitialize()); 95385649d77SJunchao Zhang 95485649d77SJunchao Zhang /* 95585649d77SJunchao Zhang Check system options and print help 95685649d77SJunchao Zhang */ 9579566063dSJacob Faibussowitsch PetscCall(PetscOptionsCheckInitial_Private(help)); 95885649d77SJunchao Zhang 959a4af0ceeSJacob Faibussowitsch /* 960a4af0ceeSJacob Faibussowitsch Initialize PetscDevice and PetscDeviceContext 961a4af0ceeSJacob Faibussowitsch 962a4af0ceeSJacob Faibussowitsch Note to any future devs thinking of moving this, proper initialization requires: 963a4af0ceeSJacob Faibussowitsch 1. MPI initialized 964a4af0ceeSJacob Faibussowitsch 2. Options DB initialized 965a4af0ceeSJacob Faibussowitsch 3. Petsc error handling initialized, specifically signal handlers. This expects to set up its own SIGSEV handler via 966a4af0ceeSJacob Faibussowitsch the push/pop interface. 967a4af0ceeSJacob Faibussowitsch */ 968a2158755SJunchao Zhang #if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP) || PetscDefined(HAVE_SYCL)) 9699566063dSJacob Faibussowitsch PetscCall(PetscDeviceInitializeFromOptions_Internal(PETSC_COMM_WORLD)); 970a4af0ceeSJacob Faibussowitsch #endif 971a4af0ceeSJacob Faibussowitsch 972a4af0ceeSJacob Faibussowitsch #if PetscDefined(HAVE_VIENNACL) 973a4af0ceeSJacob Faibussowitsch flg = PETSC_FALSE; 9749566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-log_summary", &flg)); 9759566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsHasName(NULL, NULL, "-log_view", &flg)); 9769566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsGetBool(NULL, NULL, "-viennacl_synchronize", &flg, NULL)); 977a4af0ceeSJacob Faibussowitsch PetscViennaCLSynchronize = flg; 9789566063dSJacob Faibussowitsch PetscCall(PetscViennaCLInit()); 979a4af0ceeSJacob Faibussowitsch #endif 980a4af0ceeSJacob Faibussowitsch 981a4af0ceeSJacob Faibussowitsch /* 982a4af0ceeSJacob Faibussowitsch Creates the logging data structures; this is enabled even if logging is not turned on 983a4af0ceeSJacob Faibussowitsch This is the last thing we do before returning to the user code to prevent having the 984a4af0ceeSJacob Faibussowitsch logging numbers contaminated by any startup time associated with MPI 985a4af0ceeSJacob Faibussowitsch */ 986a4af0ceeSJacob Faibussowitsch #if defined(PETSC_USE_LOG) 9879566063dSJacob Faibussowitsch PetscCall(PetscLogInitialize()); 988a4af0ceeSJacob Faibussowitsch #endif 989a4af0ceeSJacob Faibussowitsch 9909566063dSJacob Faibussowitsch PetscCall(PetscCitationsInitialize()); 99185649d77SJunchao Zhang 99285649d77SJunchao Zhang #if defined(PETSC_HAVE_SAWS) 9939566063dSJacob Faibussowitsch PetscCall(PetscInitializeSAWs(ftn ? NULL : help)); 99427104ee2SJacob Faibussowitsch flg = PETSC_FALSE; 9959566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-stack_view", &flg)); 9969566063dSJacob Faibussowitsch if (flg) PetscCall(PetscStackViewSAWs()); 99785649d77SJunchao Zhang #endif 99885649d77SJunchao Zhang 99985649d77SJunchao Zhang /* 100085649d77SJunchao Zhang Load the dynamic libraries (on machines that support them), this registers all 100185649d77SJunchao Zhang the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 100285649d77SJunchao Zhang */ 10039566063dSJacob Faibussowitsch PetscCall(PetscInitialize_DynamicLibraries()); 100485649d77SJunchao Zhang 10059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 10069566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PETSc successfully started: number of processors = %d\n", size)); 10079566063dSJacob Faibussowitsch PetscCall(PetscGetHostName(hostname, 256)); 10089566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Running on machine: %s\n", hostname)); 100985649d77SJunchao Zhang #if defined(PETSC_HAVE_OPENMP) 101085649d77SJunchao Zhang { 101185649d77SJunchao Zhang PetscBool omp_view_flag; 101285649d77SJunchao Zhang char *threads = getenv("OMP_NUM_THREADS"); 101385649d77SJunchao Zhang 101485649d77SJunchao Zhang if (threads) { 10159566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n", threads)); 101685649d77SJunchao Zhang (void)sscanf(threads, "%" PetscInt_FMT, &PetscNumOMPThreads); 101785649d77SJunchao Zhang } else { 10182f840973SStefano Zampini PetscNumOMPThreads = (PetscInt)omp_get_max_threads(); 10199566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP threads %" PetscInt_FMT " (as given by omp_get_max_threads())\n", PetscNumOMPThreads)); 102085649d77SJunchao Zhang } 1021d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "OpenMP options", "Sys"); 10229566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-omp_num_threads", "Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS", "None", PetscNumOMPThreads, &PetscNumOMPThreads, &flg)); 10239566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-omp_view", "Display OpenMP number of threads", NULL, &omp_view_flag)); 1024d0609cedSBarry Smith PetscOptionsEnd(); 102585649d77SJunchao Zhang if (flg) { 10269566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Number of OpenMP theads %" PetscInt_FMT " (given by -omp_num_threads)\n", PetscNumOMPThreads)); 102785649d77SJunchao Zhang omp_set_num_threads((int)PetscNumOMPThreads); 102885649d77SJunchao Zhang } 10299371c9d4SSatish Balay if (omp_view_flag) { PetscCall(PetscPrintf(PETSC_COMM_WORLD, "OpenMP: number of threads %" PetscInt_FMT "\n", PetscNumOMPThreads)); } 103085649d77SJunchao Zhang } 103185649d77SJunchao Zhang #endif 103285649d77SJunchao Zhang 103385649d77SJunchao Zhang #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 103485649d77SJunchao Zhang /* 103585649d77SJunchao Zhang Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 103685649d77SJunchao Zhang 103785649d77SJunchao Zhang Currently not used because it is not supported by MPICH. 103885649d77SJunchao Zhang */ 10399566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCallMPI(MPI_Register_datarep((char *)"petsc", PetscDataRep_read_conv_fn, PetscDataRep_write_conv_fn, PetscDataRep_extent_fn, NULL)); 104085649d77SJunchao Zhang #endif 104185649d77SJunchao Zhang 104285649d77SJunchao Zhang #if defined(PETSC_SERIALIZE_FUNCTIONS) 10439566063dSJacob Faibussowitsch PetscCall(PetscFPTCreate(10000)); 104485649d77SJunchao Zhang #endif 104585649d77SJunchao Zhang 104685649d77SJunchao Zhang #if defined(PETSC_HAVE_HWLOC) 104785649d77SJunchao Zhang { 104885649d77SJunchao Zhang PetscViewer viewer; 10499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-process_view", &viewer, NULL, &flg)); 105085649d77SJunchao Zhang if (flg) { 10519566063dSJacob Faibussowitsch PetscCall(PetscProcessPlacementView(viewer)); 10529566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 105385649d77SJunchao Zhang } 105485649d77SJunchao Zhang } 105585649d77SJunchao Zhang #endif 105685649d77SJunchao Zhang 105785649d77SJunchao Zhang flg = PETSC_TRUE; 10589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewfromoptions", &flg, NULL)); 10599566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscOptionsPushGetViewerOff(PETSC_TRUE)); 106085649d77SJunchao Zhang 106185649d77SJunchao Zhang #if defined(PETSC_HAVE_ADIOS) 10629566063dSJacob Faibussowitsch PetscCall(adios_init_noxml(PETSC_COMM_WORLD)); 10639566063dSJacob Faibussowitsch PetscCall(adios_declare_group(&Petsc_adios_group, "PETSc", "", adios_stat_default)); 10649566063dSJacob Faibussowitsch PetscCall(adios_select_method(Petsc_adios_group, "MPI", "", "")); 10659566063dSJacob Faibussowitsch PetscCall(adios_read_init_method(ADIOS_READ_METHOD_BP, PETSC_COMM_WORLD, "")); 106685649d77SJunchao Zhang #endif 106785649d77SJunchao Zhang 106885649d77SJunchao Zhang #if defined(__VALGRIND_H) 106985649d77SJunchao Zhang PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND ? PETSC_TRUE : PETSC_FALSE; 107085649d77SJunchao Zhang #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE) 10719566063dSJacob 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")); 107285649d77SJunchao Zhang #endif 107385649d77SJunchao Zhang #endif 107485649d77SJunchao Zhang /* 107585649d77SJunchao Zhang Set flag that we are completely initialized 107685649d77SJunchao Zhang */ 107785649d77SJunchao Zhang PetscInitializeCalled = PETSC_TRUE; 107885649d77SJunchao Zhang 10799566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-python", &flg)); 10809566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPythonInitialize(NULL, NULL)); 1081f1f2ae84SBarry Smith 1082f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1083f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerBegin()); 1084f1f2ae84SBarry 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"); 108527104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 108685649d77SJunchao Zhang } 108785649d77SJunchao Zhang 1088e5c89e4eSSatish Balay /*@C 1089e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 1090e5c89e4eSSatish Balay PetscInitialize() calls MPI_Init() if that has yet to be called, 1091e5c89e4eSSatish Balay so this routine should always be called near the beginning of 1092e5c89e4eSSatish Balay your program -- usually the very first line! 1093e5c89e4eSSatish Balay 1094e5c89e4eSSatish Balay Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set 1095e5c89e4eSSatish Balay 1096e5c89e4eSSatish Balay Input Parameters: 1097e5c89e4eSSatish Balay + argc - count of number of command line arguments 1098e5c89e4eSSatish Balay . args - the command line arguments 1099be10d61cSLisandro Dalcin . file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format. 1100be10d61cSLisandro Dalcin Use NULL or empty string to not check for code specific file. 1101be10d61cSLisandro Dalcin Also checks ~/.petscrc, .petscrc and petscrc. 1102c5b5d8d5SVaclav Hapla Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files. 11030298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message 1104e5c89e4eSSatish Balay 110505827820SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that 110605827820SBarry Smith communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a 110705827820SBarry Smith four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not, 110805827820SBarry Smith then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even 110905827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 1110e5c89e4eSSatish Balay 1111e5c89e4eSSatish Balay Options Database Keys: 11127ca660e7SBarry Smith + -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message 11137ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 1114e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 11158a690491SBarry Smith . -on_error_emacs <machinename> - causes emacsclient to jump to error file 11168a690491SBarry Smith . -on_error_abort - calls abort() when error detected (no traceback) 11178a690491SBarry Smith . -on_error_mpiabort - calls MPI_abort() when error detected 11181af3601dSBarry Smith . -error_output_stdout - prints PETSc error messages to stdout instead of the default stderr 11198a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 1120bf4d2887SBarry Smith . -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger 1121e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds) - Pauses debugger 1122e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 1123e5c89e4eSSatish Balay process and wait (-debugger_pause) seconds for attachment 112479dccf82SBarry Smith . -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug) 112579dccf82SBarry Smith . -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no) 112679dccf82SBarry Smith . -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug() 1127aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 112892f119d6SBarry 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 112992f119d6SBarry Smith . -malloc_view - show a list of all allocated memory during PetscFinalize() 113092f119d6SBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view 1131608c71bfSMatthew G. Knepley . -malloc_requested_size - malloc logging will record the requested size rather than size after alignment 113292f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 1133e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 1134e5c89e4eSSatish Balay . -shared_tmp - indicates /tmp directory is shared by all processors 1135e5c89e4eSSatish Balay . -not_shared_tmp - each processor has own /tmp 1136e5c89e4eSSatish Balay . -tmp - alternative name of /tmp directory 1137e5c89e4eSSatish Balay . -get_total_flops - returns total flops done by all processors 11380841954dSBarry Smith - -memory_view - Print memory usage at end of run 1139e5c89e4eSSatish Balay 1140c5b5d8d5SVaclav Hapla Options Database Keys for Option Database: 1141c5b5d8d5SVaclav Hapla + -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc 1142c5b5d8d5SVaclav Hapla . -options_monitor - monitor all set options to standard output for the whole program run 1143c5b5d8d5SVaclav Hapla - -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet() 1144c5b5d8d5SVaclav Hapla 1145c5b5d8d5SVaclav Hapla Options -options_monitor_{all,cancel} are 1146c5b5d8d5SVaclav Hapla position-independent and apply to all options set since the PETSc start. 1147c5b5d8d5SVaclav Hapla They can be used also in option files. 1148c5b5d8d5SVaclav Hapla 1149c5b5d8d5SVaclav Hapla See PetscOptionsMonitorSet() to do monitoring programmatically. 1150c5b5d8d5SVaclav Hapla 1151e5c89e4eSSatish Balay Options Database Keys for Profiling: 1152a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 1153fe9b927eSVaclav Hapla + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo(). 1154217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 1155217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 1156495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 1157e5c89e4eSSatish Balay hangs without running in the debugger). See PetscLogTraceBegin(). 11589a9a5d4cSBarry Smith . -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView(). 1159156b51fbSBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each event, see PetscLogView(). 1160156b51fbSBarry Smith . -log_view_gpu_time - Includes in the summary from -log_view the time used in each GPU kernel, see PetscLogView(). 11619a9a5d4cSBarry Smith . -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the 1162495fc317SBarry Smith summary is written to the file. See PetscLogView(). 1163f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 1164495fc317SBarry Smith . -log_all [filename] - Logs extensive profiling information See PetscLogDump(). 1165495fc317SBarry Smith . -log [filename] - Logs basic profiline information See PetscLogDump(). 1166c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 116787c3beb6SLisandro Dalcin . -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off 1168c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 1169495fc317SBarry Smith 11706a77f485SPierre Jolivet Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time 1171e5c89e4eSSatish Balay 1172ffbd1cfbSBarry Smith Options Database Keys for SAWs: 1173ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 1174ffbd1cfbSBarry 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 1175ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 1176ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 1177ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 1178ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 1179ffbd1cfbSBarry Smith 1180e5c89e4eSSatish Balay Environmental Variables: 1181e5c89e4eSSatish Balay + PETSC_TMP - alternative tmp directory 1182e5c89e4eSSatish Balay . PETSC_SHARED_TMP - tmp is shared by all processes 1183e5c89e4eSSatish Balay . PETSC_NOT_SHARED_TMP - each process has its own private tmp 11844a971ea4SToby Isaac . PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs 11854a971ea4SToby Isaac . PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document 1186e5c89e4eSSatish Balay . PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer 1187e5c89e4eSSatish Balay - PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to 1188e5c89e4eSSatish Balay 1189e5c89e4eSSatish Balay Level: beginner 1190e5c89e4eSSatish Balay 1191e5c89e4eSSatish Balay Notes: 1192e5c89e4eSSatish Balay If for some reason you must call MPI_Init() separately, call 1193e5c89e4eSSatish Balay it before PetscInitialize(). 1194e5c89e4eSSatish Balay 1195e5c89e4eSSatish Balay Fortran Version: 1196e5c89e4eSSatish Balay In Fortran this routine has the format 1197e5c89e4eSSatish Balay $ call PetscInitialize(file,ierr) 1198e5c89e4eSSatish Balay 1199e5c89e4eSSatish Balay + ierr - error return code 1200c5b5d8d5SVaclav Hapla - file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc. 1201c5b5d8d5SVaclav Hapla Use PETSC_NULL_CHARACTER to not check for code specific file. 1202c5b5d8d5SVaclav Hapla Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files. 1203e5c89e4eSSatish Balay 1204e5c89e4eSSatish Balay Important Fortran Note: 12050eb4c9c0SBarry Smith In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a 12060298fd71SBarry Smith null character string; you CANNOT just use NULL as 1207a7f22e61SSatish Balay in the C version. See Users-Manual: ch_fortran for details. 1208e5c89e4eSSatish Balay 120901cb0274SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after 121001cb0274SBarry Smith calling PetscInitialize(). 1211e5c89e4eSSatish Balay 1212db781477SPatrick Sanan .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`, `PetscLogGpuTime()` 1213e5c89e4eSSatish Balay 1214e5c89e4eSSatish Balay @*/ 12159371c9d4SSatish Balay PetscErrorCode PetscInitialize(int *argc, char ***args, const char file[], const char help[]) { 121685649d77SJunchao Zhang PetscMPIInt flag; 1217*21ef0414SBarry Smith const char *prog = "Unknown Name", *mpienv; 1218e5c89e4eSSatish Balay 121927104ee2SJacob Faibussowitsch PetscFunctionBegin; 122027104ee2SJacob Faibussowitsch if (PetscInitializeCalled) PetscFunctionReturn(0); 12219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Initialized(&flag)); 1222e5c89e4eSSatish Balay if (!flag) { 122339a651e2SJacob 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"); 12249566063dSJacob Faibussowitsch PetscCall(PetscPreMPIInit_Private()); 12255e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 12265e765c61SJed Brown { 122739a651e2SJacob Faibussowitsch PetscMPIInt PETSC_UNUSED provided; 12289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init_thread(argc, args, PETSC_MPI_THREAD_REQUIRED, &provided)); 12295e765c61SJed Brown } 12305e765c61SJed Brown #else 12319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Init(argc, args)); 12325e765c61SJed Brown #endif 1233*21ef0414SBarry Smith if (PetscDefined(HAVE_MPIUNI)) { 1234*21ef0414SBarry Smith mpienv = getenv("PMI_SIZE"); 1235*21ef0414SBarry Smith if (!mpienv) mpienv = getenv("OMPI_COMM_WORLD_SIZE"); 1236*21ef0414SBarry Smith if (mpienv) { 1237*21ef0414SBarry Smith PetscInt isize; 1238*21ef0414SBarry Smith PetscCall(PetscOptionsStringToInt(mpienv, &isize)); 1239*21ef0414SBarry 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"); 1240*21ef0414SBarry 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"); 1241*21ef0414SBarry Smith } 1242*21ef0414SBarry Smith } 1243e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 1244e5c89e4eSSatish Balay } 12454dfee713SSatish Balay 124685649d77SJunchao Zhang if (argc && *argc) prog = **args; 1247e5c89e4eSSatish Balay if (argc && args) { 1248e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 1249e5c89e4eSSatish Balay PetscGlobalArgs = *args; 1250e5c89e4eSSatish Balay } 12519566063dSJacob Faibussowitsch PetscCall(PetscInitialize_Common(prog, file, help, PETSC_FALSE /*C*/, PETSC_FALSE, 0)); 125227104ee2SJacob Faibussowitsch PetscFunctionReturn(0); 1253e5c89e4eSSatish Balay } 1254e5c89e4eSSatish Balay 1255a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 125695c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1257ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1258ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 125995c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 12604097062eSBarry Smith #endif 1261e5c89e4eSSatish Balay 1262008a6e76SBarry Smith /* 1263008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1264008a6e76SBarry Smith */ 12659371c9d4SSatish Balay PetscErrorCode PetscFreeMPIResources(void) { 1266008a6e76SBarry Smith PetscFunctionBegin; 1267613bf2b2SPierre Jolivet #if defined(PETSC_HAVE_REAL___FLOAT128) 12689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FLOAT128)); 12699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___COMPLEX128)); 1270613bf2b2SPierre Jolivet #if !defined(PETSC_USE_REAL___FLOAT128) 1271613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_SUM___FLOAT128)); 1272008a6e76SBarry Smith #endif 1273613bf2b2SPierre Jolivet #endif 1274613bf2b2SPierre Jolivet #if defined(PETSC_USE_REAL___FP16) 12759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU___FP16)); 1276008a6e76SBarry Smith #endif 1277008a6e76SBarry Smith 1278de272c7aSSatish Balay #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 12799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_SUM)); 1280613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MAX)); 1281613bf2b2SPierre Jolivet PetscCallMPI(MPI_Op_free(&MPIU_MIN)); 1282008a6e76SBarry Smith #endif 1283008a6e76SBarry Smith 12849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2SCALAR)); 12859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_REAL_INT)); 12869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_SCALAR_INT)); 128740df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 12889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_2INT)); 1289008a6e76SBarry Smith #endif 12909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPI_4INT)); 12919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&MPIU_4INT)); 12929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Op_free(&MPIU_MAXSUM_OP)); 1293008a6e76SBarry Smith PetscFunctionReturn(0); 1294008a6e76SBarry Smith } 1295008a6e76SBarry Smith 1296a4af0ceeSJacob Faibussowitsch #if PetscDefined(USE_LOG) 1297a4af0ceeSJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscLogFinalize(void); 1298a4af0ceeSJacob Faibussowitsch #endif 1299a4af0ceeSJacob Faibussowitsch 1300e5c89e4eSSatish Balay /*@C 1301e5c89e4eSSatish Balay PetscFinalize - Checks for options to be called at the conclusion 1302e5c89e4eSSatish Balay of the program. MPI_Finalize() is called only if the user had not 1303e5c89e4eSSatish Balay called MPI_Init() before calling PetscInitialize(). 1304e5c89e4eSSatish Balay 1305e5c89e4eSSatish Balay Collective on PETSC_COMM_WORLD 1306e5c89e4eSSatish Balay 1307e5c89e4eSSatish Balay Options Database Keys: 130826a7e8d4SBarry Smith + -options_view - Calls PetscOptionsView() 1309e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 13107eb1d149SBarry 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 1311e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 131292f119d6SBarry Smith . -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed 1313e5c89e4eSSatish Balay . -malloc_info - Prints total memory usage 131492f119d6SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and where 1315e5c89e4eSSatish Balay 1316e5c89e4eSSatish Balay Level: beginner 1317e5c89e4eSSatish Balay 1318e5c89e4eSSatish Balay Note: 1319e5c89e4eSSatish Balay See PetscInitialize() for more general runtime options. 1320e5c89e4eSSatish Balay 1321db781477SPatrick Sanan .seealso: `PetscInitialize()`, `PetscOptionsView()`, `PetscMallocDump()`, `PetscMPIDump()`, `PetscEnd()` 1322e5c89e4eSSatish Balay @*/ 13239371c9d4SSatish Balay PetscErrorCode PetscFinalize(void) { 13244bb5149bSJed Brown PetscMPIInt rank; 1325a8d2bbe5SBarry Smith PetscInt nopt; 13262bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE, flg3 = PETSC_FALSE; 1327dff31646SBarry Smith PetscBool flg; 132810463e74SBarry Smith #if defined(PETSC_USE_LOG) 132910463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 133010463e74SBarry Smith #endif 1331e5c89e4eSSatish Balay 13323cbbc5ffSBarry Smith PetscFunctionBegin; 133339a651e2SJacob Faibussowitsch PetscCheck(PetscInitializeCalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscInitialize() must be called before PetscFinalize()"); 13349566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "PetscFinalize() called\n")); 1335b022a5c1SBarry Smith 1336f1f2ae84SBarry Smith PetscCall(PetscOptionsHasName(NULL, NULL, "-mpi_linear_solver_server", &flg)); 1337f1f2ae84SBarry Smith if (PetscDefined(USE_SINGLE_LIBRARY) && flg) PetscCall(PCMPIServerEnd()); 1338f1f2ae84SBarry Smith 13399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1340a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 13419566063dSJacob Faibussowitsch PetscCall(adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE)); 13429566063dSJacob Faibussowitsch PetscCall(adios_finalize(rank)); 1343a56f64adSBarry Smith #endif 13449566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-citations", &flg)); 1345dff31646SBarry Smith if (flg) { 13461f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 13471f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 13481f817a21SBarry Smith 13499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-citations", filename, sizeof(filename), NULL)); 13509371c9d4SSatish Balay if (filename[0]) { PetscCall(PetscFOpen(PETSC_COMM_WORLD, filename, "w", &fd)); } 13519566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(PetscCitationsList, 1, &cits)); 1352dff31646SBarry Smith cits[0] = 0; 13539566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(PetscCitationsList, &cits)); 13549566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "If you publish results based on this computation please cite the following:\n")); 13559566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 13569566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%s", cits)); 13579566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "===========================================================================\n")); 13589566063dSJacob Faibussowitsch PetscCall(PetscFClose(PETSC_COMM_WORLD, fd)); 13599566063dSJacob Faibussowitsch PetscCall(PetscFree(cits)); 1360dff31646SBarry Smith } 13619566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&PetscCitationsList)); 1362dff31646SBarry Smith 1363c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER) 136404102261SBarry Smith /* TextBelt is run for testing purposes only, please do not use this feature often */ 136504102261SBarry Smith { 136604102261SBarry Smith PetscInt nmax = 2; 136704102261SBarry Smith char **buffs; 13689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &buffs)); 13699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-textbelt", buffs, &nmax, &flg1)); 137004102261SBarry Smith if (flg1) { 137128b400f6SJacob Faibussowitsch PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-textbelt requires either the phone number or number,\"message\""); 137204102261SBarry Smith if (nmax == 1) { 13739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(128, &buffs[1])); 13749566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(buffs[1], 32)); 13759566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buffs[1], " has completed")); 137604102261SBarry Smith } 13779566063dSJacob Faibussowitsch PetscCall(PetscTextBelt(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL)); 13789566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[0])); 13799566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[1])); 138004102261SBarry Smith } 13819566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs)); 138204102261SBarry Smith } 1383763ec7b1SBarry Smith { 1384763ec7b1SBarry Smith PetscInt nmax = 2; 1385763ec7b1SBarry Smith char **buffs; 13869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &buffs)); 13879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-tellmycell", buffs, &nmax, &flg1)); 1388763ec7b1SBarry Smith if (flg1) { 138928b400f6SJacob Faibussowitsch PetscCheck(nmax, PETSC_COMM_WORLD, PETSC_ERR_USER, "-tellmycell requires either the phone number or number,\"message\""); 1390763ec7b1SBarry Smith if (nmax == 1) { 13919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(128, &buffs[1])); 13929566063dSJacob Faibussowitsch PetscCall(PetscGetProgramName(buffs[1], 32)); 13939566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buffs[1], " has completed")); 1394763ec7b1SBarry Smith } 13959566063dSJacob Faibussowitsch PetscCall(PetscTellMyCell(PETSC_COMM_WORLD, buffs[0], buffs[1], NULL)); 13969566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[0])); 13979566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs[1])); 1398763ec7b1SBarry Smith } 13999566063dSJacob Faibussowitsch PetscCall(PetscFree(buffs)); 1400763ec7b1SBarry Smith } 140104102261SBarry Smith #endif 140204102261SBarry Smith 14032d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 14049566063dSJacob Faibussowitsch PetscCall(PetscFPTDestroy()); 14052d53ad75SBarry Smith #endif 14062d53ad75SBarry Smith 1407e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1408dff31646SBarry Smith flg = PETSC_FALSE; 14099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-saw_options", &flg, NULL)); 14101baa6e33SBarry Smith if (flg) PetscCall(PetscOptionsSAWsDestroy()); 1411d5649816SBarry Smith #endif 1412d5649816SBarry Smith 1413681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1414681455b2SBarry Smith flg1 = PETSC_FALSE; 14159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-x_virtual", &flg1, NULL)); 1416681455b2SBarry Smith if (flg1) { 1417681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 14189566063dSJacob Faibussowitsch PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, "pkill -9 Xvfb", "r", NULL)); 1419681455b2SBarry Smith } 1420681455b2SBarry Smith #endif 1421681455b2SBarry Smith 142267584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 14239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_info", &flg2, NULL)); 1424e5c89e4eSSatish Balay if (!flg2) { 142590d69ab7SBarry Smith flg2 = PETSC_FALSE; 14269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-memory_view", &flg2, NULL)); 1427e5c89e4eSSatish Balay } 14289371c9d4SSatish Balay if (flg2) { PetscCall(PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD, "Summary of Memory Usage in PETSc\n")); } 142967584ceeSBarry Smith #endif 1430e5c89e4eSSatish Balay 1431e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 143290d69ab7SBarry Smith flg1 = PETSC_FALSE; 14339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-get_total_flops", &flg1, NULL)); 1434e5c89e4eSSatish Balay if (flg1) { 1435e5c89e4eSSatish Balay PetscLogDouble flops = 0; 14369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&petsc_TotalFlops, &flops, 1, MPI_DOUBLE, MPI_SUM, 0, PETSC_COMM_WORLD)); 14379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total flops over all processors %g\n", flops)); 1438e5c89e4eSSatish Balay } 1439e5c89e4eSSatish Balay #endif 1440e5c89e4eSSatish Balay 1441e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 1442e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE) 1443e5c89e4eSSatish Balay mname[0] = 0; 14449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_mpe", mname, sizeof(mname), &flg1)); 1445e5c89e4eSSatish Balay if (flg1) { 14469566063dSJacob Faibussowitsch if (mname[0]) PetscCall(PetscLogMPEDump(mname)); 14479566063dSJacob Faibussowitsch else PetscCall(PetscLogMPEDump(0)); 1448e5c89e4eSSatish Balay } 1449e5c89e4eSSatish Balay #endif 1450356e5837SBarry Smith #endif 1451a297a907SKarl Rupp 1452dd710f27SBarry Smith /* 1453dd710f27SBarry Smith Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 1454dd710f27SBarry Smith */ 14559566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1456dd710f27SBarry Smith 1457356e5837SBarry Smith #if defined(PETSC_USE_LOG) 14589566063dSJacob Faibussowitsch PetscCall(PetscOptionsPushGetViewerOff(PETSC_FALSE)); 14599566063dSJacob Faibussowitsch PetscCall(PetscLogViewFromOptions()); 14609566063dSJacob Faibussowitsch PetscCall(PetscOptionsPopGetViewerOff()); 146187c3beb6SLisandro Dalcin 1462356e5837SBarry Smith mname[0] = 0; 14639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_summary", mname, sizeof(mname), &flg1)); 1464e5c89e4eSSatish Balay if (flg1) { 146591eabc43SBarry Smith PetscViewer viewer; 14669566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PETSC_COMM_WORLD, "\n\n WARNING: -log_summary is being deprecated; switch to -log_view\n\n\n")); 146791eabc43SBarry Smith if (mname[0]) { 14689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, mname, &viewer)); 14699566063dSJacob Faibussowitsch PetscCall(PetscLogView(viewer)); 14709566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 147133f85c2fSBarry Smith } else { 147233f85c2fSBarry Smith viewer = PETSC_VIEWER_STDOUT_WORLD; 14739566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT)); 14749566063dSJacob Faibussowitsch PetscCall(PetscLogView(viewer)); 14759566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 147633f85c2fSBarry Smith } 1477e5c89e4eSSatish Balay } 1478a297a907SKarl Rupp 1479dd710f27SBarry Smith /* 1480dd710f27SBarry Smith Free any objects created by the last block of code. 1481dd710f27SBarry Smith */ 14829566063dSJacob Faibussowitsch PetscCall(PetscObjectRegisterDestroyAll()); 1483dd710f27SBarry Smith 1484dd710f27SBarry Smith mname[0] = 0; 14859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_all", mname, sizeof(mname), &flg1)); 14869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log", mname, sizeof(mname), &flg2)); 14879566063dSJacob Faibussowitsch if (flg1 || flg2) PetscCall(PetscLogDump(mname)); 1488e5c89e4eSSatish Balay #endif 148910463e74SBarry Smith 149090d69ab7SBarry Smith flg1 = PETSC_FALSE; 14919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_signal_handler", &flg1, NULL)); 14929566063dSJacob Faibussowitsch if (!flg1) PetscCall(PetscPopSignalHandler()); 149390d69ab7SBarry Smith flg1 = PETSC_FALSE; 14949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-mpidump", &flg1, NULL)); 14951baa6e33SBarry Smith if (flg1) PetscCall(PetscMPIDump(stdout)); 149690d69ab7SBarry Smith flg1 = PETSC_FALSE; 149790d69ab7SBarry Smith flg2 = PETSC_FALSE; 14988bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 14999566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-malloc_dump", &flg1)); 15009566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 15019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_view", &flg2, NULL)); 1502e4c476e2SSatish Balay 1503e5c89e4eSSatish Balay if (flg2) { 1504be56827dSJed Brown PetscViewer viewer; 15059566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 15069566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 15079566063dSJacob Faibussowitsch PetscCall(PetscOptionsView(NULL, viewer)); 15089566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1509e5c89e4eSSatish Balay } 1510e5c89e4eSSatish Balay 1511e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 15129566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox", &flg1)); 15139566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-nox_warning", &flg1)); 1514e5c89e4eSSatish Balay 151533fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 15169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-options_left", &flg3, &flg1)); 15179245e749SBarry Smith if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE; 1518e5c89e4eSSatish Balay if (flg3) { 15193de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 1520be56827dSJed Brown PetscViewer viewer; 15219566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 15229566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 15239566063dSJacob Faibussowitsch PetscCall(PetscOptionsView(NULL, viewer)); 15249566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1525e5c89e4eSSatish Balay } 15269566063dSJacob Faibussowitsch PetscCall(PetscOptionsAllUsed(NULL, &nopt)); 15273de2bfdfSBarry Smith if (nopt) { 15289566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! There are options you set that were not used!\n")); 15299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING! could be spelling mistake, etc!\n")); 15303de2bfdfSBarry Smith if (nopt == 1) { 15319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There is one unused database option. It is:\n")); 1532e5c89e4eSSatish Balay } else { 15339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are %" PetscInt_FMT " unused database options. They are:\n", nopt)); 1534e5c89e4eSSatish Balay } 15353de2bfdfSBarry Smith } else if (flg3 && flg1) { 15369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "There are no unused options.\n")); 1537df12ba86SBarry Smith } 15389566063dSJacob Faibussowitsch PetscCall(PetscOptionsLeft(NULL)); 1539e5c89e4eSSatish Balay } 1540e5c89e4eSSatish Balay 1541e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1542d45a07a7SBarry Smith if (!PetscGlobalRank) { 15439566063dSJacob Faibussowitsch PetscCall(PetscStackSAWsViewOff()); 1544792fecdfSBarry Smith PetscCallSAWs(SAWs_Finalize, ()); 1545d45a07a7SBarry Smith } 1546ec957eceSBarry Smith #endif 1547ec957eceSBarry Smith 15484097062eSBarry Smith #if defined(PETSC_USE_LOG) 154910463e74SBarry Smith /* 1550dbc8283eSBarry Smith List all objects the user may have forgot to free 15512eff7a51SBarry Smith */ 155205df10baSBarry Smith if (PetscObjectsLog) { 15539566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-objects_dump", &flg1)); 1554a64a8e02SBarry Smith if (flg1) { 1555a64a8e02SBarry Smith MPI_Comm local_comm; 15567eb1d149SBarry Smith char string[64]; 1557a64a8e02SBarry Smith 15589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-objects_dump", string, sizeof(string), NULL)); 1559afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 15609566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 15619566063dSJacob Faibussowitsch PetscCall(PetscObjectsDump(stdout, (string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE)); 15629566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 15639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 15640a1571b3SBarry Smith } 156505df10baSBarry Smith } 15664097062eSBarry Smith #endif 15674097062eSBarry Smith 15684097062eSBarry Smith #if defined(PETSC_USE_LOG) 1569dbc8283eSBarry Smith PetscObjectsCounts = 0; 1570dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 15719566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscObjects)); 15724097062eSBarry Smith #endif 15732eff7a51SBarry Smith 157433f85c2fSBarry Smith /* 157533f85c2fSBarry Smith Destroy any packages that registered a finalize 157633f85c2fSBarry Smith */ 15779566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalizeAll()); 157833f85c2fSBarry Smith 1579101409b8SToby Isaac #if defined(PETSC_USE_LOG) 15809566063dSJacob Faibussowitsch PetscCall(PetscLogFinalize()); 1581101409b8SToby Isaac #endif 1582101409b8SToby Isaac 158333f85c2fSBarry Smith /* 158448dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 158548dd1dffSBarry Smith */ 15862e956fe4SStefano Zampini if (PetscPrintFunctionList) PetscCall(PetscFunctionListPrintAll()); 158737e93019SBarry Smith 15884028d114SSatish Balay if (petsc_history) { 15899566063dSJacob Faibussowitsch PetscCall(PetscCloseHistoryFile(&petsc_history)); 159002c9f0b5SLisandro Dalcin petsc_history = NULL; 1591e5c89e4eSSatish Balay } 15929566063dSJacob Faibussowitsch PetscCall(PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton)); 15939566063dSJacob Faibussowitsch PetscCall(PetscInfoDestroy()); 1594e5c89e4eSSatish Balay 159567584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 159692f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1597e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 159892f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1599e5c89e4eSSatish Balay FILE *fd; 1600ed9cf6e9SBarry Smith int err; 1601e5c89e4eSSatish Balay 1602dc92acbaSJed Brown flg2 = PETSC_FALSE; 160392f119d6SBarry Smith flg3 = PETSC_FALSE; 16049566063dSJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_test", &flg2, NULL)); 16059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-malloc_debug", &flg3, NULL)); 160692f119d6SBarry Smith fname[0] = 0; 16079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_dump", fname, sizeof(fname), &flg1)); 1608e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1609589a23caSBarry Smith PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank); 16109371c9d4SSatish Balay fd = fopen(sname, "w"); 16119371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16129566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(fd)); 1613ed9cf6e9SBarry Smith err = fclose(fd); 161428b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 161592f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1616e5c89e4eSSatish Balay MPI_Comm local_comm; 1617e5c89e4eSSatish Balay 1618afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16199566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16209566063dSJacob Faibussowitsch PetscCall(PetscMallocDump(stdout)); 16219566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1623e5c89e4eSSatish Balay } 1624e5c89e4eSSatish Balay fname[0] = 0; 16259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-malloc_view", fname, sizeof(fname), &flg1)); 1626e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1627589a23caSBarry Smith PetscSNPrintf(sname, sizeof(sname), "%s_%d", fname, rank); 16289371c9d4SSatish Balay fd = fopen(sname, "w"); 16299371c9d4SSatish Balay PetscCheck(fd, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", sname); 16309566063dSJacob Faibussowitsch PetscCall(PetscMallocView(fd)); 1631ed9cf6e9SBarry Smith err = fclose(fd); 163228b400f6SJacob Faibussowitsch PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file"); 163392f119d6SBarry Smith } else if (flg1) { 163492f119d6SBarry Smith MPI_Comm local_comm; 163592f119d6SBarry Smith 1636afd7ed4bSMatthew G. Knepley PetscCallMPI(MPI_Comm_dup(PETSC_COMM_WORLD, &local_comm)); 16379566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseBegin_Private(local_comm, 1)); 16389566063dSJacob Faibussowitsch PetscCall(PetscMallocView(stdout)); 16399566063dSJacob Faibussowitsch PetscCall(PetscSequentialPhaseEnd_Private(local_comm, 1)); 16409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&local_comm)); 1641e5c89e4eSSatish Balay } 1642e5c89e4eSSatish Balay } 164367584ceeSBarry Smith #endif 164420e2c332SMatthew G. Knepley 16455486ca60SMatthew G. Knepley /* 16465486ca60SMatthew G. Knepley Close any open dynamic libraries 16475486ca60SMatthew G. Knepley */ 16489566063dSJacob Faibussowitsch PetscCall(PetscFinalize_DynamicLibraries()); 16495486ca60SMatthew G. Knepley 1650e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 16519566063dSJacob Faibussowitsch PetscCall(PetscOptionsDestroyDefault()); 1652e5c89e4eSSatish Balay 1653e5c89e4eSSatish Balay PetscGlobalArgc = 0; 165402c9f0b5SLisandro Dalcin PetscGlobalArgs = NULL; 1655e5c89e4eSSatish Balay 1656c2b86a48SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 1657c2b86a48SJunchao Zhang if (PetscBeganKokkos) { 16589566063dSJacob Faibussowitsch PetscCall(PetscKokkosFinalize_Private()); 1659c2b86a48SJunchao Zhang PetscBeganKokkos = PETSC_FALSE; 166045639126SStefano Zampini PetscKokkosInitialized = PETSC_FALSE; 1661c2b86a48SJunchao Zhang } 1662c2b86a48SJunchao Zhang #endif 1663c2b86a48SJunchao Zhang 166471438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM) 166571438e86SJunchao Zhang if (PetscBeganNvshmem) { 16669566063dSJacob Faibussowitsch PetscCall(PetscNvshmemFinalize()); 166771438e86SJunchao Zhang PetscBeganNvshmem = PETSC_FALSE; 166871438e86SJunchao Zhang } 166971438e86SJunchao Zhang #endif 1670a0e72f99SJunchao Zhang 16719566063dSJacob Faibussowitsch PetscCall(PetscFreeMPIResources()); 1672e5c89e4eSSatish Balay 1673dbc8283eSBarry Smith /* 1674efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1675efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1676efb80d3cSBarry Smith 1677efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1678efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1679dbc8283eSBarry Smith */ 1680b770b1f6SSatish Balay { 1681dbc8283eSBarry Smith PetscCommCounter *counter; 1682dbc8283eSBarry Smith PetscMPIInt flg; 1683dbc8283eSBarry Smith MPI_Comm icomm; 16849371c9d4SSatish Balay union 16859371c9d4SSatish Balay { 16869371c9d4SSatish Balay MPI_Comm comm; 16879371c9d4SSatish Balay void *ptr; 16889371c9d4SSatish Balay } ucomm; 16899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval, &ucomm, &flg)); 1690dbc8283eSBarry Smith if (flg) { 1691265f3f35SJed Brown icomm = ucomm.comm; 16929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 169328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1694dbc8283eSBarry Smith 16959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_SELF, Petsc_InnerComm_keyval)); 16969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 16979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1698dbc8283eSBarry Smith } 16999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval, &ucomm, &flg)); 1700dbc8283eSBarry Smith if (flg) { 1701265f3f35SJed Brown icomm = ucomm.comm; 17029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_get_attr(icomm, Petsc_Counter_keyval, &counter, &flg)); 170328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_CORRUPT, "Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1704dbc8283eSBarry Smith 17059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(PETSC_COMM_WORLD, Petsc_InnerComm_keyval)); 17069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_delete_attr(icomm, Petsc_Counter_keyval)); 17079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&icomm)); 1708dbc8283eSBarry Smith } 1709b770b1f6SSatish Balay } 1710dbc8283eSBarry Smith 17119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_Counter_keyval)); 17129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_InnerComm_keyval)); 17139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_OuterComm_keyval)); 17149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free_keyval(&Petsc_ShmComm_keyval)); 1715480cf27aSJed Brown 17169566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen)); 17179566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout)); 17189566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr)); 17199566063dSJacob Faibussowitsch PetscCall(PetscSpinlockDestroy(&PetscCommSpinLock)); 1720ef19f930SBarry Smith 1721e5c89e4eSSatish Balay if (PetscBeganMPI) { 172299b1327fSBarry Smith PetscMPIInt flag; 17239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalized(&flag)); 172439a651e2SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 172539a651e2SJacob Faibussowitsch /* wait until the very last moment to disable error handling */ 172639a651e2SJacob Faibussowitsch PetscErrorHandlingInitialized = PETSC_FALSE; 17279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Finalize()); 172839a651e2SJacob Faibussowitsch } else PetscErrorHandlingInitialized = PETSC_FALSE; 172939a651e2SJacob Faibussowitsch 1730e5c89e4eSSatish Balay /* 1731e5c89e4eSSatish Balay 1732e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1733e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1734e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1735e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1736e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 17370e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1738e5c89e4eSSatish Balay memory was not freed. 1739e5c89e4eSSatish Balay 1740e5c89e4eSSatish Balay */ 17419566063dSJacob Faibussowitsch PetscCall(PetscMallocClear()); 17429566063dSJacob Faibussowitsch PetscCall(PetscStackReset()); 1743a297a907SKarl Rupp 1744e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1745e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 174656883f60SBarry Smith #if defined(PETSC_USE_GCOV) 174756883f60SBarry Smith /* 174856883f60SBarry Smith flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the 174956883f60SBarry Smith gcov files are still being added to the directories as git tries to remove the directories. 175056883f60SBarry Smith */ 175156883f60SBarry Smith __gcov_flush(); 175256883f60SBarry Smith #endif 17531724198aSStefano Zampini /* To match PetscFunctionBegin() at the beginning of this function */ 17541724198aSStefano Zampini PetscStackClearTop; 175527104ee2SJacob Faibussowitsch return 0; 1756e5c89e4eSSatish Balay } 1757e5c89e4eSSatish Balay 175843db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 17599371c9d4SSatish Balay PETSC_EXTERN int lsame_(char *a, char *b) { 176043db4dbbSBarry Smith if (*a == *b) return 1; 176143db4dbbSBarry Smith if (*a + 32 == *b) return 1; 176243db4dbbSBarry Smith if (*a - 32 == *b) return 1; 176343db4dbbSBarry Smith return 0; 176443db4dbbSBarry Smith } 1765a70650f6SBarry Smith #endif 176643db4dbbSBarry Smith 176743db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 17689371c9d4SSatish Balay PETSC_EXTERN int lsame(char *a, char *b) { 176943db4dbbSBarry Smith if (*a == *b) return 1; 177043db4dbbSBarry Smith if (*a + 32 == *b) return 1; 177143db4dbbSBarry Smith if (*a - 32 == *b) return 1; 177243db4dbbSBarry Smith return 0; 177343db4dbbSBarry Smith } 177443db4dbbSBarry Smith #endif 1775