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*/ 6022afb99SBarry Smith #include <petscvalgrind.h> 7665c2dedSJed Brown #include <petscviewer.h> 88101f56cSMatthew Knepley 9a9f03627SSatish Balay #if defined(PETSC_USE_LOG) 1095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscLogFinalize(void); 11a9f03627SSatish Balay #endif 12f2d66bcaSShri Abhyankar 132d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 1495c0884eSLisandro Dalcin PETSC_INTERN PetscFPT PetscFPTData; 152d53ad75SBarry Smith PetscFPT PetscFPTData = 0; 162d53ad75SBarry Smith #endif 172d53ad75SBarry Smith 18a6790183SBarry Smith #if defined(PETSC_HAVE_SAWS) 1916ad0300SBarry Smith #include <petscviewersaws.h> 20a6790183SBarry Smith #endif 21eb02082bSJunchao Zhang 22e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/ 23e5c89e4eSSatish Balay 2495c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history; 25e5c89e4eSSatish Balay 2695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void); 2795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void); 2895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void); 2995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int); 3095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int); 3195c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**); 320069ddf5SShri Abhyankar 33e5c89e4eSSatish Balay /* user may set this BEFORE calling PetscInitialize() */ 34e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL; 35e5c89e4eSSatish Balay 36480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval = MPI_KEYVAL_INVALID; 37480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID; 38480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID; 395f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval = MPI_KEYVAL_INVALID; 40480cf27aSJed Brown 41e5c89e4eSSatish Balay /* 42e5c89e4eSSatish Balay Declare and set all the string names of the PETSc enums 43e5c89e4eSSatish Balay */ 4402c9f0b5SLisandro Dalcin const char *const PetscBools[] = {"FALSE","TRUE","PetscBool","PETSC_",NULL}; 4502c9f0b5SLisandro Dalcin const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL}; 46e5c89e4eSSatish Balay 47ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE; 48ace3abfcSBarry Smith PetscBool PetscPreLoadingOn = PETSC_FALSE; 490f8e0872SSatish Balay 50a2f94806SJed Brown PetscInt PetscHotRegionDepth; 51a2f94806SJed Brown 52b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY) 53b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen; 54b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout; 55b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr; 56b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock; 57b22622e2STadeu Manoel #endif 58b22622e2STadeu Manoel 59e5c89e4eSSatish Balay /* 60945d1669SBarry Smith PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args 6172a42c3cSBarry Smith 6272a42c3cSBarry Smith Collective 6372a42c3cSBarry Smith 6472a42c3cSBarry Smith Level: advanced 6572a42c3cSBarry Smith 6695452b02SPatrick Sanan Notes: 67a56f64adSBarry Smith this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to 680f11a792SBarry Smith indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to 69a56f64adSBarry Smith be called multiple times from Julia without the problem of trying to initialize MPI more than once. 700f11a792SBarry Smith 71a56f64adSBarry Smith Developer Note: Turns off PETSc signal handling to allow Julia to manage signals 721ea5a559SBarry Smith 7372a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments() 740f11a792SBarry Smith */ 75945d1669SBarry Smith PetscErrorCode PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help) 7672a42c3cSBarry Smith { 7772a42c3cSBarry Smith PetscErrorCode ierr; 7872a42c3cSBarry Smith int myargc = argc; 7972a42c3cSBarry Smith char **myargs = args; 8072a42c3cSBarry Smith 8172a42c3cSBarry Smith PetscFunctionBegin; 82c52ac268SRichard Tran Mills ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr; 831ea5a559SBarry Smith ierr = PetscPopSignalHandler();CHKERRQ(ierr); 84df413903SBarry Smith PetscBeganMPI = PETSC_FALSE; 8572a42c3cSBarry Smith PetscFunctionReturn(ierr); 8672a42c3cSBarry Smith } 8772a42c3cSBarry Smith 88f0865b08SBarry Smith /* 89a56f64adSBarry Smith Used by Julia interface to get communicator 90f0865b08SBarry Smith */ 91945d1669SBarry Smith PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm) 92f0865b08SBarry Smith { 93f0865b08SBarry Smith PetscFunctionBegin; 94f0865b08SBarry Smith *comm = PETSC_COMM_SELF; 95f0865b08SBarry Smith PetscFunctionReturn(0); 96f0865b08SBarry Smith } 97f0865b08SBarry Smith 98e5c89e4eSSatish Balay /*@C 99e5c89e4eSSatish Balay PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without 100e5c89e4eSSatish Balay the command line arguments. 101e5c89e4eSSatish Balay 102e5c89e4eSSatish Balay Collective 103e5c89e4eSSatish Balay 104e5c89e4eSSatish Balay Level: advanced 105e5c89e4eSSatish Balay 106e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran() 107e5c89e4eSSatish Balay @*/ 1087087cfbeSBarry Smith PetscErrorCode PetscInitializeNoArguments(void) 109e5c89e4eSSatish Balay { 110e5c89e4eSSatish Balay PetscErrorCode ierr; 111e5c89e4eSSatish Balay int argc = 0; 11202c9f0b5SLisandro Dalcin char **args = NULL; 113e5c89e4eSSatish Balay 114e5c89e4eSSatish Balay PetscFunctionBegin; 1150298fd71SBarry Smith ierr = PetscInitialize(&argc,&args,NULL,NULL); 116e5c89e4eSSatish Balay PetscFunctionReturn(ierr); 117e5c89e4eSSatish Balay } 118e5c89e4eSSatish Balay 119e5c89e4eSSatish Balay /*@ 120e5c89e4eSSatish Balay PetscInitialized - Determine whether PETSc is initialized. 121e5c89e4eSSatish Balay 12293b6d2d1SJed Brown Level: beginner 123e5c89e4eSSatish Balay 124e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran() 125e5c89e4eSSatish Balay @*/ 1267087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized) 127e5c89e4eSSatish Balay { 128e5c89e4eSSatish Balay *isInitialized = PetscInitializeCalled; 12993b6d2d1SJed Brown return 0; 130e5c89e4eSSatish Balay } 131e5c89e4eSSatish Balay 132e5c89e4eSSatish Balay /*@ 133e5c89e4eSSatish Balay PetscFinalized - Determine whether PetscFinalize() has been called yet 134e5c89e4eSSatish Balay 135e5c89e4eSSatish Balay Level: developer 136e5c89e4eSSatish Balay 137e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran() 138e5c89e4eSSatish Balay @*/ 1397087cfbeSBarry Smith PetscErrorCode PetscFinalized(PetscBool *isFinalized) 140e5c89e4eSSatish Balay { 141e5c89e4eSSatish Balay *isFinalized = PetscFinalizeCalled; 14293b6d2d1SJed Brown return 0; 143e5c89e4eSSatish Balay } 144e5c89e4eSSatish Balay 14595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void); 146e5c89e4eSSatish Balay 147e5c89e4eSSatish Balay /* 148e5c89e4eSSatish Balay This function is the MPI reduction operation used to compute the sum of the 149e5c89e4eSSatish Balay first half of the datatype and the max of the second half. 150e5c89e4eSSatish Balay */ 151367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0; 152e5c89e4eSSatish Balay 153367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype) 154e5c89e4eSSatish Balay { 155e5c89e4eSSatish Balay PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt; 156e5c89e4eSSatish Balay 157e5c89e4eSSatish Balay PetscFunctionBegin; 158e5c89e4eSSatish Balay if (*datatype != MPIU_2INT) { 159e5c89e4eSSatish Balay (*PetscErrorPrintf)("Can only handle MPIU_2INT data types"); 16041e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 161e5c89e4eSSatish Balay } 162e5c89e4eSSatish Balay 163e5c89e4eSSatish Balay for (i=0; i<count; i++) { 164e5c89e4eSSatish Balay xout[2*i] = PetscMax(xout[2*i],xin[2*i]); 165e5c89e4eSSatish Balay xout[2*i+1] += xin[2*i+1]; 166e5c89e4eSSatish Balay } 167812af9f3SBarry Smith PetscFunctionReturnVoid(); 168e5c89e4eSSatish Balay } 169e5c89e4eSSatish Balay 170e5c89e4eSSatish Balay /* 171e5c89e4eSSatish Balay Returns the max of the first entry owned by this processor and the 172e5c89e4eSSatish Balay sum of the second entry. 173b693b147SBarry Smith 17476f543a4SJed Brown The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero 175367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths 176b693b147SBarry Smith there would be no place to store the both needed results. 177e5c89e4eSSatish Balay */ 17876ec1555SBarry Smith PetscErrorCode PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum) 179e5c89e4eSSatish Balay { 180e5c89e4eSSatish Balay PetscErrorCode ierr; 181e5c89e4eSSatish Balay 182e5c89e4eSSatish Balay PetscFunctionBegin; 183d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 184d6e4c47cSJed Brown { 185d6e4c47cSJed Brown struct {PetscInt max,sum;} work; 186367daffbSBarry Smith ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr); 187d6e4c47cSJed Brown *max = work.max; 188d6e4c47cSJed Brown *sum = work.sum; 189d6e4c47cSJed Brown } 190d6e4c47cSJed Brown #else 191d6e4c47cSJed Brown { 192d6e4c47cSJed Brown PetscMPIInt size,rank; 193d6e4c47cSJed Brown struct {PetscInt max,sum;} *work; 194e5c89e4eSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 195e5c89e4eSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 196785e854fSJed Brown ierr = PetscMalloc1(size,&work);CHKERRQ(ierr); 197367daffbSBarry Smith ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr); 1986ac3741eSJed Brown *max = work[rank].max; 1996ac3741eSJed Brown *sum = work[rank].sum; 200e5c89e4eSSatish Balay ierr = PetscFree(work);CHKERRQ(ierr); 201d6e4c47cSJed Brown } 202d6e4c47cSJed Brown #endif 203e5c89e4eSSatish Balay PetscFunctionReturn(0); 204e5c89e4eSSatish Balay } 205e5c89e4eSSatish Balay 206e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/ 207e5c89e4eSSatish Balay 208570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 20906a205a8SBarry Smith MPI_Op MPIU_SUM = 0; 210e5c89e4eSSatish Balay 2118cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 212e5c89e4eSSatish Balay { 213e5c89e4eSSatish Balay PetscInt i,count = *cnt; 214e5c89e4eSSatish Balay 215e5c89e4eSSatish Balay PetscFunctionBegin; 2167c2de775SJed Brown if (*datatype == MPIU_REAL) { 217e2e03761SBarry Smith PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 218a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] += xin[i]; 2197c2de775SJed Brown } 2207c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2217c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2227c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 223a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] += xin[i]; 2247c2de775SJed Brown } 2257c2de775SJed Brown #endif 2267c2de775SJed Brown else { 2277c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 22841e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 229e2e03761SBarry Smith } 230812af9f3SBarry Smith PetscFunctionReturnVoid(); 231e5c89e4eSSatish Balay } 232e5c89e4eSSatish Balay #endif 233e5c89e4eSSatish Balay 234570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 235d9822059SBarry Smith MPI_Op MPIU_MAX = 0; 236d9822059SBarry Smith MPI_Op MPIU_MIN = 0; 237d9822059SBarry Smith 2388cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 239d9822059SBarry Smith { 240d9822059SBarry Smith PetscInt i,count = *cnt; 241d9822059SBarry Smith 242d9822059SBarry Smith PetscFunctionBegin; 2437c2de775SJed Brown if (*datatype == MPIU_REAL) { 2448c764dc5SJose Roman PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 245a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]); 2467c2de775SJed Brown } 2477c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2487c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2497c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 2507c2de775SJed Brown for (i=0; i<count; i++) { 2517c2de775SJed Brown xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 2527c2de775SJed Brown } 2537c2de775SJed Brown } 2547c2de775SJed Brown #endif 2557c2de775SJed Brown else { 2567c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 25741e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 2588c764dc5SJose Roman } 259d9822059SBarry Smith PetscFunctionReturnVoid(); 260d9822059SBarry Smith } 261d9822059SBarry Smith 2628cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 263d9822059SBarry Smith { 264d9822059SBarry Smith PetscInt i,count = *cnt; 265d9822059SBarry Smith 266d9822059SBarry Smith PetscFunctionBegin; 2677c2de775SJed Brown if (*datatype == MPIU_REAL) { 2688c764dc5SJose Roman PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 269a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]); 2707c2de775SJed Brown } 2717c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2727c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2737c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 2747c2de775SJed Brown for (i=0; i<count; i++) { 2757c2de775SJed Brown xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 2767c2de775SJed Brown } 2777c2de775SJed Brown } 2787c2de775SJed Brown #endif 2797c2de775SJed Brown else { 2808c764dc5SJose Roman (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"); 28141e02c4dSJunchao Zhang PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG); 2828c764dc5SJose Roman } 283d9822059SBarry Smith PetscFunctionReturnVoid(); 284d9822059SBarry Smith } 285d9822059SBarry Smith #endif 286d9822059SBarry Smith 287480cf27aSJed Brown /* 288480cf27aSJed Brown Private routine to delete internal tag/name counter storage when a communicator is freed. 289480cf27aSJed Brown 290ff0e51ddSBarry 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. 291480cf27aSJed Brown 29212801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 293480cf27aSJed Brown 294480cf27aSJed Brown */ 29533779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state) 296480cf27aSJed Brown { 297480cf27aSJed Brown PetscErrorCode ierr; 29805342407SJunchao Zhang PetscCommCounter *counter=(PetscCommCounter*)count_val; 299480cf27aSJed Brown 300480cf27aSJed Brown PetscFunctionBegin; 30102c9f0b5SLisandro Dalcin ierr = PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr); 30205342407SJunchao Zhang ierr = PetscFree(counter->iflags);CHKERRMPI(ierr); 30305342407SJunchao Zhang ierr = PetscFree(counter);CHKERRMPI(ierr); 304480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 305480cf27aSJed Brown } 306480cf27aSJed Brown 307480cf27aSJed Brown /* 30847435625SJed Brown This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user 309da3039f7SJed Brown calls MPI_Comm_free(). 310da3039f7SJed Brown 311da3039f7SJed Brown This is the only entry point for breaking the links between inner and outer comms. 312480cf27aSJed Brown 313ff0e51ddSBarry Smith This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator. 314480cf27aSJed Brown 31512801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 316480cf27aSJed Brown 317480cf27aSJed Brown */ 31833779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state) 319480cf27aSJed Brown { 320480cf27aSJed Brown PetscErrorCode ierr; 32133779a13SJunchao Zhang union {MPI_Comm comm; void *ptr;} icomm; 322480cf27aSJed Brown 323480cf27aSJed Brown PetscFunctionBegin; 32412801b39SBarry Smith if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval"); 325ec4fadc2SJed Brown icomm.ptr = attr_val; 32633779a13SJunchao Zhang #if defined(PETSC_USE_DEBUG) 32733779a13SJunchao Zhang { 32833779a13SJunchao Zhang /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */ 32933779a13SJunchao Zhang PetscMPIInt flg; 33033779a13SJunchao Zhang union {MPI_Comm comm; void *ptr;} ocomm; 33147435625SJed Brown ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr); 33233779a13SJunchao Zhang if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute"); 33333779a13SJunchao 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"); 33433779a13SJunchao Zhang } 33533779a13SJunchao Zhang #endif 33647435625SJed Brown ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr); 33733779a13SJunchao Zhang ierr = PetscInfo2(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);CHKERRMPI(ierr); 338da3039f7SJed Brown PetscFunctionReturn(MPI_SUCCESS); 339b89831e5SBarry Smith } 340da3039f7SJed Brown 341da3039f7SJed Brown /* 34233779a13SJunchao 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. 343da3039f7SJed Brown */ 34433779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state) 345da3039f7SJed Brown { 346da3039f7SJed Brown PetscErrorCode ierr; 347da3039f7SJed Brown 348da3039f7SJed Brown PetscFunctionBegin; 34902c9f0b5SLisandro Dalcin ierr = PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr); 350480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 351480cf27aSJed Brown } 352480cf27aSJed Brown 35333779a13SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *); 35442218b76SBarry Smith 355951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 3568cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*); 3578cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*); 3588cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*); 359e39fd77fSBarry Smith #endif 360e39fd77fSBarry Smith 3610f9be574SPeter Hill PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE; 36212801b39SBarry Smith 363eb27c7c8SSatish Balay PETSC_INTERN int PetscGlobalArgc; 364eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs; 3656ae9a8a6SBarry Smith int PetscGlobalArgc = 0; 36602c9f0b5SLisandro Dalcin char **PetscGlobalArgs = NULL; 367dff31646SBarry Smith PetscSegBuffer PetscCitationsList; 368e5c89e4eSSatish Balay 369dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void) 370051e4cf2SJed Brown { 371051e4cf2SJed Brown PetscErrorCode ierr; 372051e4cf2SJed Brown 373051e4cf2SJed Brown PetscFunctionBegin; 374051e4cf2SJed Brown ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr); 375a1601671SBarry Smith ierr = PetscCitationsRegister("@TechReport{petsc-user-ref,\n Author = {Satish Balay and Shrirang Abhyankar and Mark F. Adams and Jed Brown \n and Peter Brune and Kris Buschelman and Lisandro Dalcin and\n Victor Eijkhout and William D. Gropp and Dmitry Karpeyev and\n Dinesh Kaushik and Matthew G. Knepley and Dave A. May and Lois Curfman McInnes\n and Richard Tran Mills and Todd Munson and Karl Rupp and Patrick Sanan\n and Barry F. Smith and Stefano Zampini and Hong Zhang and Hong Zhang},\n Title = {{PETS}c Users Manual},\n Number = {ANL-95/11 - Revision 3.11},\n Institution = {Argonne National Laboratory},\n Year = {2019}\n}\n",NULL);CHKERRQ(ierr); 376051e4cf2SJed Brown ierr = PetscCitationsRegister("@InProceedings{petsc-efficient,\n Author = {Satish Balay and William D. Gropp and Lois Curfman McInnes and Barry F. Smith},\n Title = {Efficient Management of Parallelism in Object Oriented Numerical Software Libraries},\n Booktitle = {Modern Software Tools in Scientific Computing},\n Editor = {E. Arge and A. M. Bruaset and H. P. Langtangen},\n Pages = {163--202},\n Publisher = {Birkh{\\\"{a}}user Press},\n Year = {1997}\n}\n",NULL);CHKERRQ(ierr); 377051e4cf2SJed Brown PetscFunctionReturn(0); 378051e4cf2SJed Brown } 379e5c89e4eSSatish Balay 3802d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */ 3812d747510SLisandro Dalcin 3822d747510SLisandro Dalcin PetscErrorCode PetscSetProgramName(const char name[]) 3832d747510SLisandro Dalcin { 3842d747510SLisandro Dalcin PetscErrorCode ierr; 3852d747510SLisandro Dalcin 3862d747510SLisandro Dalcin PetscFunctionBegin; 3872d747510SLisandro Dalcin ierr = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 3882d747510SLisandro Dalcin PetscFunctionReturn(0); 3892d747510SLisandro Dalcin } 3902d747510SLisandro Dalcin 3912d747510SLisandro Dalcin /*@C 3922d747510SLisandro Dalcin PetscGetProgramName - Gets the name of the running program. 3932d747510SLisandro Dalcin 3942d747510SLisandro Dalcin Not Collective 3952d747510SLisandro Dalcin 3962d747510SLisandro Dalcin Input Parameter: 3972d747510SLisandro Dalcin . len - length of the string name 3982d747510SLisandro Dalcin 3992d747510SLisandro Dalcin Output Parameter: 4002d747510SLisandro Dalcin . name - the name of the running program 4012d747510SLisandro Dalcin 4022d747510SLisandro Dalcin Level: advanced 4032d747510SLisandro Dalcin 4042d747510SLisandro Dalcin Notes: 4052d747510SLisandro Dalcin The name of the program is copied into the user-provided character 4062d747510SLisandro Dalcin array of length len. On some machines the program name includes 4072d747510SLisandro Dalcin its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN. 4082d747510SLisandro Dalcin @*/ 4092d747510SLisandro Dalcin PetscErrorCode PetscGetProgramName(char name[],size_t len) 4102d747510SLisandro Dalcin { 4112d747510SLisandro Dalcin PetscErrorCode ierr; 4122d747510SLisandro Dalcin 4132d747510SLisandro Dalcin PetscFunctionBegin; 4142d747510SLisandro Dalcin ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr); 4152d747510SLisandro Dalcin PetscFunctionReturn(0); 4162d747510SLisandro Dalcin } 4172d747510SLisandro Dalcin 418e5c89e4eSSatish Balay /*@C 419e5c89e4eSSatish Balay PetscGetArgs - Allows you to access the raw command line arguments anywhere 420e5c89e4eSSatish Balay after PetscInitialize() is called but before PetscFinalize(). 421e5c89e4eSSatish Balay 422e5c89e4eSSatish Balay Not Collective 423e5c89e4eSSatish Balay 424e5c89e4eSSatish Balay Output Parameters: 425e5c89e4eSSatish Balay + argc - count of number of command line arguments 426e5c89e4eSSatish Balay - args - the command line arguments 427e5c89e4eSSatish Balay 428e5c89e4eSSatish Balay Level: intermediate 429e5c89e4eSSatish Balay 430e5c89e4eSSatish Balay Notes: 431e5c89e4eSSatish Balay This is usually used to pass the command line arguments into other libraries 432e5c89e4eSSatish Balay that are called internally deep in PETSc or the application. 433e5c89e4eSSatish Balay 434f177e3b1SBarry Smith The first argument contains the program name as is normal for C arguments. 435f177e3b1SBarry Smith 436793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments() 437e5c89e4eSSatish Balay 438e5c89e4eSSatish Balay @*/ 4397087cfbeSBarry Smith PetscErrorCode PetscGetArgs(int *argc,char ***args) 440e5c89e4eSSatish Balay { 441e5c89e4eSSatish Balay PetscFunctionBegin; 44217186662SBarry Smith if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()"); 443e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 444e5c89e4eSSatish Balay *args = PetscGlobalArgs; 445e5c89e4eSSatish Balay PetscFunctionReturn(0); 446e5c89e4eSSatish Balay } 447e5c89e4eSSatish Balay 448793721a6SBarry Smith /*@C 449793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 450793721a6SBarry Smith after PetscInitialize() is called but before PetscFinalize(). 451793721a6SBarry Smith 452793721a6SBarry Smith Not Collective 453793721a6SBarry Smith 454793721a6SBarry Smith Output Parameters: 455793721a6SBarry Smith . args - the command line arguments 456793721a6SBarry Smith 457793721a6SBarry Smith Level: intermediate 458793721a6SBarry Smith 459793721a6SBarry Smith Notes: 460793721a6SBarry Smith This does NOT start with the program name and IS null terminated (final arg is void) 461793721a6SBarry Smith 462793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments() 463793721a6SBarry Smith 464793721a6SBarry Smith @*/ 4657087cfbeSBarry Smith PetscErrorCode PetscGetArguments(char ***args) 466793721a6SBarry Smith { 467793721a6SBarry Smith PetscInt i,argc = PetscGlobalArgc; 468793721a6SBarry Smith PetscErrorCode ierr; 469793721a6SBarry Smith 470793721a6SBarry Smith PetscFunctionBegin; 47117186662SBarry Smith if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()"); 4722d747510SLisandro Dalcin if (!argc) {*args = NULL; PetscFunctionReturn(0);} 473785e854fSJed Brown ierr = PetscMalloc1(argc,args);CHKERRQ(ierr); 474793721a6SBarry Smith for (i=0; i<argc-1; i++) { 475793721a6SBarry Smith ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr); 476793721a6SBarry Smith } 4772d747510SLisandro Dalcin (*args)[argc-1] = NULL; 478793721a6SBarry Smith PetscFunctionReturn(0); 479793721a6SBarry Smith } 480793721a6SBarry Smith 481793721a6SBarry Smith /*@C 482793721a6SBarry Smith PetscFreeArguments - Frees the memory obtained with PetscGetArguments() 483793721a6SBarry Smith 484793721a6SBarry Smith Not Collective 485793721a6SBarry Smith 486793721a6SBarry Smith Output Parameters: 487793721a6SBarry Smith . args - the command line arguments 488793721a6SBarry Smith 489793721a6SBarry Smith Level: intermediate 490793721a6SBarry Smith 491793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments() 492793721a6SBarry Smith 493793721a6SBarry Smith @*/ 4947087cfbeSBarry Smith PetscErrorCode PetscFreeArguments(char **args) 495793721a6SBarry Smith { 496793721a6SBarry Smith PetscInt i = 0; 497793721a6SBarry Smith PetscErrorCode ierr; 498793721a6SBarry Smith 499793721a6SBarry Smith PetscFunctionBegin; 500a297a907SKarl Rupp if (!args) PetscFunctionReturn(0); 501793721a6SBarry Smith while (args[i]) { 502793721a6SBarry Smith ierr = PetscFree(args[i]);CHKERRQ(ierr); 503793721a6SBarry Smith i++; 504793721a6SBarry Smith } 505793721a6SBarry Smith ierr = PetscFree(args);CHKERRQ(ierr); 506793721a6SBarry Smith PetscFunctionReturn(0); 507793721a6SBarry Smith } 508793721a6SBarry Smith 50911525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS) 51030befbd2SBarry Smith #include <petscconfiginfo.h> 51130befbd2SBarry Smith 51295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) 51311525c0dSBarry Smith { 51411525c0dSBarry Smith if (!PetscGlobalRank) { 51530befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64]; 51611525c0dSBarry Smith int port; 517ffbd1cfbSBarry Smith PetscBool flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE; 51811525c0dSBarry Smith size_t applinelen,introlen; 51911525c0dSBarry Smith PetscErrorCode ierr; 520ffbd1cfbSBarry Smith char sawsurl[256]; 52111525c0dSBarry Smith 522c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr); 52311525c0dSBarry Smith if (flg) { 52411525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 52511525c0dSBarry Smith 526c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 52711525c0dSBarry Smith if (sawslog[0]) { 52811525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog)); 52911525c0dSBarry Smith } else { 53011525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL)); 53111525c0dSBarry Smith } 53211525c0dSBarry Smith } 533c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 53411525c0dSBarry Smith if (flg) { 53511525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert)); 53611525c0dSBarry Smith } 537c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr); 538ffbd1cfbSBarry Smith if (selectport) { 539ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Get_Available_Port,(&port)); 540ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 541ffbd1cfbSBarry Smith } else { 542c5929fdfSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr); 54311525c0dSBarry Smith if (flg) { 54411525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 54511525c0dSBarry Smith } 546ffbd1cfbSBarry Smith } 547c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 54811525c0dSBarry Smith if (flg) { 54911525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr); 55011525c0dSBarry Smith ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr); 5519c1e0ce8SBarry Smith } else { 552c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr); 5539c1e0ce8SBarry Smith if (flg) { 5543c01dfcfSBarry Smith ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 5559c1e0ce8SBarry Smith PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr); 5569c1e0ce8SBarry Smith } 55711525c0dSBarry Smith } 558c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr); 55911525c0dSBarry Smith if (flg2) { 56011525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 56111525c0dSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option"); 56211525c0dSBarry Smith ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr); 56311525c0dSBarry Smith ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr); 56411525c0dSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory"); 56543da4ab2SBarry Smith PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr); 56611525c0dSBarry Smith } 56711525c0dSBarry Smith ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr); 56811525c0dSBarry Smith ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr); 56911525c0dSBarry Smith introlen = 4096 + applinelen; 57030a8c9c0SSurtai Han applinelen += 1024; 57111525c0dSBarry Smith ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr); 57211525c0dSBarry Smith ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr); 57311525c0dSBarry Smith 57411525c0dSBarry Smith if (rootlocal) { 57511525c0dSBarry Smith ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr); 57611525c0dSBarry Smith ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr); 57711525c0dSBarry Smith } 57876a34f28SBarry Smith ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr); 57911525c0dSBarry Smith if (rootlocal && help) { 580928bb9adSStefano Zampini ierr = 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);CHKERRQ(ierr); 58111525c0dSBarry Smith } else if (help) { 582928bb9adSStefano Zampini ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr); 58311525c0dSBarry Smith } else { 584928bb9adSStefano Zampini ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr); 58511525c0dSBarry Smith } 586b0bb5815SBarry Smith ierr = PetscFree(options);CHKERRQ(ierr); 58730befbd2SBarry Smith ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr); 58811525c0dSBarry Smith ierr = PetscSNPrintf(intro,introlen,"<body>\n" 589a8d69d7bSBarry Smith "<center><h2> <a href=\"https://www.mcs.anl.gov/petsc\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n" 590df62c222SBarry 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" 591928bb9adSStefano Zampini "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr); 59243da4ab2SBarry Smith PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro)); 59311525c0dSBarry Smith ierr = PetscFree(intro);CHKERRQ(ierr); 59411525c0dSBarry Smith ierr = PetscFree(appline);CHKERRQ(ierr); 595ffbd1cfbSBarry Smith if (selectport) { 596aa573868SBarry Smith PetscBool silent; 5977d812c46SBarry Smith 5987d812c46SBarry Smith ierr = SAWs_Initialize(); 5997d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 6007d812c46SBarry Smith while (ierr) { 6017d812c46SBarry Smith PetscStackCallSAWs(SAWs_Get_Available_Port,(&port)); 6027d812c46SBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 6037d812c46SBarry Smith ierr = SAWs_Initialize(); 6047d812c46SBarry Smith } 6057d812c46SBarry Smith 606aa573868SBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr); 607aa573868SBarry Smith if (!silent) { 608ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl)); 609ffbd1cfbSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr); 610ffbd1cfbSBarry Smith } 6117d812c46SBarry Smith } else { 6127d812c46SBarry Smith PetscStackCallSAWs(SAWs_Initialize,()); 613aa573868SBarry Smith } 6140af79b04SBarry Smith ierr = PetscCitationsRegister("@TechReport{ saws,\n" 6150af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 6160af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 6170af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 6180af79b04SBarry Smith " Year = 2013\n}\n",NULL);CHKERRQ(ierr); 61911525c0dSBarry Smith } 62011525c0dSBarry Smith PetscFunctionReturn(0); 62111525c0dSBarry Smith } 62211525c0dSBarry Smith #endif 62311525c0dSBarry Smith 6244dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 6254dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) 6264dfee713SSatish Balay { 6274dfee713SSatish Balay PetscFunctionBegin; 6284dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 6294dfee713SSatish Balay /* see MPI.py for details on this bug */ 6304dfee713SSatish Balay (void) setenv("HWLOC_COMPONENTS","-x86",1); 6314dfee713SSatish Balay #endif 6324dfee713SSatish Balay PetscFunctionReturn(0); 6334dfee713SSatish Balay } 6344dfee713SSatish Balay 635a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 636a56f64adSBarry Smith #include <adios.h> 63722580e64SBarry Smith #include <adios_read.h> 6387b56e58cSSatish Balay int64_t Petsc_adios_group; 639a56f64adSBarry Smith #endif 64055d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2) 64155d657eeSBarry Smith #include <adios2_c.h> 64255d657eeSBarry Smith #endif 643cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP) 644cd1b99a6SBarry Smith #include <omp.h> 645f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 646cd1b99a6SBarry Smith #endif 647a56f64adSBarry Smith 648bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLFCN_H) 649bc8a24ecSBarry Smith #include <dlfcn.h> 650bc8a24ecSBarry Smith #endif 651bc8a24ecSBarry Smith 652e5c89e4eSSatish Balay /*@C 653e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 654e5c89e4eSSatish Balay PetscInitialize() calls MPI_Init() if that has yet to be called, 655e5c89e4eSSatish Balay so this routine should always be called near the beginning of 656e5c89e4eSSatish Balay your program -- usually the very first line! 657e5c89e4eSSatish Balay 658e5c89e4eSSatish Balay Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set 659e5c89e4eSSatish Balay 660e5c89e4eSSatish Balay Input Parameters: 661e5c89e4eSSatish Balay + argc - count of number of command line arguments 662e5c89e4eSSatish Balay . args - the command line arguments 6630298fd71SBarry Smith . file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for 664fc2bca9aSBarry Smith code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files 6650298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message 666e5c89e4eSSatish Balay 66705827820SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that 66805827820SBarry Smith communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a 66905827820SBarry Smith four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not, 67005827820SBarry Smith then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even 67105827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 672e5c89e4eSSatish Balay 673e5c89e4eSSatish Balay Options Database Keys: 6747ca660e7SBarry Smith + -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message 6757ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 676e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 6778a690491SBarry Smith . -on_error_emacs <machinename> - causes emacsclient to jump to error file 6788a690491SBarry Smith . -on_error_abort - calls abort() when error detected (no traceback) 6798a690491SBarry Smith . -on_error_mpiabort - calls MPI_abort() when error detected 6808a690491SBarry Smith . -error_output_stderr - prints error messages to stderr instead of the default stdout 6818a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 682e5c89e4eSSatish Balay . -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger 683e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds) - Pauses debugger 684e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 685e5c89e4eSSatish Balay process and wait (-debugger_pause) seconds for attachment 68679dccf82SBarry Smith . -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug) 68779dccf82SBarry Smith . -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no) 68879dccf82SBarry Smith . -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug() 689aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 69092f119d6SBarry 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 69192f119d6SBarry Smith . -malloc_view - show a list of all allocated memory during PetscFinalize() 69292f119d6SBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view 69392f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 694e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 695e5c89e4eSSatish Balay . -shared_tmp - indicates /tmp directory is shared by all processors 696e5c89e4eSSatish Balay . -not_shared_tmp - each processor has own /tmp 697e5c89e4eSSatish Balay . -tmp - alternative name of /tmp directory 698e5c89e4eSSatish Balay . -get_total_flops - returns total flops done by all processors 6990841954dSBarry Smith - -memory_view - Print memory usage at end of run 700e5c89e4eSSatish Balay 701e5c89e4eSSatish Balay Options Database Keys for Profiling: 702a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 703*fe9b927eSVaclav Hapla + -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo(). 704217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 705217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 706495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 707e5c89e4eSSatish Balay hangs without running in the debugger). See PetscLogTraceBegin(). 7089a9a5d4cSBarry Smith . -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView(). 70979dccf82SBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView(). 7109a9a5d4cSBarry Smith . -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the 711495fc317SBarry Smith summary is written to the file. See PetscLogView(). 712f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 713495fc317SBarry Smith . -log_all [filename] - Logs extensive profiling information See PetscLogDump(). 714495fc317SBarry Smith . -log [filename] - Logs basic profiline information See PetscLogDump(). 715c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 71687c3beb6SLisandro Dalcin . -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off 717c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 718495fc317SBarry Smith 719609bdbeeSBarry Smith Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time 720e5c89e4eSSatish Balay 721ffbd1cfbSBarry Smith Options Database Keys for SAWs: 722ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 723ffbd1cfbSBarry 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 724ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 725ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 726ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 727ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 728ffbd1cfbSBarry Smith 729e5c89e4eSSatish Balay Environmental Variables: 730e5c89e4eSSatish Balay + PETSC_TMP - alternative tmp directory 731e5c89e4eSSatish Balay . PETSC_SHARED_TMP - tmp is shared by all processes 732e5c89e4eSSatish Balay . PETSC_NOT_SHARED_TMP - each process has its own private tmp 733e5c89e4eSSatish Balay . PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer 734e5c89e4eSSatish Balay - PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to 735e5c89e4eSSatish Balay 736e5c89e4eSSatish Balay 737e5c89e4eSSatish Balay Level: beginner 738e5c89e4eSSatish Balay 739e5c89e4eSSatish Balay Notes: 740e5c89e4eSSatish Balay If for some reason you must call MPI_Init() separately, call 741e5c89e4eSSatish Balay it before PetscInitialize(). 742e5c89e4eSSatish Balay 743e5c89e4eSSatish Balay Fortran Version: 744e5c89e4eSSatish Balay In Fortran this routine has the format 745e5c89e4eSSatish Balay $ call PetscInitialize(file,ierr) 746e5c89e4eSSatish Balay 747e5c89e4eSSatish Balay + ierr - error return code 7480eb4c9c0SBarry Smith - file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for 749fc2bca9aSBarry Smith code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files 750e5c89e4eSSatish Balay 751e5c89e4eSSatish Balay Important Fortran Note: 7520eb4c9c0SBarry Smith In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a 7530298fd71SBarry Smith null character string; you CANNOT just use NULL as 754a7f22e61SSatish Balay in the C version. See Users-Manual: ch_fortran for details. 755e5c89e4eSSatish Balay 75601cb0274SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after 75701cb0274SBarry Smith calling PetscInitialize(). 758e5c89e4eSSatish Balay 75901cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments() 760e5c89e4eSSatish Balay 761e5c89e4eSSatish Balay @*/ 7627087cfbeSBarry Smith PetscErrorCode PetscInitialize(int *argc,char ***args,const char file[],const char help[]) 763e5c89e4eSSatish Balay { 764e5c89e4eSSatish Balay PetscErrorCode ierr; 7654bb5149bSJed Brown PetscMPIInt flag, size; 76619c5658aSBarry Smith PetscBool flg = PETSC_TRUE; 767e5c89e4eSSatish Balay char hostname[256]; 768e5c89e4eSSatish Balay 769e5c89e4eSSatish Balay PetscFunctionBegin; 770e5c89e4eSSatish Balay if (PetscInitializeCalled) PetscFunctionReturn(0); 7713d96e996SBarry Smith /* 7723d96e996SBarry Smith The checking over compatible runtime libraries is complicated by the MPI ABI initiative 7733d96e996SBarry Smith https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 7743d96e996SBarry Smith MPICH v3.1 (Released Feburary 2014) 7753d96e996SBarry Smith IBM MPI v2.1 (December 2014) 7763d96e996SBarry Smith Intel® MPI Library v5.0 (2014) 7773d96e996SBarry Smith Cray MPT v7.0.0 (June 2014) 7783d96e996SBarry Smith As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 7793d96e996SBarry Smith listed above and since that time are compatible. 7803d96e996SBarry Smith 7813d96e996SBarry Smith Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 7823d96e996SBarry Smith at compile time or runtime. Thus we will need to systematically track the allowed versions 7833d96e996SBarry Smith and how they are represented in the mpi.h and MPI_Get_library_version() output in order 7843d96e996SBarry Smith to perform the checking. 7853d96e996SBarry Smith 7863d96e996SBarry Smith Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 7873d96e996SBarry Smith 7883d96e996SBarry Smith Questions: 7893d96e996SBarry Smith 7903d96e996SBarry Smith Should the checks for ABI incompatibility be only on the major version number below? 7913d96e996SBarry Smith Presumably the output to stderr will be removed before a release. 7923d96e996SBarry Smith */ 7933d96e996SBarry Smith 79419c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 79519c5658aSBarry Smith { 79619c5658aSBarry Smith char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 79719c5658aSBarry Smith PetscMPIInt mpilibraryversionlength; 79819c5658aSBarry Smith ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr; 7993d96e996SBarry Smith /* check for MPICH versions before MPI ABI initiative */ 80019c5658aSBarry Smith #if defined(MPICH_VERSION) 8013d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000 80219c5658aSBarry Smith { 80319c5658aSBarry Smith char *ver,*lf; 80419c5658aSBarry Smith flg = PETSC_FALSE; 80519c5658aSBarry Smith ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr; 80619c5658aSBarry Smith if (ver) { 80719c5658aSBarry Smith ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr; 80819c5658aSBarry Smith if (lf) { 80919c5658aSBarry Smith *lf = 0; 81019c5658aSBarry Smith ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr; 81119c5658aSBarry Smith } 81219c5658aSBarry Smith } 81319c5658aSBarry Smith if (!flg) { 81419c5658aSBarry Smith fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION); 8153d96e996SBarry Smith return PETSC_ERR_MPI_LIB_INCOMP; 81619c5658aSBarry Smith } 81719c5658aSBarry Smith } 8183d96e996SBarry Smith #endif 8193d96e996SBarry Smith /* check for OpenMPI version, it is not part of the MPI ABI initiative (is it part of another initiative that needs to be handled?) */ 82019c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION) 82119c5658aSBarry Smith { 82219c5658aSBarry Smith char *ver,bs[32],*bsf; 82319c5658aSBarry Smith flg = PETSC_FALSE; 82419c5658aSBarry Smith ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr; 82519c5658aSBarry Smith if (ver) { 8262e924ca5SSatish Balay PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION); 82719c5658aSBarry Smith ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr; 82819c5658aSBarry Smith if (bsf) flg = PETSC_TRUE; 82919c5658aSBarry Smith } 83019c5658aSBarry Smith if (!flg) { 83119c5658aSBarry Smith fprintf(stderr,"PETSc Error --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d, aborting\n",mpilibraryversion,OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION); 8323d96e996SBarry Smith return PETSC_ERR_MPI_LIB_INCOMP; 83319c5658aSBarry Smith } 83419c5658aSBarry Smith } 83519c5658aSBarry Smith #endif 83619c5658aSBarry Smith } 83719c5658aSBarry Smith #endif 83819c5658aSBarry Smith 839bc8a24ecSBarry Smith #if defined(PETSC_HAVE_DLSYM) 840bc8a24ecSBarry Smith { 841bc8a24ecSBarry Smith PetscInt cnt = 0; 842bc8a24ecSBarry Smith /* 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 */ 843bc8a24ecSBarry Smith if (dlsym(RTLD_DEFAULT,"ompi_mpi_init")) cnt++; 844bc8a24ecSBarry Smith if (dlsym(RTLD_DEFAULT,"MPL_exit")) cnt++; 845bc8a24ecSBarry Smith if (cnt > 1) { 846bc8a24ecSBarry Smith fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n"); 847bc8a24ecSBarry Smith return PETSC_ERR_MPI_LIB_INCOMP; 848bc8a24ecSBarry Smith } 849bc8a24ecSBarry Smith } 850bc8a24ecSBarry Smith #endif 851bc8a24ecSBarry Smith 852ae9b4142SLisandro Dalcin /* these must be initialized in a routine, not as a constant declaration*/ 853d89683f4Sbcordonn PETSC_STDOUT = stdout; 854ae9b4142SLisandro Dalcin PETSC_STDERR = stderr; 855e5c89e4eSSatish Balay 8560c30907bSSatish Balay /* on Windows - set printf to default to printing 2 digit exponents */ 8570c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT) 8580c30907bSSatish Balay _set_output_format(_TWO_DIGIT_EXPONENT); 8590c30907bSSatish Balay #endif 8600c30907bSSatish Balay 8614416b707SBarry Smith ierr = PetscOptionsCreateDefault();CHKERRQ(ierr); 862e5c89e4eSSatish Balay 863e5c89e4eSSatish Balay /* 864e5c89e4eSSatish Balay We initialize the program name here (before MPI_Init()) because MPICH has a bug in 865e5c89e4eSSatish Balay it that it sets args[0] on all processors to be args[0] on the first processor. 866e5c89e4eSSatish Balay */ 867e5c89e4eSSatish Balay if (argc && *argc) { 868e5c89e4eSSatish Balay ierr = PetscSetProgramName(**args);CHKERRQ(ierr); 869e5c89e4eSSatish Balay } else { 870e5c89e4eSSatish Balay ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr); 871e5c89e4eSSatish Balay } 872e5c89e4eSSatish Balay 873e5c89e4eSSatish Balay ierr = MPI_Initialized(&flag);CHKERRQ(ierr); 874e5c89e4eSSatish Balay if (!flag) { 875e32f2f54SBarry Smith if (PETSC_COMM_WORLD != MPI_COMM_NULL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You cannot set PETSC_COMM_WORLD if you have not initialized MPI first"); 8764dfee713SSatish Balay ierr = PetscPreMPIInit_Private();CHKERRQ(ierr); 8775e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 8785e765c61SJed Brown { 8795e765c61SJed Brown PetscMPIInt provided; 8805e765c61SJed Brown ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr); 8815e765c61SJed Brown } 8825e765c61SJed Brown #else 883e5c89e4eSSatish Balay ierr = MPI_Init(argc,args);CHKERRQ(ierr); 8845e765c61SJed Brown #endif 885e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 886e5c89e4eSSatish Balay } 8874dfee713SSatish Balay 888e5c89e4eSSatish Balay if (argc && args) { 889e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 890e5c89e4eSSatish Balay PetscGlobalArgs = *args; 891e5c89e4eSSatish Balay } 892e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_FALSE; 8935ad9ad5bSBarry Smith ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr); 8945ad9ad5bSBarry Smith ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr); 8955ad9ad5bSBarry Smith ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr); 896ef19f930SBarry Smith ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr); 897e5c89e4eSSatish Balay 898a297a907SKarl Rupp if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 899d54338ecSKarl Rupp ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr); 900e5c89e4eSSatish Balay 9010f9be574SPeter Hill if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) { 90212801b39SBarry Smith ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr); 90312801b39SBarry Smith ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr); 9040f9be574SPeter Hill } 90512801b39SBarry Smith 906e5c89e4eSSatish Balay /* Done after init due to a bug in MPICH-GM? */ 907e5c89e4eSSatish Balay ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr); 908e5c89e4eSSatish Balay 909e5c89e4eSSatish Balay ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr); 910e5c89e4eSSatish Balay ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr); 911e5c89e4eSSatish Balay 9128ad47952SJed Brown MPIU_BOOL = MPI_INT; 9138ad47952SJed Brown MPIU_ENUM = MPI_INT; 9147cdaf61dSJed Brown MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64; 915e316c87fSJed Brown if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 916e316c87fSJed Brown else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 917e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG) 918e316c87fSJed Brown else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 919e316c87fSJed Brown #endif 920e316c87fSJed Brown else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;} 9218ad47952SJed Brown 922e5c89e4eSSatish Balay /* 923e5c89e4eSSatish Balay Initialized the global complex variable; this is because with 924e5c89e4eSSatish Balay shared libraries the constructors for global variables 925e5c89e4eSSatish Balay are not called; at least on IRIX. 926e5c89e4eSSatish Balay */ 927886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX) 928e5c89e4eSSatish Balay { 929252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 93050f81f78SJed Brown PetscComplex ic(0.0,1.0); 931e5c89e4eSSatish Balay PETSC_i = ic; 932252ecd31SSatish Balay #else 93350f81f78SJed Brown PETSC_i = _Complex_I; 934b7940d39SSatish Balay #endif 935762437b8SSatish Balay } 936762437b8SSatish Balay 9372c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX) 938e69cd0e6SSatish Balay ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr); 939500d8756SSatish Balay ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr); 940500d8756SSatish Balay ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr); 941500d8756SSatish Balay ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr); 9422c876bd9SBarry Smith #endif 943886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */ 944e5c89e4eSSatish Balay 945e5c89e4eSSatish Balay /* 946e5c89e4eSSatish Balay Create the PETSc MPI reduction operator that sums of the first 947e5c89e4eSSatish Balay half of the entries and maxes the second half. 948e5c89e4eSSatish Balay */ 949367daffbSBarry Smith ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr); 950e5c89e4eSSatish Balay 951ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 952c90a1750SBarry Smith ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr); 953c90a1750SBarry Smith ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr); 9547c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 9558c764dc5SJose Roman ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr); 9568c764dc5SJose Roman ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr); 9578c764dc5SJose Roman #endif 958d9822059SBarry Smith ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr); 959d9822059SBarry Smith ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr); 960570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 961570b7f6dSBarry Smith ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr); 962570b7f6dSBarry Smith ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr); 963570b7f6dSBarry Smith ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr); 964570b7f6dSBarry Smith ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr); 965c90a1750SBarry Smith #endif 966c90a1750SBarry Smith 967570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 968cca4cb22SSatish Balay ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr); 969cca4cb22SSatish Balay #endif 970cca4cb22SSatish Balay 971e5c89e4eSSatish Balay ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr); 972e5c89e4eSSatish Balay ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr); 973e5c89e4eSSatish Balay 97440df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 975e5c89e4eSSatish Balay ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr); 976e5c89e4eSSatish Balay ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr); 97744041f26SJed Brown #endif 978e5c89e4eSSatish Balay 979ec957eceSBarry Smith 980e5c89e4eSSatish Balay /* 981480cf27aSJed Brown Attributes to be set on PETSc communicators 982480cf27aSJed Brown */ 98333779a13SJunchao Zhang ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr); 98433779a13SJunchao Zhang ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr); 98533779a13SJunchao Zhang ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr); 98633779a13SJunchao Zhang ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr); 987480cf27aSJed Brown 988480cf27aSJed Brown /* 989e8fb0fc0SBarry Smith Build the options database 990e5c89e4eSSatish Balay */ 991c5929fdfSBarry Smith ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr); 992e5c89e4eSSatish Balay 993703f3eceSBarry Smith /* call a second time so it can look in the options database */ 994703f3eceSBarry Smith ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr); 9956dc8fec2Sbcordonn 996e5c89e4eSSatish Balay /* 997e5c89e4eSSatish Balay Print main application help message 998e5c89e4eSSatish Balay */ 9992d747510SLisandro Dalcin ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr); 1000e5c89e4eSSatish Balay if (help && flg) { 1001e5c89e4eSSatish Balay ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr); 1002e5c89e4eSSatish Balay } 1003e5c89e4eSSatish Balay ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr); 1004e5c89e4eSSatish Balay 1005d45a07a7SBarry Smith ierr = PetscCitationsInitialize();CHKERRQ(ierr); 1006d45a07a7SBarry Smith 1007e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 100811525c0dSBarry Smith ierr = PetscInitializeSAWs(help);CHKERRQ(ierr); 1009f4202a44SBarry Smith #endif 1010f4202a44SBarry Smith 1011e5c89e4eSSatish Balay /* 1012e5c89e4eSSatish Balay Load the dynamic libraries (on machines that support them), this registers all 1013e5c89e4eSSatish Balay the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 1014e5c89e4eSSatish Balay */ 1015e5c89e4eSSatish Balay ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr); 1016e5c89e4eSSatish Balay 1017e5c89e4eSSatish Balay ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 101802c9f0b5SLisandro Dalcin ierr = PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr); 1019e5c89e4eSSatish Balay ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr); 102002c9f0b5SLisandro Dalcin ierr = PetscInfo1(NULL,"Running on machine: %s\n",hostname);CHKERRQ(ierr); 1021cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP) 1022101491d0SBarry Smith { 1023101491d0SBarry Smith PetscBool omp_view_flag; 1024cd1b99a6SBarry Smith char *threads = getenv("OMP_NUM_THREADS"); 1025e5c89e4eSSatish Balay 1026cd1b99a6SBarry Smith if (threads) { 102702c9f0b5SLisandro Dalcin ierr = PetscInfo1(NULL,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr); 1028f90b075cSBarry Smith (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads); 1029cd1b99a6SBarry Smith } else { 1030cd1b99a6SBarry Smith #define NMAX 10000 1031101491d0SBarry Smith int i; 1032cd1b99a6SBarry Smith PetscScalar *x; 1033cd1b99a6SBarry Smith ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr); 1034cd1b99a6SBarry Smith #pragma omp parallel for 1035cd1b99a6SBarry Smith for (i=0; i<NMAX; i++) { 1036cd1b99a6SBarry Smith x[i] = 0.0; 1037f90b075cSBarry Smith PetscNumOMPThreads = (PetscInt) omp_get_num_threads(); 1038cd1b99a6SBarry Smith } 1039cd1b99a6SBarry Smith ierr = PetscFree(x);CHKERRQ(ierr); 104002c9f0b5SLisandro Dalcin ierr = PetscInfo1(NULL,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr); 1041101491d0SBarry Smith } 1042101491d0SBarry Smith ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr); 1043f90b075cSBarry Smith ierr = PetscOptionsInt("-omp_num_threads","Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS","None",PetscNumOMPThreads,&PetscNumOMPThreads,&flg);CHKERRQ(ierr); 1044101491d0SBarry Smith ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr); 1045101491d0SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 1046101491d0SBarry Smith if (flg) { 104702c9f0b5SLisandro Dalcin ierr = PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr); 1048f90b075cSBarry Smith omp_set_num_threads((int)PetscNumOMPThreads); 1049101491d0SBarry Smith } 1050101491d0SBarry Smith if (omp_view_flag) { 1051f90b075cSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr); 1052cd1b99a6SBarry Smith } 1053cd1b99a6SBarry Smith } 1054cd1b99a6SBarry Smith #endif 1055ef6c6fedSBoyana Norris /* Check the options database for options related to the options database itself */ 1056c5929fdfSBarry Smith ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr); 1057ef6c6fedSBoyana Norris 1058951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 1059e39fd77fSBarry Smith /* 1060e39fd77fSBarry Smith Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 1061e39fd77fSBarry Smith 1062e39fd77fSBarry Smith Currently not used because it is not supported by MPICH. 1063e39fd77fSBarry Smith */ 106430815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) { 10650298fd71SBarry Smith ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr); 106630815ce0SLisandro Dalcin } 1067951e3c8eSBarry Smith #endif 1068e39fd77fSBarry Smith 106941c0b4b3SShri Abhyankar /* 107041c0b4b3SShri Abhyankar Setup building of stack frames for all function calls 107141c0b4b3SShri Abhyankar */ 1072ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1073e1167bb9SShri Abhyankar ierr = PetscStackCreate();CHKERRQ(ierr); 1074e1167bb9SShri Abhyankar #endif 1075e1167bb9SShri Abhyankar 10762d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 10772d53ad75SBarry Smith ierr = PetscFPTCreate(10000);CHKERRQ(ierr); 10782d53ad75SBarry Smith #endif 10792d53ad75SBarry Smith 10805e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC) 108187c3beb6SLisandro Dalcin { 108287c3beb6SLisandro Dalcin PetscViewer viewer; 108322e7e69cSBarry Smith ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr); 10845e71baefSBarry Smith if (flg) { 10855e71baefSBarry Smith ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr); 108687c3beb6SLisandro Dalcin ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 108787c3beb6SLisandro Dalcin } 10885e71baefSBarry Smith } 10895e71baefSBarry Smith #endif 1090dff31646SBarry Smith 109187c3beb6SLisandro Dalcin flg = PETSC_TRUE; 109287c3beb6SLisandro Dalcin ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr); 109387c3beb6SLisandro Dalcin if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);} 109487c3beb6SLisandro Dalcin 1095a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 1096a56f64adSBarry Smith ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr); 10977b56e58cSSatish Balay ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr); 1098a56f64adSBarry Smith ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr); 10999fc7e16cSBarry Smith ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr); 1100a56f64adSBarry Smith #endif 110155d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2) 110255d657eeSBarry Smith #endif 1103a56f64adSBarry Smith 1104301d30feSBarry Smith /* 1105301d30feSBarry Smith Once we are completedly initialized then we can set this variables 1106301d30feSBarry Smith */ 1107301d30feSBarry Smith PetscInitializeCalled = PETSC_TRUE; 11082db0e300SLisandro Dalcin 11092db0e300SLisandro Dalcin ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr); 11102db0e300SLisandro Dalcin if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);} 1111301d30feSBarry Smith PetscFunctionReturn(0); 1112e5c89e4eSSatish Balay } 1113e5c89e4eSSatish Balay 11144097062eSBarry Smith #if defined(PETSC_USE_LOG) 111595c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1116ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1117ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 111895c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 11194097062eSBarry Smith #endif 1120e5c89e4eSSatish Balay 1121008a6e76SBarry Smith /* 1122008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1123008a6e76SBarry Smith */ 1124008a6e76SBarry Smith PetscErrorCode PetscFreeMPIResources(void) 1125008a6e76SBarry Smith { 1126008a6e76SBarry Smith PetscErrorCode ierr; 1127008a6e76SBarry Smith 1128008a6e76SBarry Smith PetscFunctionBegin; 1129008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 1130008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr); 1131008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1132008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr); 1133008a6e76SBarry Smith #endif 1134008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr); 1135008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr); 1136008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16) 1137008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr); 1138008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr); 1139008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr); 1140008a6e76SBarry Smith #endif 1141008a6e76SBarry Smith 1142008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1143008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX) 1144008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr); 1145008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr); 1146008a6e76SBarry Smith #endif 1147008a6e76SBarry Smith #endif 1148008a6e76SBarry Smith 1149008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 1150008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr); 1151008a6e76SBarry Smith #endif 1152008a6e76SBarry Smith 1153008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr); 115440df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 1155008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr); 1156008a6e76SBarry Smith #endif 1157008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr); 1158008a6e76SBarry Smith PetscFunctionReturn(0); 1159008a6e76SBarry Smith } 1160008a6e76SBarry Smith 1161e5c89e4eSSatish Balay /*@C 1162e5c89e4eSSatish Balay PetscFinalize - Checks for options to be called at the conclusion 1163e5c89e4eSSatish Balay of the program. MPI_Finalize() is called only if the user had not 1164e5c89e4eSSatish Balay called MPI_Init() before calling PetscInitialize(). 1165e5c89e4eSSatish Balay 1166e5c89e4eSSatish Balay Collective on PETSC_COMM_WORLD 1167e5c89e4eSSatish Balay 1168e5c89e4eSSatish Balay Options Database Keys: 116926a7e8d4SBarry Smith + -options_view - Calls PetscOptionsView() 1170e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 11717eb1d149SBarry 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 1172e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 117392f119d6SBarry Smith . -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed 1174e5c89e4eSSatish Balay . -malloc_info - Prints total memory usage 117592f119d6SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and where 1176e5c89e4eSSatish Balay 1177e5c89e4eSSatish Balay Level: beginner 1178e5c89e4eSSatish Balay 1179e5c89e4eSSatish Balay Note: 1180e5c89e4eSSatish Balay See PetscInitialize() for more general runtime options. 1181e5c89e4eSSatish Balay 118288c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd() 1183e5c89e4eSSatish Balay @*/ 11847087cfbeSBarry Smith PetscErrorCode PetscFinalize(void) 1185e5c89e4eSSatish Balay { 1186e5c89e4eSSatish Balay PetscErrorCode ierr; 11874bb5149bSJed Brown PetscMPIInt rank; 1188a8d2bbe5SBarry Smith PetscInt nopt; 11892bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE; 1190dff31646SBarry Smith PetscBool flg; 119110463e74SBarry Smith #if defined(PETSC_USE_LOG) 119210463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 119310463e74SBarry Smith #endif 1194e5c89e4eSSatish Balay 1195e5c89e4eSSatish Balay if (!PetscInitializeCalled) { 11964b09e917SBarry Smith printf("PetscInitialize() must be called before PetscFinalize()\n"); 11973cbbc5ffSBarry Smith return(PETSC_ERR_ARG_WRONGSTATE); 1198e5c89e4eSSatish Balay } 11993cbbc5ffSBarry Smith PetscFunctionBegin; 12000298fd71SBarry Smith ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr); 1201b022a5c1SBarry Smith 12021f817a21SBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 1203a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 120422580e64SBarry Smith ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr); 1205a56f64adSBarry Smith ierr = adios_finalize(rank);CHKERRQ(ierr); 1206a56f64adSBarry Smith #endif 120755d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2) 120855d657eeSBarry Smith #endif 1209c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr); 1210dff31646SBarry Smith if (flg) { 12111f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 12121f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 12131f817a21SBarry Smith 1214c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 12151f817a21SBarry Smith if (filename[0]) { 12161f817a21SBarry Smith ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr); 12171f817a21SBarry Smith } 1218dff31646SBarry Smith ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr); 1219dff31646SBarry Smith cits[0] = 0; 1220dff31646SBarry Smith ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr); 12211f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr); 12221f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr); 12231f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr); 12241f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr); 12251f817a21SBarry Smith ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr); 1226dff31646SBarry Smith ierr = PetscFree(cits);CHKERRQ(ierr); 1227dff31646SBarry Smith } 1228dff31646SBarry Smith ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr); 1229dff31646SBarry Smith 1230c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER) 123104102261SBarry Smith /* TextBelt is run for testing purposes only, please do not use this feature often */ 123204102261SBarry Smith { 123304102261SBarry Smith PetscInt nmax = 2; 123404102261SBarry Smith char **buffs; 123504102261SBarry Smith ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr); 1236c5929fdfSBarry Smith ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr); 123704102261SBarry Smith if (flg1) { 123804102261SBarry Smith if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\""); 123904102261SBarry Smith if (nmax == 1) { 124004102261SBarry Smith ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr); 124104102261SBarry Smith ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr); 124204102261SBarry Smith ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr); 124304102261SBarry Smith } 124404102261SBarry Smith ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr); 124504102261SBarry Smith ierr = PetscFree(buffs[0]);CHKERRQ(ierr); 124604102261SBarry Smith ierr = PetscFree(buffs[1]);CHKERRQ(ierr); 124704102261SBarry Smith } 124804102261SBarry Smith ierr = PetscFree(buffs);CHKERRQ(ierr); 124904102261SBarry Smith } 1250763ec7b1SBarry Smith { 1251763ec7b1SBarry Smith PetscInt nmax = 2; 1252763ec7b1SBarry Smith char **buffs; 1253763ec7b1SBarry Smith ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr); 1254763ec7b1SBarry Smith ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr); 1255763ec7b1SBarry Smith if (flg1) { 1256763ec7b1SBarry Smith if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\""); 1257763ec7b1SBarry Smith if (nmax == 1) { 1258763ec7b1SBarry Smith ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr); 1259763ec7b1SBarry Smith ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr); 1260763ec7b1SBarry Smith ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr); 1261763ec7b1SBarry Smith } 1262763ec7b1SBarry Smith ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr); 1263763ec7b1SBarry Smith ierr = PetscFree(buffs[0]);CHKERRQ(ierr); 1264763ec7b1SBarry Smith ierr = PetscFree(buffs[1]);CHKERRQ(ierr); 1265763ec7b1SBarry Smith } 1266763ec7b1SBarry Smith ierr = PetscFree(buffs);CHKERRQ(ierr); 1267763ec7b1SBarry Smith } 126804102261SBarry Smith #endif 126967234432SDmitry Karpeev /* 127067234432SDmitry Karpeev It should be safe to cancel the options monitors, since we don't expect to be setting options 127167234432SDmitry Karpeev here (at least that are worth monitoring). Monitors ought to be released so that they release 127267234432SDmitry Karpeev whatever memory was allocated there before -malloc_dump reports unfreed memory. 127367234432SDmitry Karpeev */ 127467234432SDmitry Karpeev ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr); 127504102261SBarry Smith 12762d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 12772d53ad75SBarry Smith ierr = PetscFPTDestroy();CHKERRQ(ierr); 12782d53ad75SBarry Smith #endif 12792d53ad75SBarry Smith 12802d53ad75SBarry Smith 1281e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1282dff31646SBarry Smith flg = PETSC_FALSE; 1283c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr); 1284d5649816SBarry Smith if (flg) { 1285e04113cfSBarry Smith ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr); 1286d5649816SBarry Smith } 1287d5649816SBarry Smith #endif 1288d5649816SBarry Smith 1289681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1290681455b2SBarry Smith flg1 = PETSC_FALSE; 1291c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr); 1292681455b2SBarry Smith if (flg1) { 1293681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 1294681455b2SBarry Smith ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr); 1295681455b2SBarry Smith } 1296681455b2SBarry Smith #endif 1297681455b2SBarry Smith 129867584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 1299c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr); 1300e5c89e4eSSatish Balay if (!flg2) { 130190d69ab7SBarry Smith flg2 = PETSC_FALSE; 1302c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr); 1303e5c89e4eSSatish Balay } 1304e5c89e4eSSatish Balay if (flg2) { 13050841954dSBarry Smith ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr); 1306e5c89e4eSSatish Balay } 130767584ceeSBarry Smith #endif 1308e5c89e4eSSatish Balay 1309e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 131090d69ab7SBarry Smith flg1 = PETSC_FALSE; 1311c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr); 1312e5c89e4eSSatish Balay if (flg1) { 1313e5c89e4eSSatish Balay PetscLogDouble flops = 0; 1314205a32c2SJed Brown ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr); 1315e5c89e4eSSatish Balay ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr); 1316e5c89e4eSSatish Balay } 1317e5c89e4eSSatish Balay #endif 1318e5c89e4eSSatish Balay 1319e5c89e4eSSatish Balay 1320e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 1321e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE) 1322e5c89e4eSSatish Balay mname[0] = 0; 1323c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr); 1324e5c89e4eSSatish Balay if (flg1) { 1325e5c89e4eSSatish Balay if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);} 1326e5c89e4eSSatish Balay else {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);} 1327e5c89e4eSSatish Balay } 1328e5c89e4eSSatish Balay #endif 1329356e5837SBarry Smith #endif 1330a297a907SKarl Rupp 1331dd710f27SBarry Smith /* 1332dd710f27SBarry Smith Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 1333dd710f27SBarry Smith */ 1334dd710f27SBarry Smith ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr); 1335dd710f27SBarry Smith 1336356e5837SBarry Smith #if defined(PETSC_USE_LOG) 133787c3beb6SLisandro Dalcin ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr); 1338f14045dbSBarry Smith ierr = PetscLogViewFromOptions();CHKERRQ(ierr); 133987c3beb6SLisandro Dalcin ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr); 134087c3beb6SLisandro Dalcin 1341356e5837SBarry Smith mname[0] = 0; 1342c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr); 1343e5c89e4eSSatish Balay if (flg1) { 134491eabc43SBarry Smith PetscViewer viewer; 134520a8bfc3SBarry Smith ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING: -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr); 134691eabc43SBarry Smith if (mname[0]) { 134791eabc43SBarry Smith ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr); 134891eabc43SBarry Smith ierr = PetscLogView(viewer);CHKERRQ(ierr); 13496bf464f9SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 135033f85c2fSBarry Smith } else { 135133f85c2fSBarry Smith viewer = PETSC_VIEWER_STDOUT_WORLD; 13529a9a5d4cSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 135333f85c2fSBarry Smith ierr = PetscLogView(viewer);CHKERRQ(ierr); 13549a9a5d4cSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 135533f85c2fSBarry Smith } 1356e5c89e4eSSatish Balay } 1357a297a907SKarl Rupp 1358dd710f27SBarry Smith /* 1359dd710f27SBarry Smith Free any objects created by the last block of code. 1360dd710f27SBarry Smith */ 1361dd710f27SBarry Smith ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr); 1362dd710f27SBarry Smith 1363dd710f27SBarry Smith mname[0] = 0; 1364c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr); 1365c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr); 13667ff663adSLisandro Dalcin if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);} 1367e5c89e4eSSatish Balay #endif 136810463e74SBarry Smith 136933f85c2fSBarry Smith ierr = PetscStackDestroy();CHKERRQ(ierr); 137010463e74SBarry Smith 137190d69ab7SBarry Smith flg1 = PETSC_FALSE; 1372c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr); 1373e5c89e4eSSatish Balay if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);} 137490d69ab7SBarry Smith flg1 = PETSC_FALSE; 1375c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr); 1376e5c89e4eSSatish Balay if (flg1) { 1377e5c89e4eSSatish Balay ierr = PetscMPIDump(stdout);CHKERRQ(ierr); 1378e5c89e4eSSatish Balay } 137990d69ab7SBarry Smith flg1 = PETSC_FALSE; 138090d69ab7SBarry Smith flg2 = PETSC_FALSE; 13818bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 1382c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr); 1383c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr); 1384c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr); 1385e4c476e2SSatish Balay 1386e5c89e4eSSatish Balay if (flg2) { 1387be56827dSJed Brown PetscViewer viewer; 138802ba9f54SBarry Smith ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr); 138902ba9f54SBarry Smith ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr); 1390c5929fdfSBarry Smith ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr); 1391be56827dSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1392e5c89e4eSSatish Balay } 1393e5c89e4eSSatish Balay 1394e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 1395c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr); 1396c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr); 1397e5c89e4eSSatish Balay 139833fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 1399c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr); 14003de2bfdfSBarry Smith #if defined(PETSC_USE_DEBUG) 14013de2bfdfSBarry Smith if (!flg1) flg3 = PETSC_TRUE; 14023de2bfdfSBarry Smith #endif 1403e5c89e4eSSatish Balay if (flg3) { 14043de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 1405be56827dSJed Brown PetscViewer viewer; 140602ba9f54SBarry Smith ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr); 140702ba9f54SBarry Smith ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr); 1408c5929fdfSBarry Smith ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr); 1409be56827dSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1410e5c89e4eSSatish Balay } 14113de2bfdfSBarry Smith ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr); 14123de2bfdfSBarry Smith if (nopt) { 14133de2bfdfSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr); 14143de2bfdfSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr); 14153de2bfdfSBarry Smith if (nopt == 1) { 1416e5c89e4eSSatish Balay ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr); 1417e5c89e4eSSatish Balay } else { 14187582186dSLisandro Dalcin ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr); 1419e5c89e4eSSatish Balay } 14203de2bfdfSBarry Smith } else if (flg3 && flg1) { 14213de2bfdfSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr); 1422df12ba86SBarry Smith } 1423c5929fdfSBarry Smith ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr); 1424e5c89e4eSSatish Balay } 1425e5c89e4eSSatish Balay 1426e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1427d45a07a7SBarry Smith if (!PetscGlobalRank) { 142887f587eeSBarry Smith ierr = PetscStackSAWsViewOff();CHKERRQ(ierr); 142916ad0300SBarry Smith PetscStackCallSAWs(SAWs_Finalize,()); 1430d45a07a7SBarry Smith } 1431ec957eceSBarry Smith #endif 1432ec957eceSBarry Smith 14334097062eSBarry Smith #if defined(PETSC_USE_LOG) 143410463e74SBarry Smith /* 1435dbc8283eSBarry Smith List all objects the user may have forgot to free 14362eff7a51SBarry Smith */ 143705df10baSBarry Smith if (PetscObjectsLog) { 1438c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr); 1439a64a8e02SBarry Smith if (flg1) { 1440a64a8e02SBarry Smith MPI_Comm local_comm; 14417eb1d149SBarry Smith char string[64]; 1442a64a8e02SBarry Smith 1443c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr); 1444a64a8e02SBarry Smith ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr); 1445a64a8e02SBarry Smith ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr); 14467eb1d149SBarry Smith ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 1447a64a8e02SBarry Smith ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr); 1448a64a8e02SBarry Smith ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr); 14490a1571b3SBarry Smith } 145005df10baSBarry Smith } 14514097062eSBarry Smith #endif 14524097062eSBarry Smith 14534097062eSBarry Smith #if defined(PETSC_USE_LOG) 1454dbc8283eSBarry Smith PetscObjectsCounts = 0; 1455dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 1456a297a907SKarl Rupp ierr = PetscFree(PetscObjects);CHKERRQ(ierr); 14574097062eSBarry Smith #endif 14582eff7a51SBarry Smith 145933f85c2fSBarry Smith /* 146033f85c2fSBarry Smith Destroy any packages that registered a finalize 146133f85c2fSBarry Smith */ 146233f85c2fSBarry Smith ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr); 146333f85c2fSBarry Smith 1464101409b8SToby Isaac #if defined(PETSC_USE_LOG) 1465fa2bb9feSLisandro Dalcin ierr = PetscLogFinalize();CHKERRQ(ierr); 1466101409b8SToby Isaac #endif 1467101409b8SToby Isaac 146833f85c2fSBarry Smith /* 146948dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 147048dd1dffSBarry Smith 147137e93019SBarry Smith ierr = PetscFunctionListPrintAll();CHKERRQ(ierr); 147248dd1dffSBarry Smith */ 147337e93019SBarry Smith 14744028d114SSatish Balay if (petsc_history) { 1475f3dea69dSBarry Smith ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr); 147602c9f0b5SLisandro Dalcin petsc_history = NULL; 1477e5c89e4eSSatish Balay } 14789de0f6ecSBarry Smith ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr); 1479e94e781bSJacob Faibussowitsch ierr = PetscInfoDestroy();CHKERRQ(ierr); 1480e5c89e4eSSatish Balay 148167584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 148292f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1483e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 148492f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1485e5c89e4eSSatish Balay FILE *fd; 1486ed9cf6e9SBarry Smith int err; 1487e5c89e4eSSatish Balay 1488dc92acbaSJed Brown flg2 = PETSC_FALSE; 148992f119d6SBarry Smith flg3 = PETSC_FALSE; 14908bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG) 149192f119d6SBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr); 1492dc92acbaSJed Brown #endif 149392f119d6SBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr); 149492f119d6SBarry Smith fname[0] = 0; 149592f119d6SBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr); 1496e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1497e5c89e4eSSatish Balay 14982e924ca5SSatish Balay PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank); 1499e32f2f54SBarry Smith fd = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname); 1500e5c89e4eSSatish Balay ierr = PetscMallocDump(fd);CHKERRQ(ierr); 1501ed9cf6e9SBarry Smith err = fclose(fd); 1502e32f2f54SBarry Smith if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 150392f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1504e5c89e4eSSatish Balay MPI_Comm local_comm; 1505e5c89e4eSSatish Balay 1506e5c89e4eSSatish Balay ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr); 1507e5c89e4eSSatish Balay ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr); 1508e5c89e4eSSatish Balay ierr = PetscMallocDump(stdout);CHKERRQ(ierr); 1509e5c89e4eSSatish Balay ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr); 1510e5c89e4eSSatish Balay ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr); 1511e5c89e4eSSatish Balay } 1512e5c89e4eSSatish Balay fname[0] = 0; 151392f119d6SBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,250,&flg1);CHKERRQ(ierr); 1514e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1515e5c89e4eSSatish Balay 151692f119d6SBarry Smith PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank); 151792f119d6SBarry Smith fd = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname); 151892f119d6SBarry Smith ierr = PetscMallocView(fd);CHKERRQ(ierr); 1519ed9cf6e9SBarry Smith err = fclose(fd); 1520e32f2f54SBarry Smith if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 152192f119d6SBarry Smith } else if (flg1) { 152292f119d6SBarry Smith MPI_Comm local_comm; 152392f119d6SBarry Smith 152492f119d6SBarry Smith ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr); 152592f119d6SBarry Smith ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr); 152692f119d6SBarry Smith ierr = PetscMallocView(stdout);CHKERRQ(ierr); 152792f119d6SBarry Smith ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr); 152892f119d6SBarry Smith ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr); 1529e5c89e4eSSatish Balay } 1530e5c89e4eSSatish Balay } 153167584ceeSBarry Smith #endif 153220e2c332SMatthew G. Knepley 15335486ca60SMatthew G. Knepley /* 15345486ca60SMatthew G. Knepley Close any open dynamic libraries 15355486ca60SMatthew G. Knepley */ 15365486ca60SMatthew G. Knepley ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr); 15375486ca60SMatthew G. Knepley 1538e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 15394416b707SBarry Smith ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr); 1540e5c89e4eSSatish Balay 1541e5c89e4eSSatish Balay PetscGlobalArgc = 0; 154202c9f0b5SLisandro Dalcin PetscGlobalArgs = NULL; 1543e5c89e4eSSatish Balay 1544008a6e76SBarry Smith ierr = PetscFreeMPIResources();CHKERRQ(ierr); 1545e5c89e4eSSatish Balay 1546dbc8283eSBarry Smith /* 1547efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1548efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1549efb80d3cSBarry Smith 1550efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1551efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1552dbc8283eSBarry Smith */ 1553b770b1f6SSatish Balay { 1554dbc8283eSBarry Smith PetscCommCounter *counter; 1555dbc8283eSBarry Smith PetscMPIInt flg; 1556dbc8283eSBarry Smith MPI_Comm icomm; 1557265f3f35SJed Brown union {MPI_Comm comm; void *ptr;} ucomm; 155847435625SJed Brown ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr); 1559dbc8283eSBarry Smith if (flg) { 1560265f3f35SJed Brown icomm = ucomm.comm; 156147435625SJed Brown ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr); 1562dbc8283eSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1563dbc8283eSBarry Smith 156447435625SJed Brown ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr); 156547435625SJed Brown ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr); 1566efb80d3cSBarry Smith ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr); 1567dbc8283eSBarry Smith } 156847435625SJed Brown ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr); 1569dbc8283eSBarry Smith if (flg) { 1570265f3f35SJed Brown icomm = ucomm.comm; 157147435625SJed Brown ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr); 1572dbc8283eSBarry Smith if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory"); 1573dbc8283eSBarry Smith 157447435625SJed Brown ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr); 157547435625SJed Brown ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr); 1576efb80d3cSBarry Smith ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr); 1577dbc8283eSBarry Smith } 1578b770b1f6SSatish Balay } 1579dbc8283eSBarry Smith 158047435625SJed Brown ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr); 158147435625SJed Brown ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr); 158247435625SJed Brown ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr); 15835f7487a0SJunchao Zhang ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr); 1584480cf27aSJed Brown 15855ad9ad5bSBarry Smith ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr); 15865ad9ad5bSBarry Smith ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr); 15875ad9ad5bSBarry Smith ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr); 1588ef19f930SBarry Smith ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr); 1589ef19f930SBarry Smith 1590e5c89e4eSSatish Balay if (PetscBeganMPI) { 159199608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED) 159299b1327fSBarry Smith PetscMPIInt flag; 159399b1327fSBarry Smith ierr = MPI_Finalized(&flag);CHKERRQ(ierr); 1594e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 159599608316SBarry Smith #endif 1596e5c89e4eSSatish Balay ierr = MPI_Finalize();CHKERRQ(ierr); 1597e5c89e4eSSatish Balay } 1598e5c89e4eSSatish Balay /* 1599e5c89e4eSSatish Balay 1600e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1601e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1602e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1603e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1604e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 16050e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1606e5c89e4eSSatish Balay memory was not freed. 1607e5c89e4eSSatish Balay 1608e5c89e4eSSatish Balay */ 16091d1a0024SBarry Smith ierr = PetscMallocClear();CHKERRQ(ierr); 1610a297a907SKarl Rupp 1611e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1612e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 16133db9a53dSBarry Smith PetscFunctionReturn(0); 1614e5c89e4eSSatish Balay } 1615e5c89e4eSSatish Balay 161643db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 16178cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b) 161843db4dbbSBarry Smith { 161943db4dbbSBarry Smith if (*a == *b) return 1; 162043db4dbbSBarry Smith if (*a + 32 == *b) return 1; 162143db4dbbSBarry Smith if (*a - 32 == *b) return 1; 162243db4dbbSBarry Smith return 0; 162343db4dbbSBarry Smith } 1624a70650f6SBarry Smith #endif 162543db4dbbSBarry Smith 162643db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 16278cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b) 162843db4dbbSBarry Smith { 162943db4dbbSBarry Smith if (*a == *b) return 1; 163043db4dbbSBarry Smith if (*a + 32 == *b) return 1; 163143db4dbbSBarry Smith if (*a - 32 == *b) return 1; 163243db4dbbSBarry Smith return 0; 163343db4dbbSBarry Smith } 163443db4dbbSBarry Smith #endif 1635