17d0a6c19SBarry Smith 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 21e5c89e4eSSatish Balay /* -----------------------------------------------------------------------------------------*/ 22e5c89e4eSSatish Balay 2395c0884eSLisandro Dalcin PETSC_INTERN FILE *petsc_history; 24e5c89e4eSSatish Balay 2595c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void); 2695c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void); 2795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void); 2895c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int); 2995c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int); 3095c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**); 310069ddf5SShri Abhyankar 32e5c89e4eSSatish Balay /* user may set this BEFORE calling PetscInitialize() */ 33e8373e55SMatthew Knepley MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL; 34e5c89e4eSSatish Balay 35480cf27aSJed Brown PetscMPIInt Petsc_Counter_keyval = MPI_KEYVAL_INVALID; 36480cf27aSJed Brown PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID; 37480cf27aSJed Brown PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID; 385f7487a0SJunchao Zhang PetscMPIInt Petsc_ShmComm_keyval = MPI_KEYVAL_INVALID; 39480cf27aSJed Brown 40e5c89e4eSSatish Balay /* 41e5c89e4eSSatish Balay Declare and set all the string names of the PETSc enums 42e5c89e4eSSatish Balay */ 436a6fc655SJed Brown const char *const PetscBools[] = {"FALSE","TRUE","PetscBool","PETSC_",0}; 446a6fc655SJed Brown const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",0}; 45e5c89e4eSSatish Balay 46ace3abfcSBarry Smith PetscBool PetscPreLoadingUsed = PETSC_FALSE; 47ace3abfcSBarry Smith PetscBool PetscPreLoadingOn = PETSC_FALSE; 480f8e0872SSatish Balay 49a2f94806SJed Brown PetscInt PetscHotRegionDepth; 50a2f94806SJed Brown 51b22622e2STadeu Manoel #if defined(PETSC_HAVE_THREADSAFETY) 52b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockOpen; 53b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStdout; 54b22622e2STadeu Manoel PetscSpinlock PetscViewerASCIISpinLockStderr; 55b22622e2STadeu Manoel PetscSpinlock PetscCommSpinLock; 56b22622e2STadeu Manoel #endif 57b22622e2STadeu Manoel 58e5c89e4eSSatish Balay /* 59e5c89e4eSSatish Balay Checks the options database for initializations related to the 60e5c89e4eSSatish Balay PETSc components 61e5c89e4eSSatish Balay */ 6295c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Components(void) 63e5c89e4eSSatish Balay { 642d747510SLisandro Dalcin PetscBool flg; 65e5c89e4eSSatish Balay PetscErrorCode ierr; 66e5c89e4eSSatish Balay 67e5c89e4eSSatish Balay PetscFunctionBegin; 682d747510SLisandro Dalcin ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr); 692d747510SLisandro Dalcin if (flg) { 70e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 71e8e7597cSSatish Balay MPI_Comm comm = PETSC_COMM_WORLD; 72e5c89e4eSSatish Balay ierr = (*PetscHelpPrintf)(comm,"------Additional PETSc component options--------\n");CHKERRQ(ierr); 738e81d068SLisandro Dalcin ierr = (*PetscHelpPrintf)(comm," -log_exclude: <vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr); 748e81d068SLisandro Dalcin ierr = (*PetscHelpPrintf)(comm," -info_exclude: <null,vec,mat,pc,ksp,snes,tao,ts>\n");CHKERRQ(ierr); 75e5c89e4eSSatish Balay ierr = (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");CHKERRQ(ierr); 76e5c89e4eSSatish Balay #endif 77e5c89e4eSSatish Balay } 78e5c89e4eSSatish Balay PetscFunctionReturn(0); 79e5c89e4eSSatish Balay } 80e5c89e4eSSatish Balay 810f11a792SBarry Smith /* 82945d1669SBarry Smith PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args 8372a42c3cSBarry Smith 8472a42c3cSBarry Smith Collective 8572a42c3cSBarry Smith 8672a42c3cSBarry Smith Level: advanced 8772a42c3cSBarry Smith 8895452b02SPatrick Sanan Notes: 89a56f64adSBarry Smith this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to 900f11a792SBarry Smith indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to 91a56f64adSBarry Smith be called multiple times from Julia without the problem of trying to initialize MPI more than once. 920f11a792SBarry Smith 93a56f64adSBarry Smith Developer Note: Turns off PETSc signal handling to allow Julia to manage signals 941ea5a559SBarry Smith 9572a42c3cSBarry Smith .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments() 960f11a792SBarry Smith */ 97945d1669SBarry Smith PetscErrorCode PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help) 9872a42c3cSBarry Smith { 9972a42c3cSBarry Smith PetscErrorCode ierr; 10072a42c3cSBarry Smith int myargc = argc; 10172a42c3cSBarry Smith char **myargs = args; 10272a42c3cSBarry Smith 10372a42c3cSBarry Smith PetscFunctionBegin; 104c52ac268SRichard Tran Mills ierr = PetscInitialize(&myargc,&myargs,filename,help);if (ierr) return ierr; 1051ea5a559SBarry Smith ierr = PetscPopSignalHandler();CHKERRQ(ierr); 106df413903SBarry Smith PetscBeganMPI = PETSC_FALSE; 10772a42c3cSBarry Smith PetscFunctionReturn(ierr); 10872a42c3cSBarry Smith } 10972a42c3cSBarry Smith 110f0865b08SBarry Smith /* 111a56f64adSBarry Smith Used by Julia interface to get communicator 112f0865b08SBarry Smith */ 113945d1669SBarry Smith PetscErrorCode PetscGetPETSC_COMM_SELF(MPI_Comm *comm) 114f0865b08SBarry Smith { 115f0865b08SBarry Smith PetscFunctionBegin; 116f0865b08SBarry Smith *comm = PETSC_COMM_SELF; 117f0865b08SBarry Smith PetscFunctionReturn(0); 118f0865b08SBarry Smith } 119f0865b08SBarry Smith 120e5c89e4eSSatish Balay /*@C 121e5c89e4eSSatish Balay PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without 122e5c89e4eSSatish Balay the command line arguments. 123e5c89e4eSSatish Balay 124e5c89e4eSSatish Balay Collective 125e5c89e4eSSatish Balay 126e5c89e4eSSatish Balay Level: advanced 127e5c89e4eSSatish Balay 128e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeFortran() 129e5c89e4eSSatish Balay @*/ 1307087cfbeSBarry Smith PetscErrorCode PetscInitializeNoArguments(void) 131e5c89e4eSSatish Balay { 132e5c89e4eSSatish Balay PetscErrorCode ierr; 133e5c89e4eSSatish Balay int argc = 0; 134e5c89e4eSSatish Balay char **args = 0; 135e5c89e4eSSatish Balay 136e5c89e4eSSatish Balay PetscFunctionBegin; 1370298fd71SBarry Smith ierr = PetscInitialize(&argc,&args,NULL,NULL); 138e5c89e4eSSatish Balay PetscFunctionReturn(ierr); 139e5c89e4eSSatish Balay } 140e5c89e4eSSatish Balay 141e5c89e4eSSatish Balay /*@ 142e5c89e4eSSatish Balay PetscInitialized - Determine whether PETSc is initialized. 143e5c89e4eSSatish Balay 14493b6d2d1SJed Brown Level: beginner 145e5c89e4eSSatish Balay 146e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran() 147e5c89e4eSSatish Balay @*/ 1487087cfbeSBarry Smith PetscErrorCode PetscInitialized(PetscBool *isInitialized) 149e5c89e4eSSatish Balay { 150e5c89e4eSSatish Balay *isInitialized = PetscInitializeCalled; 15193b6d2d1SJed Brown return 0; 152e5c89e4eSSatish Balay } 153e5c89e4eSSatish Balay 154e5c89e4eSSatish Balay /*@ 155e5c89e4eSSatish Balay PetscFinalized - Determine whether PetscFinalize() has been called yet 156e5c89e4eSSatish Balay 157e5c89e4eSSatish Balay Level: developer 158e5c89e4eSSatish Balay 159e5c89e4eSSatish Balay .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran() 160e5c89e4eSSatish Balay @*/ 1617087cfbeSBarry Smith PetscErrorCode PetscFinalized(PetscBool *isFinalized) 162e5c89e4eSSatish Balay { 163e5c89e4eSSatish Balay *isFinalized = PetscFinalizeCalled; 16493b6d2d1SJed Brown return 0; 165e5c89e4eSSatish Balay } 166e5c89e4eSSatish Balay 16795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(void); 168e5c89e4eSSatish Balay 169e5c89e4eSSatish Balay /* 170e5c89e4eSSatish Balay This function is the MPI reduction operation used to compute the sum of the 171e5c89e4eSSatish Balay first half of the datatype and the max of the second half. 172e5c89e4eSSatish Balay */ 173367daffbSBarry Smith MPI_Op MPIU_MAXSUM_OP = 0; 174e5c89e4eSSatish Balay 175367daffbSBarry Smith PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype) 176e5c89e4eSSatish Balay { 177e5c89e4eSSatish Balay PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt; 178e5c89e4eSSatish Balay 179e5c89e4eSSatish Balay PetscFunctionBegin; 180e5c89e4eSSatish Balay if (*datatype != MPIU_2INT) { 181e5c89e4eSSatish Balay (*PetscErrorPrintf)("Can only handle MPIU_2INT data types"); 182e5c89e4eSSatish Balay MPI_Abort(MPI_COMM_WORLD,1); 183e5c89e4eSSatish Balay } 184e5c89e4eSSatish Balay 185e5c89e4eSSatish Balay for (i=0; i<count; i++) { 186e5c89e4eSSatish Balay xout[2*i] = PetscMax(xout[2*i],xin[2*i]); 187e5c89e4eSSatish Balay xout[2*i+1] += xin[2*i+1]; 188e5c89e4eSSatish Balay } 189812af9f3SBarry Smith PetscFunctionReturnVoid(); 190e5c89e4eSSatish Balay } 191e5c89e4eSSatish Balay 192e5c89e4eSSatish Balay /* 193e5c89e4eSSatish Balay Returns the max of the first entry owned by this processor and the 194e5c89e4eSSatish Balay sum of the second entry. 195b693b147SBarry Smith 19676f543a4SJed Brown The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero 197367daffbSBarry Smith is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths 198b693b147SBarry Smith there would be no place to store the both needed results. 199e5c89e4eSSatish Balay */ 20076ec1555SBarry Smith PetscErrorCode PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum) 201e5c89e4eSSatish Balay { 202e5c89e4eSSatish Balay PetscErrorCode ierr; 203e5c89e4eSSatish Balay 204e5c89e4eSSatish Balay PetscFunctionBegin; 205d6e4c47cSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 206d6e4c47cSJed Brown { 207d6e4c47cSJed Brown struct {PetscInt max,sum;} work; 208367daffbSBarry Smith ierr = MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr); 209d6e4c47cSJed Brown *max = work.max; 210d6e4c47cSJed Brown *sum = work.sum; 211d6e4c47cSJed Brown } 212d6e4c47cSJed Brown #else 213d6e4c47cSJed Brown { 214d6e4c47cSJed Brown PetscMPIInt size,rank; 215d6e4c47cSJed Brown struct {PetscInt max,sum;} *work; 216e5c89e4eSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 217e5c89e4eSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 218785e854fSJed Brown ierr = PetscMalloc1(size,&work);CHKERRQ(ierr); 219367daffbSBarry Smith ierr = MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);CHKERRQ(ierr); 2206ac3741eSJed Brown *max = work[rank].max; 2216ac3741eSJed Brown *sum = work[rank].sum; 222e5c89e4eSSatish Balay ierr = PetscFree(work);CHKERRQ(ierr); 223d6e4c47cSJed Brown } 224d6e4c47cSJed Brown #endif 225e5c89e4eSSatish Balay PetscFunctionReturn(0); 226e5c89e4eSSatish Balay } 227e5c89e4eSSatish Balay 228e5c89e4eSSatish Balay /* ----------------------------------------------------------------------------*/ 229e5c89e4eSSatish Balay 230570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 23106a205a8SBarry Smith MPI_Op MPIU_SUM = 0; 232e5c89e4eSSatish Balay 2338cc058d9SJed Brown PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 234e5c89e4eSSatish Balay { 235e5c89e4eSSatish Balay PetscInt i,count = *cnt; 236e5c89e4eSSatish Balay 237e5c89e4eSSatish Balay PetscFunctionBegin; 2387c2de775SJed Brown if (*datatype == MPIU_REAL) { 239e2e03761SBarry Smith PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 240a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] += xin[i]; 2417c2de775SJed Brown } 2427c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2437c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2447c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 245a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] += xin[i]; 2467c2de775SJed Brown } 2477c2de775SJed Brown #endif 2487c2de775SJed Brown else { 2497c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 250e2e03761SBarry Smith MPI_Abort(MPI_COMM_WORLD,1); 251e2e03761SBarry Smith } 252812af9f3SBarry Smith PetscFunctionReturnVoid(); 253e5c89e4eSSatish Balay } 254e5c89e4eSSatish Balay #endif 255e5c89e4eSSatish Balay 256570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 257d9822059SBarry Smith MPI_Op MPIU_MAX = 0; 258d9822059SBarry Smith MPI_Op MPIU_MIN = 0; 259d9822059SBarry Smith 2608cc058d9SJed Brown PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 261d9822059SBarry Smith { 262d9822059SBarry Smith PetscInt i,count = *cnt; 263d9822059SBarry Smith 264d9822059SBarry Smith PetscFunctionBegin; 2657c2de775SJed Brown if (*datatype == MPIU_REAL) { 2668c764dc5SJose Roman PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 267a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]); 2687c2de775SJed Brown } 2697c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2707c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2717c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 2727c2de775SJed Brown for (i=0; i<count; i++) { 2737c2de775SJed Brown xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 2747c2de775SJed Brown } 2757c2de775SJed Brown } 2767c2de775SJed Brown #endif 2777c2de775SJed Brown else { 2787c2de775SJed Brown (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types"); 2798c764dc5SJose Roman MPI_Abort(MPI_COMM_WORLD,1); 2808c764dc5SJose Roman } 281d9822059SBarry Smith PetscFunctionReturnVoid(); 282d9822059SBarry Smith } 283d9822059SBarry Smith 2848cc058d9SJed Brown PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype) 285d9822059SBarry Smith { 286d9822059SBarry Smith PetscInt i,count = *cnt; 287d9822059SBarry Smith 288d9822059SBarry Smith PetscFunctionBegin; 2897c2de775SJed Brown if (*datatype == MPIU_REAL) { 2908c764dc5SJose Roman PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out; 291a297a907SKarl Rupp for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]); 2927c2de775SJed Brown } 2937c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 2947c2de775SJed Brown else if (*datatype == MPIU_COMPLEX) { 2957c2de775SJed Brown PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out; 2967c2de775SJed Brown for (i=0; i<count; i++) { 2977c2de775SJed Brown xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i]; 2987c2de775SJed Brown } 2997c2de775SJed Brown } 3007c2de775SJed Brown #endif 3017c2de775SJed Brown else { 3028c764dc5SJose Roman (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types"); 3038c764dc5SJose Roman MPI_Abort(MPI_COMM_WORLD,1); 3048c764dc5SJose Roman } 305d9822059SBarry Smith PetscFunctionReturnVoid(); 306d9822059SBarry Smith } 307d9822059SBarry Smith #endif 308d9822059SBarry Smith 309480cf27aSJed Brown /* 310480cf27aSJed Brown Private routine to delete internal tag/name counter storage when a communicator is freed. 311480cf27aSJed Brown 312ff0e51ddSBarry 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. 313480cf27aSJed Brown 31412801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 315480cf27aSJed Brown 316480cf27aSJed Brown */ 3178cc058d9SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state) 318480cf27aSJed Brown { 319480cf27aSJed Brown PetscErrorCode ierr; 320480cf27aSJed Brown 321480cf27aSJed Brown PetscFunctionBegin; 32212801b39SBarry Smith ierr = PetscInfo1(0,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr); 32312801b39SBarry Smith ierr = PetscFree(count_val);CHKERRMPI(ierr); 324480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 325480cf27aSJed Brown } 326480cf27aSJed Brown 327480cf27aSJed Brown /* 32847435625SJed Brown This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user 329da3039f7SJed Brown calls MPI_Comm_free(). 330da3039f7SJed Brown 331da3039f7SJed Brown This is the only entry point for breaking the links between inner and outer comms. 332480cf27aSJed Brown 333ff0e51ddSBarry Smith This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator. 334480cf27aSJed Brown 33512801b39SBarry Smith Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval() 336480cf27aSJed Brown 337480cf27aSJed Brown */ 338da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state) 339480cf27aSJed Brown { 340480cf27aSJed Brown PetscErrorCode ierr; 341b89831e5SBarry Smith PetscMPIInt flg; 342265f3f35SJed Brown union {MPI_Comm comm; void *ptr;} icomm,ocomm; 343480cf27aSJed Brown 344480cf27aSJed Brown PetscFunctionBegin; 34512801b39SBarry Smith if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval"); 346ec4fadc2SJed Brown icomm.ptr = attr_val; 347da3039f7SJed Brown 34847435625SJed Brown ierr = MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);CHKERRMPI(ierr); 34912801b39SBarry Smith if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm"); 35012801b39SBarry Smith if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm"); 35147435625SJed Brown ierr = MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);CHKERRMPI(ierr); 35212801b39SBarry Smith ierr = PetscInfo1(0,"User MPI_Comm %ld is being freed after removing reference from inner PETSc comm to this outer comm\n",(long)comm);CHKERRMPI(ierr); 353da3039f7SJed Brown PetscFunctionReturn(MPI_SUCCESS); 354b89831e5SBarry Smith } 355da3039f7SJed Brown 356da3039f7SJed Brown /* 35747435625SJed Brown * This is invoked on the inner comm when Petsc_DelComm_Outer calls MPI_Comm_delete_attr. It should not be reached any other way. 358da3039f7SJed Brown */ 359da3039f7SJed Brown PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state) 360da3039f7SJed Brown { 361da3039f7SJed Brown PetscErrorCode ierr; 362da3039f7SJed Brown 363da3039f7SJed Brown PetscFunctionBegin; 36412801b39SBarry Smith ierr = PetscInfo1(0,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);CHKERRMPI(ierr); 365480cf27aSJed Brown PetscFunctionReturn(MPI_SUCCESS); 366480cf27aSJed Brown } 367480cf27aSJed Brown 3685f7487a0SJunchao Zhang PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Shm(MPI_Comm,PetscMPIInt,void *,void *); 36942218b76SBarry Smith 370951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 3718cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*); 3728cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*); 3738cc058d9SJed Brown PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*); 374e39fd77fSBarry Smith #endif 375e39fd77fSBarry Smith 37612801b39SBarry Smith PetscMPIInt PETSC_MPI_ERROR_CLASS,PETSC_MPI_ERROR_CODE; 37712801b39SBarry Smith 378eb27c7c8SSatish Balay PETSC_INTERN int PetscGlobalArgc; 379eb27c7c8SSatish Balay PETSC_INTERN char **PetscGlobalArgs; 3806ae9a8a6SBarry Smith int PetscGlobalArgc = 0; 3816ae9a8a6SBarry Smith char **PetscGlobalArgs = 0; 382dff31646SBarry Smith PetscSegBuffer PetscCitationsList; 383e5c89e4eSSatish Balay 384dd63322aSSatish Balay PetscErrorCode PetscCitationsInitialize(void) 385051e4cf2SJed Brown { 386051e4cf2SJed Brown PetscErrorCode ierr; 387051e4cf2SJed Brown 388051e4cf2SJed Brown PetscFunctionBegin; 389051e4cf2SJed Brown ierr = PetscSegBufferCreate(1,10000,&PetscCitationsList);CHKERRQ(ierr); 390a1601671SBarry 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); 391051e4cf2SJed 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); 392051e4cf2SJed Brown PetscFunctionReturn(0); 393051e4cf2SJed Brown } 394e5c89e4eSSatish Balay 3952d747510SLisandro Dalcin static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */ 3962d747510SLisandro Dalcin 3972d747510SLisandro Dalcin PetscErrorCode PetscSetProgramName(const char name[]) 3982d747510SLisandro Dalcin { 3992d747510SLisandro Dalcin PetscErrorCode ierr; 4002d747510SLisandro Dalcin 4012d747510SLisandro Dalcin PetscFunctionBegin; 4022d747510SLisandro Dalcin ierr = PetscStrncpy(programname,name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 4032d747510SLisandro Dalcin PetscFunctionReturn(0); 4042d747510SLisandro Dalcin } 4052d747510SLisandro Dalcin 4062d747510SLisandro Dalcin /*@C 4072d747510SLisandro Dalcin PetscGetProgramName - Gets the name of the running program. 4082d747510SLisandro Dalcin 4092d747510SLisandro Dalcin Not Collective 4102d747510SLisandro Dalcin 4112d747510SLisandro Dalcin Input Parameter: 4122d747510SLisandro Dalcin . len - length of the string name 4132d747510SLisandro Dalcin 4142d747510SLisandro Dalcin Output Parameter: 4152d747510SLisandro Dalcin . name - the name of the running program 4162d747510SLisandro Dalcin 4172d747510SLisandro Dalcin Level: advanced 4182d747510SLisandro Dalcin 4192d747510SLisandro Dalcin Notes: 4202d747510SLisandro Dalcin The name of the program is copied into the user-provided character 4212d747510SLisandro Dalcin array of length len. On some machines the program name includes 4222d747510SLisandro Dalcin its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN. 4232d747510SLisandro Dalcin @*/ 4242d747510SLisandro Dalcin PetscErrorCode PetscGetProgramName(char name[],size_t len) 4252d747510SLisandro Dalcin { 4262d747510SLisandro Dalcin PetscErrorCode ierr; 4272d747510SLisandro Dalcin 4282d747510SLisandro Dalcin PetscFunctionBegin; 4292d747510SLisandro Dalcin ierr = PetscStrncpy(name,programname,len);CHKERRQ(ierr); 4302d747510SLisandro Dalcin PetscFunctionReturn(0); 4312d747510SLisandro Dalcin } 4322d747510SLisandro Dalcin 433e5c89e4eSSatish Balay /*@C 434e5c89e4eSSatish Balay PetscGetArgs - Allows you to access the raw command line arguments anywhere 435e5c89e4eSSatish Balay after PetscInitialize() is called but before PetscFinalize(). 436e5c89e4eSSatish Balay 437e5c89e4eSSatish Balay Not Collective 438e5c89e4eSSatish Balay 439e5c89e4eSSatish Balay Output Parameters: 440e5c89e4eSSatish Balay + argc - count of number of command line arguments 441e5c89e4eSSatish Balay - args - the command line arguments 442e5c89e4eSSatish Balay 443e5c89e4eSSatish Balay Level: intermediate 444e5c89e4eSSatish Balay 445e5c89e4eSSatish Balay Notes: 446e5c89e4eSSatish Balay This is usually used to pass the command line arguments into other libraries 447e5c89e4eSSatish Balay that are called internally deep in PETSc or the application. 448e5c89e4eSSatish Balay 449f177e3b1SBarry Smith The first argument contains the program name as is normal for C arguments. 450f177e3b1SBarry Smith 451793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments() 452e5c89e4eSSatish Balay 453e5c89e4eSSatish Balay @*/ 4547087cfbeSBarry Smith PetscErrorCode PetscGetArgs(int *argc,char ***args) 455e5c89e4eSSatish Balay { 456e5c89e4eSSatish Balay PetscFunctionBegin; 45717186662SBarry Smith if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()"); 458e5c89e4eSSatish Balay *argc = PetscGlobalArgc; 459e5c89e4eSSatish Balay *args = PetscGlobalArgs; 460e5c89e4eSSatish Balay PetscFunctionReturn(0); 461e5c89e4eSSatish Balay } 462e5c89e4eSSatish Balay 463793721a6SBarry Smith /*@C 464793721a6SBarry Smith PetscGetArguments - Allows you to access the command line arguments anywhere 465793721a6SBarry Smith after PetscInitialize() is called but before PetscFinalize(). 466793721a6SBarry Smith 467793721a6SBarry Smith Not Collective 468793721a6SBarry Smith 469793721a6SBarry Smith Output Parameters: 470793721a6SBarry Smith . args - the command line arguments 471793721a6SBarry Smith 472793721a6SBarry Smith Level: intermediate 473793721a6SBarry Smith 474793721a6SBarry Smith Notes: 475793721a6SBarry Smith This does NOT start with the program name and IS null terminated (final arg is void) 476793721a6SBarry Smith 477793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments() 478793721a6SBarry Smith 479793721a6SBarry Smith @*/ 4807087cfbeSBarry Smith PetscErrorCode PetscGetArguments(char ***args) 481793721a6SBarry Smith { 482793721a6SBarry Smith PetscInt i,argc = PetscGlobalArgc; 483793721a6SBarry Smith PetscErrorCode ierr; 484793721a6SBarry Smith 485793721a6SBarry Smith PetscFunctionBegin; 48617186662SBarry Smith if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()"); 4872d747510SLisandro Dalcin if (!argc) {*args = NULL; PetscFunctionReturn(0);} 488785e854fSJed Brown ierr = PetscMalloc1(argc,args);CHKERRQ(ierr); 489793721a6SBarry Smith for (i=0; i<argc-1; i++) { 490793721a6SBarry Smith ierr = PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);CHKERRQ(ierr); 491793721a6SBarry Smith } 4922d747510SLisandro Dalcin (*args)[argc-1] = NULL; 493793721a6SBarry Smith PetscFunctionReturn(0); 494793721a6SBarry Smith } 495793721a6SBarry Smith 496793721a6SBarry Smith /*@C 497793721a6SBarry Smith PetscFreeArguments - Frees the memory obtained with PetscGetArguments() 498793721a6SBarry Smith 499793721a6SBarry Smith Not Collective 500793721a6SBarry Smith 501793721a6SBarry Smith Output Parameters: 502793721a6SBarry Smith . args - the command line arguments 503793721a6SBarry Smith 504793721a6SBarry Smith Level: intermediate 505793721a6SBarry Smith 506793721a6SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments() 507793721a6SBarry Smith 508793721a6SBarry Smith @*/ 5097087cfbeSBarry Smith PetscErrorCode PetscFreeArguments(char **args) 510793721a6SBarry Smith { 511793721a6SBarry Smith PetscInt i = 0; 512793721a6SBarry Smith PetscErrorCode ierr; 513793721a6SBarry Smith 514793721a6SBarry Smith PetscFunctionBegin; 515a297a907SKarl Rupp if (!args) PetscFunctionReturn(0); 516793721a6SBarry Smith while (args[i]) { 517793721a6SBarry Smith ierr = PetscFree(args[i]);CHKERRQ(ierr); 518793721a6SBarry Smith i++; 519793721a6SBarry Smith } 520793721a6SBarry Smith ierr = PetscFree(args);CHKERRQ(ierr); 521793721a6SBarry Smith PetscFunctionReturn(0); 522793721a6SBarry Smith } 523793721a6SBarry Smith 52411525c0dSBarry Smith #if defined(PETSC_HAVE_SAWS) 52530befbd2SBarry Smith #include <petscconfiginfo.h> 52630befbd2SBarry Smith 52795c0884eSLisandro Dalcin PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[]) 52811525c0dSBarry Smith { 52911525c0dSBarry Smith if (!PetscGlobalRank) { 53030befbd2SBarry Smith char cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64]; 53111525c0dSBarry Smith int port; 532ffbd1cfbSBarry Smith PetscBool flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE; 53311525c0dSBarry Smith size_t applinelen,introlen; 53411525c0dSBarry Smith PetscErrorCode ierr; 535ffbd1cfbSBarry Smith char sawsurl[256]; 53611525c0dSBarry Smith 537c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);CHKERRQ(ierr); 53811525c0dSBarry Smith if (flg) { 53911525c0dSBarry Smith char sawslog[PETSC_MAX_PATH_LEN]; 54011525c0dSBarry Smith 541c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 54211525c0dSBarry Smith if (sawslog[0]) { 54311525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog)); 54411525c0dSBarry Smith } else { 54511525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL)); 54611525c0dSBarry Smith } 54711525c0dSBarry Smith } 548c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 54911525c0dSBarry Smith if (flg) { 55011525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert)); 55111525c0dSBarry Smith } 552c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);CHKERRQ(ierr); 553ffbd1cfbSBarry Smith if (selectport) { 554ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Get_Available_Port,(&port)); 555ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 556ffbd1cfbSBarry Smith } else { 557c5929fdfSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);CHKERRQ(ierr); 55811525c0dSBarry Smith if (flg) { 55911525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 56011525c0dSBarry Smith } 561ffbd1cfbSBarry Smith } 562c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 56311525c0dSBarry Smith if (flg) { 56411525c0dSBarry Smith PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr); 56511525c0dSBarry Smith ierr = PetscStrcmp(root,".",&rootlocal);CHKERRQ(ierr); 5669c1e0ce8SBarry Smith } else { 567c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);CHKERRQ(ierr); 5689c1e0ce8SBarry Smith if (flg) { 5693c01dfcfSBarry Smith ierr = PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 5709c1e0ce8SBarry Smith PetscStackCallSAWs(SAWs_Set_Document_Root,(root));CHKERRQ(ierr); 5719c1e0ce8SBarry Smith } 57211525c0dSBarry Smith } 573c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);CHKERRQ(ierr); 57411525c0dSBarry Smith if (flg2) { 57511525c0dSBarry Smith char jsdir[PETSC_MAX_PATH_LEN]; 57611525c0dSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option"); 57711525c0dSBarry Smith ierr = PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);CHKERRQ(ierr); 57811525c0dSBarry Smith ierr = PetscTestDirectory(jsdir,'r',&flg);CHKERRQ(ierr); 57911525c0dSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory"); 58043da4ab2SBarry Smith PetscStackCallSAWs(SAWs_Push_Local_Header,());CHKERRQ(ierr); 58111525c0dSBarry Smith } 58211525c0dSBarry Smith ierr = PetscGetProgramName(programname,64);CHKERRQ(ierr); 58311525c0dSBarry Smith ierr = PetscStrlen(help,&applinelen);CHKERRQ(ierr); 58411525c0dSBarry Smith introlen = 4096 + applinelen; 58530a8c9c0SSurtai Han applinelen += 1024; 58611525c0dSBarry Smith ierr = PetscMalloc(applinelen,&appline);CHKERRQ(ierr); 58711525c0dSBarry Smith ierr = PetscMalloc(introlen,&intro);CHKERRQ(ierr); 58811525c0dSBarry Smith 58911525c0dSBarry Smith if (rootlocal) { 59011525c0dSBarry Smith ierr = PetscSNPrintf(appline,applinelen,"%s.c.html",programname);CHKERRQ(ierr); 59111525c0dSBarry Smith ierr = PetscTestFile(appline,'r',&rootlocal);CHKERRQ(ierr); 59211525c0dSBarry Smith } 59376a34f28SBarry Smith ierr = PetscOptionsGetAll(NULL,&options);CHKERRQ(ierr); 59411525c0dSBarry Smith if (rootlocal && help) { 595928bb9adSStefano 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); 59611525c0dSBarry Smith } else if (help) { 597928bb9adSStefano Zampini ierr = PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);CHKERRQ(ierr); 59811525c0dSBarry Smith } else { 599928bb9adSStefano Zampini ierr = PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);CHKERRQ(ierr); 60011525c0dSBarry Smith } 601b0bb5815SBarry Smith ierr = PetscFree(options);CHKERRQ(ierr); 60230befbd2SBarry Smith ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr); 60311525c0dSBarry Smith ierr = PetscSNPrintf(intro,introlen,"<body>\n" 604a8d69d7bSBarry 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" 605df62c222SBarry 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" 606928bb9adSStefano Zampini "%s",version,petscconfigureoptions,appline);CHKERRQ(ierr); 60743da4ab2SBarry Smith PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro)); 60811525c0dSBarry Smith ierr = PetscFree(intro);CHKERRQ(ierr); 60911525c0dSBarry Smith ierr = PetscFree(appline);CHKERRQ(ierr); 610ffbd1cfbSBarry Smith if (selectport) { 611aa573868SBarry Smith PetscBool silent; 6127d812c46SBarry Smith 6137d812c46SBarry Smith ierr = SAWs_Initialize(); 6147d812c46SBarry Smith /* another process may have grabbed the port so keep trying */ 6157d812c46SBarry Smith while (ierr) { 6167d812c46SBarry Smith PetscStackCallSAWs(SAWs_Get_Available_Port,(&port)); 6177d812c46SBarry Smith PetscStackCallSAWs(SAWs_Set_Port,(port)); 6187d812c46SBarry Smith ierr = SAWs_Initialize(); 6197d812c46SBarry Smith } 6207d812c46SBarry Smith 621aa573868SBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);CHKERRQ(ierr); 622aa573868SBarry Smith if (!silent) { 623ffbd1cfbSBarry Smith PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl)); 624ffbd1cfbSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);CHKERRQ(ierr); 625ffbd1cfbSBarry Smith } 6267d812c46SBarry Smith } else { 6277d812c46SBarry Smith PetscStackCallSAWs(SAWs_Initialize,()); 628aa573868SBarry Smith } 6290af79b04SBarry Smith ierr = PetscCitationsRegister("@TechReport{ saws,\n" 6300af79b04SBarry Smith " Author = {Matt Otten and Jed Brown and Barry Smith},\n" 6310af79b04SBarry Smith " Title = {Scientific Application Web Server (SAWs) Users Manual},\n" 6320af79b04SBarry Smith " Institution = {Argonne National Laboratory},\n" 6330af79b04SBarry Smith " Year = 2013\n}\n",NULL);CHKERRQ(ierr); 63411525c0dSBarry Smith } 63511525c0dSBarry Smith PetscFunctionReturn(0); 63611525c0dSBarry Smith } 63711525c0dSBarry Smith #endif 63811525c0dSBarry Smith 6394dfee713SSatish Balay /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */ 6404dfee713SSatish Balay PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void) 6414dfee713SSatish Balay { 6424dfee713SSatish Balay PetscFunctionBegin; 6434dfee713SSatish Balay #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG) 6444dfee713SSatish Balay /* see MPI.py for details on this bug */ 6454dfee713SSatish Balay (void) setenv("HWLOC_COMPONENTS","-x86",1); 6464dfee713SSatish Balay #endif 6474dfee713SSatish Balay PetscFunctionReturn(0); 6484dfee713SSatish Balay } 6494dfee713SSatish Balay 650a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 651a56f64adSBarry Smith #include <adios.h> 65222580e64SBarry Smith #include <adios_read.h> 6537b56e58cSSatish Balay int64_t Petsc_adios_group; 654a56f64adSBarry Smith #endif 65555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2) 65655d657eeSBarry Smith #include <adios2_c.h> 65755d657eeSBarry Smith #endif 658cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP) 659cd1b99a6SBarry Smith #include <omp.h> 660f90b075cSBarry Smith PetscInt PetscNumOMPThreads; 661cd1b99a6SBarry Smith #endif 662a56f64adSBarry Smith 663e5c89e4eSSatish Balay /*@C 664e5c89e4eSSatish Balay PetscInitialize - Initializes the PETSc database and MPI. 665e5c89e4eSSatish Balay PetscInitialize() calls MPI_Init() if that has yet to be called, 666e5c89e4eSSatish Balay so this routine should always be called near the beginning of 667e5c89e4eSSatish Balay your program -- usually the very first line! 668e5c89e4eSSatish Balay 669e5c89e4eSSatish Balay Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set 670e5c89e4eSSatish Balay 671e5c89e4eSSatish Balay Input Parameters: 672e5c89e4eSSatish Balay + argc - count of number of command line arguments 673e5c89e4eSSatish Balay . args - the command line arguments 6740298fd71SBarry Smith . file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for 675fc2bca9aSBarry Smith code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files 6760298fd71SBarry Smith - help - [optional] Help message to print, use NULL for no message 677e5c89e4eSSatish Balay 67805827820SBarry Smith If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that 67905827820SBarry Smith communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a 68005827820SBarry Smith four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not, 68105827820SBarry Smith then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even 68205827820SBarry Smith if different subcommunicators of the job are doing different things with PETSc. 683e5c89e4eSSatish Balay 684e5c89e4eSSatish Balay Options Database Keys: 6857ca660e7SBarry Smith + -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message 6867ca660e7SBarry Smith . -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger 687e5c89e4eSSatish Balay . -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected 6888a690491SBarry Smith . -on_error_emacs <machinename> - causes emacsclient to jump to error file 6898a690491SBarry Smith . -on_error_abort - calls abort() when error detected (no traceback) 6908a690491SBarry Smith . -on_error_mpiabort - calls MPI_abort() when error detected 6918a690491SBarry Smith . -error_output_stderr - prints error messages to stderr instead of the default stdout 6928a690491SBarry Smith . -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called) 693e5c89e4eSSatish Balay . -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger 694e5c89e4eSSatish Balay . -debugger_pause [sleeptime] (in seconds) - Pauses debugger 695e5c89e4eSSatish Balay . -stop_for_debugger - Print message on how to attach debugger manually to 696e5c89e4eSSatish Balay process and wait (-debugger_pause) seconds for attachment 697*79dccf82SBarry Smith . -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug) 698*79dccf82SBarry Smith . -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no) 699*79dccf82SBarry Smith . -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug() 700aee23540SBarry Smith . -malloc_dump - prints a list of all unfreed memory at the end of the run 70192f119d6SBarry 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 70292f119d6SBarry Smith . -malloc_view - show a list of all allocated memory during PetscFinalize() 70392f119d6SBarry Smith . -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view 70492f119d6SBarry Smith . -fp_trap - Stops on floating point exceptions 705e5c89e4eSSatish Balay . -no_signal_handler - Indicates not to trap error signals 706e5c89e4eSSatish Balay . -shared_tmp - indicates /tmp directory is shared by all processors 707e5c89e4eSSatish Balay . -not_shared_tmp - each processor has own /tmp 708e5c89e4eSSatish Balay . -tmp - alternative name of /tmp directory 709e5c89e4eSSatish Balay . -get_total_flops - returns total flops done by all processors 7100841954dSBarry Smith - -memory_view - Print memory usage at end of run 711e5c89e4eSSatish Balay 712e5c89e4eSSatish Balay Options Database Keys for Profiling: 713a7f22e61SSatish Balay See Users-Manual: ch_profiling for details. 714495fc317SBarry Smith + -info <optional filename> - Prints verbose information to the screen 715495fc317SBarry Smith . -info_exclude <null,vec,mat,pc,ksp,snes,ts> - Excludes some of the verbose messages 716217044c2SLisandro Dalcin . -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event, 717217044c2SLisandro Dalcin however it slows things down and gives a distorted view of the overall runtime. 718495fc317SBarry Smith . -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program 719e5c89e4eSSatish Balay hangs without running in the debugger). See PetscLogTraceBegin(). 7209a9a5d4cSBarry Smith . -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView(). 721*79dccf82SBarry Smith . -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView(). 7229a9a5d4cSBarry Smith . -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the 723495fc317SBarry Smith summary is written to the file. See PetscLogView(). 724f5d6ab90SLisandro Dalcin . -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging 725495fc317SBarry Smith . -log_all [filename] - Logs extensive profiling information See PetscLogDump(). 726495fc317SBarry Smith . -log [filename] - Logs basic profiline information See PetscLogDump(). 727c2f74817SBarry Smith . -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution) 72887c3beb6SLisandro Dalcin . -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off 729c2f74817SBarry Smith - -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code 730495fc317SBarry Smith 731609bdbeeSBarry Smith Only one of -log_trace, -log_view, -log_view, -log_all, -log, or -log_mpe may be used at a time 732e5c89e4eSSatish Balay 733ffbd1cfbSBarry Smith Options Database Keys for SAWs: 734ffbd1cfbSBarry Smith + -saws_port <portnumber> - port number to publish SAWs data, default is 8080 735ffbd1cfbSBarry 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 736ffbd1cfbSBarry Smith this is useful when you are running many jobs that utilize SAWs at the same time 737ffbd1cfbSBarry Smith . -saws_log <filename> - save a log of all SAWs communication 738ffbd1cfbSBarry Smith . -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP 739ffbd1cfbSBarry Smith - -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files 740ffbd1cfbSBarry Smith 741e5c89e4eSSatish Balay Environmental Variables: 742e5c89e4eSSatish Balay + PETSC_TMP - alternative tmp directory 743e5c89e4eSSatish Balay . PETSC_SHARED_TMP - tmp is shared by all processes 744e5c89e4eSSatish Balay . PETSC_NOT_SHARED_TMP - each process has its own private tmp 745e5c89e4eSSatish Balay . PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer 746e5c89e4eSSatish Balay - PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to 747e5c89e4eSSatish Balay 748e5c89e4eSSatish Balay 749e5c89e4eSSatish Balay Level: beginner 750e5c89e4eSSatish Balay 751e5c89e4eSSatish Balay Notes: 752e5c89e4eSSatish Balay If for some reason you must call MPI_Init() separately, call 753e5c89e4eSSatish Balay it before PetscInitialize(). 754e5c89e4eSSatish Balay 755e5c89e4eSSatish Balay Fortran Version: 756e5c89e4eSSatish Balay In Fortran this routine has the format 757e5c89e4eSSatish Balay $ call PetscInitialize(file,ierr) 758e5c89e4eSSatish Balay 759e5c89e4eSSatish Balay + ierr - error return code 7600eb4c9c0SBarry Smith - file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for 761fc2bca9aSBarry Smith code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files 762e5c89e4eSSatish Balay 763e5c89e4eSSatish Balay Important Fortran Note: 7640eb4c9c0SBarry Smith In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a 7650298fd71SBarry Smith null character string; you CANNOT just use NULL as 766a7f22e61SSatish Balay in the C version. See Users-Manual: ch_fortran for details. 767e5c89e4eSSatish Balay 76801cb0274SBarry Smith If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after 76901cb0274SBarry Smith calling PetscInitialize(). 770e5c89e4eSSatish Balay 77101cb0274SBarry Smith .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments() 772e5c89e4eSSatish Balay 773e5c89e4eSSatish Balay @*/ 7747087cfbeSBarry Smith PetscErrorCode PetscInitialize(int *argc,char ***args,const char file[],const char help[]) 775e5c89e4eSSatish Balay { 776e5c89e4eSSatish Balay PetscErrorCode ierr; 7774bb5149bSJed Brown PetscMPIInt flag, size; 77819c5658aSBarry Smith PetscBool flg = PETSC_TRUE; 779e5c89e4eSSatish Balay char hostname[256]; 780e5c89e4eSSatish Balay 781e5c89e4eSSatish Balay PetscFunctionBegin; 782e5c89e4eSSatish Balay if (PetscInitializeCalled) PetscFunctionReturn(0); 7833d96e996SBarry Smith /* 7843d96e996SBarry Smith The checking over compatible runtime libraries is complicated by the MPI ABI initiative 7853d96e996SBarry Smith https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with 7863d96e996SBarry Smith MPICH v3.1 (Released Feburary 2014) 7873d96e996SBarry Smith IBM MPI v2.1 (December 2014) 7883d96e996SBarry Smith Intel® MPI Library v5.0 (2014) 7893d96e996SBarry Smith Cray MPT v7.0.0 (June 2014) 7903d96e996SBarry Smith As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions 7913d96e996SBarry Smith listed above and since that time are compatible. 7923d96e996SBarry Smith 7933d96e996SBarry Smith Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number 7943d96e996SBarry Smith at compile time or runtime. Thus we will need to systematically track the allowed versions 7953d96e996SBarry Smith and how they are represented in the mpi.h and MPI_Get_library_version() output in order 7963d96e996SBarry Smith to perform the checking. 7973d96e996SBarry Smith 7983d96e996SBarry Smith Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI). 7993d96e996SBarry Smith 8003d96e996SBarry Smith Questions: 8013d96e996SBarry Smith 8023d96e996SBarry Smith Should the checks for ABI incompatibility be only on the major version number below? 8033d96e996SBarry Smith Presumably the output to stderr will be removed before a release. 8043d96e996SBarry Smith */ 8053d96e996SBarry Smith 80619c5658aSBarry Smith #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION) 80719c5658aSBarry Smith { 80819c5658aSBarry Smith char mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING]; 80919c5658aSBarry Smith PetscMPIInt mpilibraryversionlength; 81019c5658aSBarry Smith ierr = MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);if (ierr) return ierr; 8113d96e996SBarry Smith /* check for MPICH versions before MPI ABI initiative */ 81219c5658aSBarry Smith #if defined(MPICH_VERSION) 8133d96e996SBarry Smith #if MPICH_NUMVERSION < 30100000 81419c5658aSBarry Smith { 81519c5658aSBarry Smith char *ver,*lf; 81619c5658aSBarry Smith flg = PETSC_FALSE; 81719c5658aSBarry Smith ierr = PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);if (ierr) return ierr; 81819c5658aSBarry Smith if (ver) { 81919c5658aSBarry Smith ierr = PetscStrchr(ver,'\n',&lf);if (ierr) return ierr; 82019c5658aSBarry Smith if (lf) { 82119c5658aSBarry Smith *lf = 0; 82219c5658aSBarry Smith ierr = PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) return ierr; 82319c5658aSBarry Smith } 82419c5658aSBarry Smith } 82519c5658aSBarry Smith if (!flg) { 82619c5658aSBarry Smith fprintf(stderr,"PETSc Error --- MPICH library version \n%s does not match what PETSc was compiled with %s, aborting\n",mpilibraryversion,MPICH_VERSION); 8273d96e996SBarry Smith return PETSC_ERR_MPI_LIB_INCOMP; 82819c5658aSBarry Smith } 82919c5658aSBarry Smith } 8303d96e996SBarry Smith #endif 8313d96e996SBarry 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?) */ 83219c5658aSBarry Smith #elif defined(OMPI_MAJOR_VERSION) 83319c5658aSBarry Smith { 83419c5658aSBarry Smith char *ver,bs[32],*bsf; 83519c5658aSBarry Smith flg = PETSC_FALSE; 83619c5658aSBarry Smith ierr = PetscStrstr(mpilibraryversion,"Open MPI",&ver);if (ierr) return ierr; 83719c5658aSBarry Smith if (ver) { 8382e924ca5SSatish Balay PetscSNPrintf(bs,32,"v%d.%d",OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION); 83919c5658aSBarry Smith ierr = PetscStrstr(ver,bs,&bsf);if (ierr) return ierr; 84019c5658aSBarry Smith if (bsf) flg = PETSC_TRUE; 84119c5658aSBarry Smith } 84219c5658aSBarry Smith if (!flg) { 84319c5658aSBarry 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); 8443d96e996SBarry Smith return PETSC_ERR_MPI_LIB_INCOMP; 84519c5658aSBarry Smith } 84619c5658aSBarry Smith } 84719c5658aSBarry Smith #endif 84819c5658aSBarry Smith } 84919c5658aSBarry Smith #endif 85019c5658aSBarry Smith 851ae9b4142SLisandro Dalcin /* these must be initialized in a routine, not as a constant declaration*/ 852d89683f4Sbcordonn PETSC_STDOUT = stdout; 853ae9b4142SLisandro Dalcin PETSC_STDERR = stderr; 854e5c89e4eSSatish Balay 8550c30907bSSatish Balay /* on Windows - set printf to default to printing 2 digit exponents */ 8560c30907bSSatish Balay #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT) 8570c30907bSSatish Balay _set_output_format(_TWO_DIGIT_EXPONENT); 8580c30907bSSatish Balay #endif 8590c30907bSSatish Balay 8604416b707SBarry Smith ierr = PetscOptionsCreateDefault();CHKERRQ(ierr); 861e5c89e4eSSatish Balay 862e5c89e4eSSatish Balay /* 863e5c89e4eSSatish Balay We initialize the program name here (before MPI_Init()) because MPICH has a bug in 864e5c89e4eSSatish Balay it that it sets args[0] on all processors to be args[0] on the first processor. 865e5c89e4eSSatish Balay */ 866e5c89e4eSSatish Balay if (argc && *argc) { 867e5c89e4eSSatish Balay ierr = PetscSetProgramName(**args);CHKERRQ(ierr); 868e5c89e4eSSatish Balay } else { 869e5c89e4eSSatish Balay ierr = PetscSetProgramName("Unknown Name");CHKERRQ(ierr); 870e5c89e4eSSatish Balay } 871e5c89e4eSSatish Balay 872e5c89e4eSSatish Balay ierr = MPI_Initialized(&flag);CHKERRQ(ierr); 873e5c89e4eSSatish Balay if (!flag) { 874e32f2f54SBarry 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"); 8754dfee713SSatish Balay ierr = PetscPreMPIInit_Private();CHKERRQ(ierr); 8765e765c61SJed Brown #if defined(PETSC_HAVE_MPI_INIT_THREAD) 8775e765c61SJed Brown { 8785e765c61SJed Brown PetscMPIInt provided; 8795e765c61SJed Brown ierr = MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);CHKERRQ(ierr); 8805e765c61SJed Brown } 8815e765c61SJed Brown #else 882e5c89e4eSSatish Balay ierr = MPI_Init(argc,args);CHKERRQ(ierr); 8835e765c61SJed Brown #endif 884e5c89e4eSSatish Balay PetscBeganMPI = PETSC_TRUE; 885e5c89e4eSSatish Balay } 8864dfee713SSatish Balay 887e5c89e4eSSatish Balay if (argc && args) { 888e5c89e4eSSatish Balay PetscGlobalArgc = *argc; 889e5c89e4eSSatish Balay PetscGlobalArgs = *args; 890e5c89e4eSSatish Balay } 891e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_FALSE; 8925ad9ad5bSBarry Smith ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr); 8935ad9ad5bSBarry Smith ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr); 8945ad9ad5bSBarry Smith ierr = PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr); 895ef19f930SBarry Smith ierr = PetscSpinlockCreate(&PetscCommSpinLock);CHKERRQ(ierr); 896e5c89e4eSSatish Balay 897a297a907SKarl Rupp if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD; 898d54338ecSKarl Rupp ierr = MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);CHKERRQ(ierr); 899e5c89e4eSSatish Balay 90012801b39SBarry Smith ierr = MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);CHKERRQ(ierr); 90112801b39SBarry Smith ierr = MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);CHKERRQ(ierr); 90212801b39SBarry Smith 903e5c89e4eSSatish Balay /* Done after init due to a bug in MPICH-GM? */ 904e5c89e4eSSatish Balay ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr); 905e5c89e4eSSatish Balay 906e5c89e4eSSatish Balay ierr = MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);CHKERRQ(ierr); 907e5c89e4eSSatish Balay ierr = MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);CHKERRQ(ierr); 908e5c89e4eSSatish Balay 9098ad47952SJed Brown MPIU_BOOL = MPI_INT; 9108ad47952SJed Brown MPIU_ENUM = MPI_INT; 9117cdaf61dSJed Brown MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64; 912e316c87fSJed Brown if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED; 913e316c87fSJed Brown else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG; 914e316c87fSJed Brown #if defined(PETSC_SIZEOF_LONG_LONG) 915e316c87fSJed Brown else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG; 916e316c87fSJed Brown #endif 917e316c87fSJed Brown else {(*PetscErrorPrintf)("PetscInitialize: Could not find MPI type for size_t\n"); return PETSC_ERR_SUP_SYS;} 9188ad47952SJed Brown 919e5c89e4eSSatish Balay /* 920e5c89e4eSSatish Balay Initialized the global complex variable; this is because with 921e5c89e4eSSatish Balay shared libraries the constructors for global variables 922e5c89e4eSSatish Balay are not called; at least on IRIX. 923e5c89e4eSSatish Balay */ 924886cfec0SSatish Balay #if defined(PETSC_HAVE_COMPLEX) 925e5c89e4eSSatish Balay { 926252ecd31SSatish Balay #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128) 92750f81f78SJed Brown PetscComplex ic(0.0,1.0); 928e5c89e4eSSatish Balay PETSC_i = ic; 929252ecd31SSatish Balay #else 93050f81f78SJed Brown PETSC_i = _Complex_I; 931b7940d39SSatish Balay #endif 932762437b8SSatish Balay } 933762437b8SSatish Balay 9342c876bd9SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX) 935e69cd0e6SSatish Balay ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr); 936500d8756SSatish Balay ierr = MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr); 937500d8756SSatish Balay ierr = MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);CHKERRQ(ierr); 938500d8756SSatish Balay ierr = MPI_Type_commit(&MPIU_C_COMPLEX);CHKERRQ(ierr); 9392c876bd9SBarry Smith #endif 940886cfec0SSatish Balay #endif /* PETSC_HAVE_COMPLEX */ 941e5c89e4eSSatish Balay 942e5c89e4eSSatish Balay /* 943e5c89e4eSSatish Balay Create the PETSc MPI reduction operator that sums of the first 944e5c89e4eSSatish Balay half of the entries and maxes the second half. 945e5c89e4eSSatish Balay */ 946367daffbSBarry Smith ierr = MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);CHKERRQ(ierr); 947e5c89e4eSSatish Balay 948ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 949c90a1750SBarry Smith ierr = MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);CHKERRQ(ierr); 950c90a1750SBarry Smith ierr = MPI_Type_commit(&MPIU___FLOAT128);CHKERRQ(ierr); 9517c2de775SJed Brown #if defined(PETSC_HAVE_COMPLEX) 9528c764dc5SJose Roman ierr = MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);CHKERRQ(ierr); 9538c764dc5SJose Roman ierr = MPI_Type_commit(&MPIU___COMPLEX128);CHKERRQ(ierr); 9548c764dc5SJose Roman #endif 955d9822059SBarry Smith ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr); 956d9822059SBarry Smith ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr); 957570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 958570b7f6dSBarry Smith ierr = MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);CHKERRQ(ierr); 959570b7f6dSBarry Smith ierr = MPI_Type_commit(&MPIU___FP16);CHKERRQ(ierr); 960570b7f6dSBarry Smith ierr = MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);CHKERRQ(ierr); 961570b7f6dSBarry Smith ierr = MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);CHKERRQ(ierr); 962c90a1750SBarry Smith #endif 963c90a1750SBarry Smith 964570b7f6dSBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 965cca4cb22SSatish Balay ierr = MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);CHKERRQ(ierr); 966cca4cb22SSatish Balay #endif 967cca4cb22SSatish Balay 968e5c89e4eSSatish Balay ierr = MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);CHKERRQ(ierr); 969e5c89e4eSSatish Balay ierr = MPI_Type_commit(&MPIU_2SCALAR);CHKERRQ(ierr); 970e5c89e4eSSatish Balay 97140df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 972e5c89e4eSSatish Balay ierr = MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);CHKERRQ(ierr); 973e5c89e4eSSatish Balay ierr = MPI_Type_commit(&MPIU_2INT);CHKERRQ(ierr); 97444041f26SJed Brown #endif 975e5c89e4eSSatish Balay 976ec957eceSBarry Smith 977e5c89e4eSSatish Balay /* 978480cf27aSJed Brown Attributes to be set on PETSc communicators 979480cf27aSJed Brown */ 98012801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);CHKERRQ(ierr); 98112801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);CHKERRQ(ierr); 98212801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);CHKERRQ(ierr); 9835f7487a0SJunchao Zhang ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_DelComm_Shm,&Petsc_ShmComm_keyval,(void*)0);CHKERRQ(ierr); 984480cf27aSJed Brown 985480cf27aSJed Brown /* 986e8fb0fc0SBarry Smith Build the options database 987e5c89e4eSSatish Balay */ 988c5929fdfSBarry Smith ierr = PetscOptionsInsert(NULL,argc,args,file);CHKERRQ(ierr); 989e5c89e4eSSatish Balay 990703f3eceSBarry Smith /* call a second time so it can look in the options database */ 991703f3eceSBarry Smith ierr = PetscErrorPrintfInitialize();CHKERRQ(ierr); 9926dc8fec2Sbcordonn 993e5c89e4eSSatish Balay /* 994e5c89e4eSSatish Balay Print main application help message 995e5c89e4eSSatish Balay */ 9962d747510SLisandro Dalcin ierr = PetscOptionsHasHelp(NULL,&flg);CHKERRQ(ierr); 997e5c89e4eSSatish Balay if (help && flg) { 998e5c89e4eSSatish Balay ierr = PetscPrintf(PETSC_COMM_WORLD,help);CHKERRQ(ierr); 999e5c89e4eSSatish Balay } 1000e5c89e4eSSatish Balay ierr = PetscOptionsCheckInitial_Private();CHKERRQ(ierr); 1001e5c89e4eSSatish Balay 1002d45a07a7SBarry Smith ierr = PetscCitationsInitialize();CHKERRQ(ierr); 1003d45a07a7SBarry Smith 1004e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 100511525c0dSBarry Smith ierr = PetscInitializeSAWs(help);CHKERRQ(ierr); 1006f4202a44SBarry Smith #endif 1007f4202a44SBarry Smith 1008e5c89e4eSSatish Balay /* 1009e5c89e4eSSatish Balay Load the dynamic libraries (on machines that support them), this registers all 1010e5c89e4eSSatish Balay the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes) 1011e5c89e4eSSatish Balay */ 1012e5c89e4eSSatish Balay ierr = PetscInitialize_DynamicLibraries();CHKERRQ(ierr); 1013e5c89e4eSSatish Balay 1014e5c89e4eSSatish Balay ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 1015ae15b995SBarry Smith ierr = PetscInfo1(0,"PETSc successfully started: number of processors = %d\n",size);CHKERRQ(ierr); 1016e5c89e4eSSatish Balay ierr = PetscGetHostName(hostname,256);CHKERRQ(ierr); 1017ae15b995SBarry Smith ierr = PetscInfo1(0,"Running on machine: %s\n",hostname);CHKERRQ(ierr); 1018cd1b99a6SBarry Smith #if defined(PETSC_HAVE_OPENMP) 1019101491d0SBarry Smith { 1020101491d0SBarry Smith PetscBool omp_view_flag; 1021cd1b99a6SBarry Smith char *threads = getenv("OMP_NUM_THREADS"); 1022e5c89e4eSSatish Balay 1023cd1b99a6SBarry Smith if (threads) { 1024101491d0SBarry Smith ierr = PetscInfo1(0,"Number of OpenMP threads %s (given by OMP_NUM_THREADS)\n",threads);CHKERRQ(ierr); 1025f90b075cSBarry Smith (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads); 1026cd1b99a6SBarry Smith } else { 1027cd1b99a6SBarry Smith #define NMAX 10000 1028101491d0SBarry Smith int i; 1029cd1b99a6SBarry Smith PetscScalar *x; 1030cd1b99a6SBarry Smith ierr = PetscMalloc1(NMAX,&x);CHKERRQ(ierr); 1031cd1b99a6SBarry Smith #pragma omp parallel for 1032cd1b99a6SBarry Smith for (i=0; i<NMAX; i++) { 1033cd1b99a6SBarry Smith x[i] = 0.0; 1034f90b075cSBarry Smith PetscNumOMPThreads = (PetscInt) omp_get_num_threads(); 1035cd1b99a6SBarry Smith } 1036cd1b99a6SBarry Smith ierr = PetscFree(x);CHKERRQ(ierr); 1037f90b075cSBarry Smith ierr = PetscInfo1(0,"Number of OpenMP threads %D (number not set with OMP_NUM_THREADS, chosen by system)\n",PetscNumOMPThreads);CHKERRQ(ierr); 1038101491d0SBarry Smith } 1039101491d0SBarry Smith ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");CHKERRQ(ierr); 1040f90b075cSBarry 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); 1041101491d0SBarry Smith ierr = PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);CHKERRQ(ierr); 1042101491d0SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 1043101491d0SBarry Smith if (flg) { 1044f90b075cSBarry Smith ierr = PetscInfo1(0,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);CHKERRQ(ierr); 1045f90b075cSBarry Smith omp_set_num_threads((int)PetscNumOMPThreads); 1046101491d0SBarry Smith } 1047101491d0SBarry Smith if (omp_view_flag) { 1048f90b075cSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);CHKERRQ(ierr); 1049cd1b99a6SBarry Smith } 1050cd1b99a6SBarry Smith } 1051cd1b99a6SBarry Smith #endif 1052e5c89e4eSSatish Balay ierr = PetscOptionsCheckInitial_Components();CHKERRQ(ierr); 1053ef6c6fedSBoyana Norris /* Check the options database for options related to the options database itself */ 1054c5929fdfSBarry Smith ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr); 1055ef6c6fedSBoyana Norris 1056951e3c8eSBarry Smith #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32) 1057e39fd77fSBarry Smith /* 1058e39fd77fSBarry Smith Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI 1059e39fd77fSBarry Smith 1060e39fd77fSBarry Smith Currently not used because it is not supported by MPICH. 1061e39fd77fSBarry Smith */ 106230815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) { 10630298fd71SBarry Smith ierr = MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);CHKERRQ(ierr); 106430815ce0SLisandro Dalcin } 1065951e3c8eSBarry Smith #endif 1066e39fd77fSBarry Smith 106741c0b4b3SShri Abhyankar /* 106841c0b4b3SShri Abhyankar Setup building of stack frames for all function calls 106941c0b4b3SShri Abhyankar */ 1070ef19f930SBarry Smith #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1071e1167bb9SShri Abhyankar ierr = PetscStackCreate();CHKERRQ(ierr); 1072e1167bb9SShri Abhyankar #endif 1073e1167bb9SShri Abhyankar 10742d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 10752d53ad75SBarry Smith ierr = PetscFPTCreate(10000);CHKERRQ(ierr); 10762d53ad75SBarry Smith #endif 10772d53ad75SBarry Smith 10785e71baefSBarry Smith #if defined(PETSC_HAVE_HWLOC) 107987c3beb6SLisandro Dalcin { 108087c3beb6SLisandro Dalcin PetscViewer viewer; 108122e7e69cSBarry Smith ierr = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);CHKERRQ(ierr); 10825e71baefSBarry Smith if (flg) { 10835e71baefSBarry Smith ierr = PetscProcessPlacementView(viewer);CHKERRQ(ierr); 108487c3beb6SLisandro Dalcin ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 108587c3beb6SLisandro Dalcin } 10865e71baefSBarry Smith } 10875e71baefSBarry Smith #endif 1088dff31646SBarry Smith 108987c3beb6SLisandro Dalcin flg = PETSC_TRUE; 109087c3beb6SLisandro Dalcin ierr = PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);CHKERRQ(ierr); 109187c3beb6SLisandro Dalcin if (!flg) {ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE); CHKERRQ(ierr);} 109287c3beb6SLisandro Dalcin 1093a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 1094a56f64adSBarry Smith ierr = adios_init_noxml(PETSC_COMM_WORLD);CHKERRQ(ierr); 10957b56e58cSSatish Balay ierr = adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);CHKERRQ(ierr); 1096a56f64adSBarry Smith ierr = adios_select_method(Petsc_adios_group,"MPI","","");CHKERRQ(ierr); 10979fc7e16cSBarry Smith ierr = adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");CHKERRQ(ierr); 1098a56f64adSBarry Smith #endif 109955d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2) 110055d657eeSBarry Smith #endif 1101a56f64adSBarry Smith 1102301d30feSBarry Smith /* 1103301d30feSBarry Smith Once we are completedly initialized then we can set this variables 1104301d30feSBarry Smith */ 1105301d30feSBarry Smith PetscInitializeCalled = PETSC_TRUE; 11062db0e300SLisandro Dalcin 11072db0e300SLisandro Dalcin ierr = PetscOptionsHasName(NULL,NULL,"-python",&flg);CHKERRQ(ierr); 11082db0e300SLisandro Dalcin if (flg) {ierr = PetscPythonInitialize(NULL,NULL);CHKERRQ(ierr);} 1109301d30feSBarry Smith PetscFunctionReturn(0); 1110e5c89e4eSSatish Balay } 1111e5c89e4eSSatish Balay 11124097062eSBarry Smith #if defined(PETSC_USE_LOG) 111395c0884eSLisandro Dalcin PETSC_INTERN PetscObject *PetscObjects; 1114ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsCounts; 1115ecfe9a72SLisandro Dalcin PETSC_INTERN PetscInt PetscObjectsMaxCounts; 111695c0884eSLisandro Dalcin PETSC_INTERN PetscBool PetscObjectsLog; 11174097062eSBarry Smith #endif 1118e5c89e4eSSatish Balay 1119008a6e76SBarry Smith /* 1120008a6e76SBarry Smith Frees all the MPI types and operations that PETSc may have created 1121008a6e76SBarry Smith */ 1122008a6e76SBarry Smith PetscErrorCode PetscFreeMPIResources(void) 1123008a6e76SBarry Smith { 1124008a6e76SBarry Smith PetscErrorCode ierr; 1125008a6e76SBarry Smith 1126008a6e76SBarry Smith PetscFunctionBegin; 1127008a6e76SBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 1128008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU___FLOAT128);CHKERRQ(ierr); 1129008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1130008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU___COMPLEX128);CHKERRQ(ierr); 1131008a6e76SBarry Smith #endif 1132008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr); 1133008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr); 1134008a6e76SBarry Smith #elif defined(PETSC_USE_REAL___FP16) 1135008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU___FP16);CHKERRQ(ierr); 1136008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MAX);CHKERRQ(ierr); 1137008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MIN);CHKERRQ(ierr); 1138008a6e76SBarry Smith #endif 1139008a6e76SBarry Smith 1140008a6e76SBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1141008a6e76SBarry Smith #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX) 1142008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);CHKERRQ(ierr); 1143008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_C_COMPLEX);CHKERRQ(ierr); 1144008a6e76SBarry Smith #endif 1145008a6e76SBarry Smith #endif 1146008a6e76SBarry Smith 1147008a6e76SBarry Smith #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 1148008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_SUM);CHKERRQ(ierr); 1149008a6e76SBarry Smith #endif 1150008a6e76SBarry Smith 1151008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_2SCALAR);CHKERRQ(ierr); 115240df0d72SSatish Balay #if defined(PETSC_USE_64BIT_INDICES) 1153008a6e76SBarry Smith ierr = MPI_Type_free(&MPIU_2INT);CHKERRQ(ierr); 1154008a6e76SBarry Smith #endif 1155008a6e76SBarry Smith ierr = MPI_Op_free(&MPIU_MAXSUM_OP);CHKERRQ(ierr); 1156008a6e76SBarry Smith PetscFunctionReturn(0); 1157008a6e76SBarry Smith } 1158008a6e76SBarry Smith 1159e5c89e4eSSatish Balay /*@C 1160e5c89e4eSSatish Balay PetscFinalize - Checks for options to be called at the conclusion 1161e5c89e4eSSatish Balay of the program. MPI_Finalize() is called only if the user had not 1162e5c89e4eSSatish Balay called MPI_Init() before calling PetscInitialize(). 1163e5c89e4eSSatish Balay 1164e5c89e4eSSatish Balay Collective on PETSC_COMM_WORLD 1165e5c89e4eSSatish Balay 1166e5c89e4eSSatish Balay Options Database Keys: 116726a7e8d4SBarry Smith + -options_view - Calls PetscOptionsView() 1168e5c89e4eSSatish Balay . -options_left - Prints unused options that remain in the database 11697eb1d149SBarry 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 1170e5c89e4eSSatish Balay . -mpidump - Calls PetscMPIDump() 117192f119d6SBarry Smith . -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed 1172e5c89e4eSSatish Balay . -malloc_info - Prints total memory usage 117392f119d6SBarry Smith - -malloc_view <optional filename> - Prints list of all memory allocated and where 1174e5c89e4eSSatish Balay 1175e5c89e4eSSatish Balay Level: beginner 1176e5c89e4eSSatish Balay 1177e5c89e4eSSatish Balay Note: 1178e5c89e4eSSatish Balay See PetscInitialize() for more general runtime options. 1179e5c89e4eSSatish Balay 118088c29154SBarry Smith .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd() 1181e5c89e4eSSatish Balay @*/ 11827087cfbeSBarry Smith PetscErrorCode PetscFinalize(void) 1183e5c89e4eSSatish Balay { 1184e5c89e4eSSatish Balay PetscErrorCode ierr; 11854bb5149bSJed Brown PetscMPIInt rank; 1186a8d2bbe5SBarry Smith PetscInt nopt; 11872bf49c77SBarry Smith PetscBool flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE; 1188dff31646SBarry Smith PetscBool flg; 118910463e74SBarry Smith #if defined(PETSC_USE_LOG) 119010463e74SBarry Smith char mname[PETSC_MAX_PATH_LEN]; 119110463e74SBarry Smith #endif 1192e5c89e4eSSatish Balay 1193e5c89e4eSSatish Balay if (!PetscInitializeCalled) { 11944b09e917SBarry Smith printf("PetscInitialize() must be called before PetscFinalize()\n"); 11953cbbc5ffSBarry Smith return(PETSC_ERR_ARG_WRONGSTATE); 1196e5c89e4eSSatish Balay } 11973cbbc5ffSBarry Smith PetscFunctionBegin; 11980298fd71SBarry Smith ierr = PetscInfo(NULL,"PetscFinalize() called\n");CHKERRQ(ierr); 1199b022a5c1SBarry Smith 12001f817a21SBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 1201a56f64adSBarry Smith #if defined(PETSC_HAVE_ADIOS) 120222580e64SBarry Smith ierr = adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);CHKERRQ(ierr); 1203a56f64adSBarry Smith ierr = adios_finalize(rank);CHKERRQ(ierr); 1204a56f64adSBarry Smith #endif 120555d657eeSBarry Smith #if defined(PETSC_HAVE_ADIOS2) 120655d657eeSBarry Smith #endif 1207c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-citations",&flg);CHKERRQ(ierr); 1208dff31646SBarry Smith if (flg) { 12091f817a21SBarry Smith char *cits, filename[PETSC_MAX_PATH_LEN]; 12101f817a21SBarry Smith FILE *fd = PETSC_STDOUT; 12111f817a21SBarry Smith 1212c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 12131f817a21SBarry Smith if (filename[0]) { 12141f817a21SBarry Smith ierr = PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);CHKERRQ(ierr); 12151f817a21SBarry Smith } 1216dff31646SBarry Smith ierr = PetscSegBufferGet(PetscCitationsList,1,&cits);CHKERRQ(ierr); 1217dff31646SBarry Smith cits[0] = 0; 1218dff31646SBarry Smith ierr = PetscSegBufferExtractAlloc(PetscCitationsList,&cits);CHKERRQ(ierr); 12191f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");CHKERRQ(ierr); 12201f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr); 12211f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);CHKERRQ(ierr); 12221f817a21SBarry Smith ierr = PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");CHKERRQ(ierr); 12231f817a21SBarry Smith ierr = PetscFClose(PETSC_COMM_WORLD,fd);CHKERRQ(ierr); 1224dff31646SBarry Smith ierr = PetscFree(cits);CHKERRQ(ierr); 1225dff31646SBarry Smith } 1226dff31646SBarry Smith ierr = PetscSegBufferDestroy(&PetscCitationsList);CHKERRQ(ierr); 1227dff31646SBarry Smith 1228c2a97968SBarry Smith #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER) 122904102261SBarry Smith /* TextBelt is run for testing purposes only, please do not use this feature often */ 123004102261SBarry Smith { 123104102261SBarry Smith PetscInt nmax = 2; 123204102261SBarry Smith char **buffs; 123304102261SBarry Smith ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr); 1234c5929fdfSBarry Smith ierr = PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);CHKERRQ(ierr); 123504102261SBarry Smith if (flg1) { 123604102261SBarry Smith if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\""); 123704102261SBarry Smith if (nmax == 1) { 123804102261SBarry Smith ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr); 123904102261SBarry Smith ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr); 124004102261SBarry Smith ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr); 124104102261SBarry Smith } 124204102261SBarry Smith ierr = PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr); 124304102261SBarry Smith ierr = PetscFree(buffs[0]);CHKERRQ(ierr); 124404102261SBarry Smith ierr = PetscFree(buffs[1]);CHKERRQ(ierr); 124504102261SBarry Smith } 124604102261SBarry Smith ierr = PetscFree(buffs);CHKERRQ(ierr); 124704102261SBarry Smith } 1248763ec7b1SBarry Smith { 1249763ec7b1SBarry Smith PetscInt nmax = 2; 1250763ec7b1SBarry Smith char **buffs; 1251763ec7b1SBarry Smith ierr = PetscMalloc1(2,&buffs);CHKERRQ(ierr); 1252763ec7b1SBarry Smith ierr = PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);CHKERRQ(ierr); 1253763ec7b1SBarry Smith if (flg1) { 1254763ec7b1SBarry Smith if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\""); 1255763ec7b1SBarry Smith if (nmax == 1) { 1256763ec7b1SBarry Smith ierr = PetscMalloc1(128,&buffs[1]);CHKERRQ(ierr); 1257763ec7b1SBarry Smith ierr = PetscGetProgramName(buffs[1],32);CHKERRQ(ierr); 1258763ec7b1SBarry Smith ierr = PetscStrcat(buffs[1]," has completed");CHKERRQ(ierr); 1259763ec7b1SBarry Smith } 1260763ec7b1SBarry Smith ierr = PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);CHKERRQ(ierr); 1261763ec7b1SBarry Smith ierr = PetscFree(buffs[0]);CHKERRQ(ierr); 1262763ec7b1SBarry Smith ierr = PetscFree(buffs[1]);CHKERRQ(ierr); 1263763ec7b1SBarry Smith } 1264763ec7b1SBarry Smith ierr = PetscFree(buffs);CHKERRQ(ierr); 1265763ec7b1SBarry Smith } 126604102261SBarry Smith #endif 126767234432SDmitry Karpeev /* 126867234432SDmitry Karpeev It should be safe to cancel the options monitors, since we don't expect to be setting options 126967234432SDmitry Karpeev here (at least that are worth monitoring). Monitors ought to be released so that they release 127067234432SDmitry Karpeev whatever memory was allocated there before -malloc_dump reports unfreed memory. 127167234432SDmitry Karpeev */ 127267234432SDmitry Karpeev ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr); 127304102261SBarry Smith 12742d53ad75SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 12752d53ad75SBarry Smith ierr = PetscFPTDestroy();CHKERRQ(ierr); 12762d53ad75SBarry Smith #endif 12772d53ad75SBarry Smith 12782d53ad75SBarry Smith 1279e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1280dff31646SBarry Smith flg = PETSC_FALSE; 1281c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);CHKERRQ(ierr); 1282d5649816SBarry Smith if (flg) { 1283e04113cfSBarry Smith ierr = PetscOptionsSAWsDestroy();CHKERRQ(ierr); 1284d5649816SBarry Smith } 1285d5649816SBarry Smith #endif 1286d5649816SBarry Smith 1287681455b2SBarry Smith #if defined(PETSC_HAVE_X) 1288681455b2SBarry Smith flg1 = PETSC_FALSE; 1289c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);CHKERRQ(ierr); 1290681455b2SBarry Smith if (flg1) { 1291681455b2SBarry Smith /* this is a crude hack, but better than nothing */ 1292681455b2SBarry Smith ierr = PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);CHKERRQ(ierr); 1293681455b2SBarry Smith } 1294681455b2SBarry Smith #endif 1295681455b2SBarry Smith 129667584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 1297c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);CHKERRQ(ierr); 1298e5c89e4eSSatish Balay if (!flg2) { 129990d69ab7SBarry Smith flg2 = PETSC_FALSE; 1300c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);CHKERRQ(ierr); 1301e5c89e4eSSatish Balay } 1302e5c89e4eSSatish Balay if (flg2) { 13030841954dSBarry Smith ierr = PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");CHKERRQ(ierr); 1304e5c89e4eSSatish Balay } 130567584ceeSBarry Smith #endif 1306e5c89e4eSSatish Balay 1307e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 130890d69ab7SBarry Smith flg1 = PETSC_FALSE; 1309c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);CHKERRQ(ierr); 1310e5c89e4eSSatish Balay if (flg1) { 1311e5c89e4eSSatish Balay PetscLogDouble flops = 0; 1312205a32c2SJed Brown ierr = MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr); 1313e5c89e4eSSatish Balay ierr = PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);CHKERRQ(ierr); 1314e5c89e4eSSatish Balay } 1315e5c89e4eSSatish Balay #endif 1316e5c89e4eSSatish Balay 1317e5c89e4eSSatish Balay 1318e5c89e4eSSatish Balay #if defined(PETSC_USE_LOG) 1319e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MPE) 1320e5c89e4eSSatish Balay mname[0] = 0; 1321c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr); 1322e5c89e4eSSatish Balay if (flg1) { 1323e5c89e4eSSatish Balay if (mname[0]) {ierr = PetscLogMPEDump(mname);CHKERRQ(ierr);} 1324e5c89e4eSSatish Balay else {ierr = PetscLogMPEDump(0);CHKERRQ(ierr);} 1325e5c89e4eSSatish Balay } 1326e5c89e4eSSatish Balay #endif 1327356e5837SBarry Smith #endif 1328a297a907SKarl Rupp 1329dd710f27SBarry Smith /* 1330dd710f27SBarry Smith Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_(). 1331dd710f27SBarry Smith */ 1332dd710f27SBarry Smith ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr); 1333dd710f27SBarry Smith 1334356e5837SBarry Smith #if defined(PETSC_USE_LOG) 133587c3beb6SLisandro Dalcin ierr = PetscOptionsPushGetViewerOff(PETSC_FALSE);CHKERRQ(ierr); 1336f14045dbSBarry Smith ierr = PetscLogViewFromOptions();CHKERRQ(ierr); 133787c3beb6SLisandro Dalcin ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr); 133887c3beb6SLisandro Dalcin 1339356e5837SBarry Smith mname[0] = 0; 1340c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr); 1341e5c89e4eSSatish Balay if (flg1) { 134291eabc43SBarry Smith PetscViewer viewer; 134320a8bfc3SBarry Smith ierr = (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING: -log_summary is being deprecated; switch to -log_view\n\n\n");CHKERRQ(ierr); 134491eabc43SBarry Smith if (mname[0]) { 134591eabc43SBarry Smith ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);CHKERRQ(ierr); 134691eabc43SBarry Smith ierr = PetscLogView(viewer);CHKERRQ(ierr); 13476bf464f9SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 134833f85c2fSBarry Smith } else { 134933f85c2fSBarry Smith viewer = PETSC_VIEWER_STDOUT_WORLD; 13509a9a5d4cSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 135133f85c2fSBarry Smith ierr = PetscLogView(viewer);CHKERRQ(ierr); 13529a9a5d4cSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 135333f85c2fSBarry Smith } 1354e5c89e4eSSatish Balay } 1355a297a907SKarl Rupp 1356dd710f27SBarry Smith /* 1357dd710f27SBarry Smith Free any objects created by the last block of code. 1358dd710f27SBarry Smith */ 1359dd710f27SBarry Smith ierr = PetscObjectRegisterDestroyAll();CHKERRQ(ierr); 1360dd710f27SBarry Smith 1361dd710f27SBarry Smith mname[0] = 0; 1362c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);CHKERRQ(ierr); 1363c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);CHKERRQ(ierr); 13647ff663adSLisandro Dalcin if (flg1 || flg2) {ierr = PetscLogDump(mname);CHKERRQ(ierr);} 1365e5c89e4eSSatish Balay #endif 136610463e74SBarry Smith 136733f85c2fSBarry Smith ierr = PetscStackDestroy();CHKERRQ(ierr); 136810463e74SBarry Smith 136990d69ab7SBarry Smith flg1 = PETSC_FALSE; 1370c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);CHKERRQ(ierr); 1371e5c89e4eSSatish Balay if (!flg1) { ierr = PetscPopSignalHandler();CHKERRQ(ierr);} 137290d69ab7SBarry Smith flg1 = PETSC_FALSE; 1373c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);CHKERRQ(ierr); 1374e5c89e4eSSatish Balay if (flg1) { 1375e5c89e4eSSatish Balay ierr = PetscMPIDump(stdout);CHKERRQ(ierr); 1376e5c89e4eSSatish Balay } 137790d69ab7SBarry Smith flg1 = PETSC_FALSE; 137890d69ab7SBarry Smith flg2 = PETSC_FALSE; 13798bb29257SSatish Balay /* preemptive call to avoid listing this option in options table as unused */ 1380c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);CHKERRQ(ierr); 1381c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr); 1382c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);CHKERRQ(ierr); 1383e4c476e2SSatish Balay 1384e5c89e4eSSatish Balay if (flg2) { 1385be56827dSJed Brown PetscViewer viewer; 138602ba9f54SBarry Smith ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr); 138702ba9f54SBarry Smith ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr); 1388c5929fdfSBarry Smith ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr); 1389be56827dSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1390e5c89e4eSSatish Balay } 1391e5c89e4eSSatish Balay 1392e5c89e4eSSatish Balay /* to prevent PETSc -options_left from warning */ 1393c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-nox",&flg1);CHKERRQ(ierr); 1394c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);CHKERRQ(ierr); 1395e5c89e4eSSatish Balay 139633fc4174SSatish Balay flg3 = PETSC_FALSE; /* default value is required */ 1397c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);CHKERRQ(ierr); 13983de2bfdfSBarry Smith #if defined(PETSC_USE_DEBUG) 13993de2bfdfSBarry Smith if (!flg1) flg3 = PETSC_TRUE; 14003de2bfdfSBarry Smith #endif 1401e5c89e4eSSatish Balay if (flg3) { 14023de2bfdfSBarry Smith if (!flg2 && flg1) { /* have not yet printed the options */ 1403be56827dSJed Brown PetscViewer viewer; 140402ba9f54SBarry Smith ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr); 140502ba9f54SBarry Smith ierr = PetscViewerSetType(viewer,PETSCVIEWERASCII);CHKERRQ(ierr); 1406c5929fdfSBarry Smith ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr); 1407be56827dSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1408e5c89e4eSSatish Balay } 14093de2bfdfSBarry Smith ierr = PetscOptionsAllUsed(NULL,&nopt);CHKERRQ(ierr); 14103de2bfdfSBarry Smith if (nopt) { 14113de2bfdfSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");CHKERRQ(ierr); 14123de2bfdfSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");CHKERRQ(ierr); 14133de2bfdfSBarry Smith if (nopt == 1) { 1414e5c89e4eSSatish Balay ierr = PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");CHKERRQ(ierr); 1415e5c89e4eSSatish Balay } else { 14167582186dSLisandro Dalcin ierr = PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);CHKERRQ(ierr); 1417e5c89e4eSSatish Balay } 14183de2bfdfSBarry Smith } else if (flg3 && flg1) { 14193de2bfdfSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");CHKERRQ(ierr); 1420df12ba86SBarry Smith } 1421c5929fdfSBarry Smith ierr = PetscOptionsLeft(NULL);CHKERRQ(ierr); 1422e5c89e4eSSatish Balay } 1423e5c89e4eSSatish Balay 1424e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1425d45a07a7SBarry Smith if (!PetscGlobalRank) { 142687f587eeSBarry Smith ierr = PetscStackSAWsViewOff();CHKERRQ(ierr); 142716ad0300SBarry Smith PetscStackCallSAWs(SAWs_Finalize,()); 1428d45a07a7SBarry Smith } 1429ec957eceSBarry Smith #endif 1430ec957eceSBarry Smith 14314097062eSBarry Smith #if defined(PETSC_USE_LOG) 143210463e74SBarry Smith /* 1433dbc8283eSBarry Smith List all objects the user may have forgot to free 14342eff7a51SBarry Smith */ 143505df10baSBarry Smith if (PetscObjectsLog) { 1436c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr); 1437a64a8e02SBarry Smith if (flg1) { 1438a64a8e02SBarry Smith MPI_Comm local_comm; 14397eb1d149SBarry Smith char string[64]; 1440a64a8e02SBarry Smith 1441c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr); 1442a64a8e02SBarry Smith ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr); 1443a64a8e02SBarry Smith ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr); 14447eb1d149SBarry Smith ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 1445a64a8e02SBarry Smith ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr); 1446a64a8e02SBarry Smith ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr); 14470a1571b3SBarry Smith } 144805df10baSBarry Smith } 14494097062eSBarry Smith #endif 14504097062eSBarry Smith 14514097062eSBarry Smith #if defined(PETSC_USE_LOG) 1452dbc8283eSBarry Smith PetscObjectsCounts = 0; 1453dbc8283eSBarry Smith PetscObjectsMaxCounts = 0; 1454a297a907SKarl Rupp ierr = PetscFree(PetscObjects);CHKERRQ(ierr); 14554097062eSBarry Smith #endif 14562eff7a51SBarry Smith 145733f85c2fSBarry Smith /* 145833f85c2fSBarry Smith Destroy any packages that registered a finalize 145933f85c2fSBarry Smith */ 146033f85c2fSBarry Smith ierr = PetscRegisterFinalizeAll();CHKERRQ(ierr); 146133f85c2fSBarry Smith 1462101409b8SToby Isaac #if defined(PETSC_USE_LOG) 1463fa2bb9feSLisandro Dalcin ierr = PetscLogFinalize();CHKERRQ(ierr); 1464101409b8SToby Isaac #endif 1465101409b8SToby Isaac 146633f85c2fSBarry Smith /* 146748dd1dffSBarry Smith Print PetscFunctionLists that have not been properly freed 146848dd1dffSBarry Smith 146937e93019SBarry Smith ierr = PetscFunctionListPrintAll();CHKERRQ(ierr); 147048dd1dffSBarry Smith */ 147137e93019SBarry Smith 14724028d114SSatish Balay if (petsc_history) { 1473f3dea69dSBarry Smith ierr = PetscCloseHistoryFile(&petsc_history);CHKERRQ(ierr); 1474e5c89e4eSSatish Balay petsc_history = 0; 1475e5c89e4eSSatish Balay } 14769de0f6ecSBarry Smith ierr = PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr); 1477e5c89e4eSSatish Balay 14780298fd71SBarry Smith ierr = PetscInfoAllow(PETSC_FALSE,NULL);CHKERRQ(ierr); 1479e5c89e4eSSatish Balay 148067584ceeSBarry Smith #if !defined(PETSC_HAVE_THREADSAFETY) 148192f119d6SBarry Smith if (!(PETSC_RUNNING_ON_VALGRIND)) { 1482e5c89e4eSSatish Balay char fname[PETSC_MAX_PATH_LEN]; 148392f119d6SBarry Smith char sname[PETSC_MAX_PATH_LEN]; 1484e5c89e4eSSatish Balay FILE *fd; 1485ed9cf6e9SBarry Smith int err; 1486e5c89e4eSSatish Balay 1487dc92acbaSJed Brown flg2 = PETSC_FALSE; 148892f119d6SBarry Smith flg3 = PETSC_FALSE; 14898bf1f09cSShri Abhyankar #if defined(PETSC_USE_DEBUG) 149092f119d6SBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);CHKERRQ(ierr); 1491dc92acbaSJed Brown #endif 149292f119d6SBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);CHKERRQ(ierr); 149392f119d6SBarry Smith fname[0] = 0; 149492f119d6SBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,250,&flg1);CHKERRQ(ierr); 1495e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1496e5c89e4eSSatish Balay 14972e924ca5SSatish Balay PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank); 1498e32f2f54SBarry Smith fd = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname); 1499e5c89e4eSSatish Balay ierr = PetscMallocDump(fd);CHKERRQ(ierr); 1500ed9cf6e9SBarry Smith err = fclose(fd); 1501e32f2f54SBarry Smith if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 150292f119d6SBarry Smith } else if (flg1 || flg2 || flg3) { 1503e5c89e4eSSatish Balay MPI_Comm local_comm; 1504e5c89e4eSSatish Balay 1505e5c89e4eSSatish Balay ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr); 1506e5c89e4eSSatish Balay ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr); 1507e5c89e4eSSatish Balay ierr = PetscMallocDump(stdout);CHKERRQ(ierr); 1508e5c89e4eSSatish Balay ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr); 1509e5c89e4eSSatish Balay ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr); 1510e5c89e4eSSatish Balay } 1511e5c89e4eSSatish Balay fname[0] = 0; 151292f119d6SBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,250,&flg1);CHKERRQ(ierr); 1513e5c89e4eSSatish Balay if (flg1 && fname[0]) { 1514e5c89e4eSSatish Balay 151592f119d6SBarry Smith PetscSNPrintf(sname,PETSC_MAX_PATH_LEN,"%s_%d",fname,rank); 151692f119d6SBarry Smith fd = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname); 151792f119d6SBarry Smith ierr = PetscMallocView(fd);CHKERRQ(ierr); 1518ed9cf6e9SBarry Smith err = fclose(fd); 1519e32f2f54SBarry Smith if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 152092f119d6SBarry Smith } else if (flg1) { 152192f119d6SBarry Smith MPI_Comm local_comm; 152292f119d6SBarry Smith 152392f119d6SBarry Smith ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr); 152492f119d6SBarry Smith ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr); 152592f119d6SBarry Smith ierr = PetscMallocView(stdout);CHKERRQ(ierr); 152692f119d6SBarry Smith ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr); 152792f119d6SBarry Smith ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr); 1528e5c89e4eSSatish Balay } 1529e5c89e4eSSatish Balay } 153067584ceeSBarry Smith #endif 153120e2c332SMatthew G. Knepley 15325486ca60SMatthew G. Knepley /* 15335486ca60SMatthew G. Knepley Close any open dynamic libraries 15345486ca60SMatthew G. Knepley */ 15355486ca60SMatthew G. Knepley ierr = PetscFinalize_DynamicLibraries();CHKERRQ(ierr); 15365486ca60SMatthew G. Knepley 1537e5c89e4eSSatish Balay /* Can be destroyed only after all the options are used */ 15384416b707SBarry Smith ierr = PetscOptionsDestroyDefault();CHKERRQ(ierr); 1539e5c89e4eSSatish Balay 1540e5c89e4eSSatish Balay PetscGlobalArgc = 0; 1541e5c89e4eSSatish Balay PetscGlobalArgs = 0; 1542e5c89e4eSSatish Balay 1543008a6e76SBarry Smith ierr = PetscFreeMPIResources();CHKERRQ(ierr); 1544e5c89e4eSSatish Balay 1545dbc8283eSBarry Smith /* 1546efb80d3cSBarry Smith Destroy any known inner MPI_Comm's and attributes pointing to them 1547efb80d3cSBarry Smith Note this will not destroy any new communicators the user has created. 1548efb80d3cSBarry Smith 1549efb80d3cSBarry Smith If all PETSc objects were not destroyed those left over objects will have hanging references to 1550efb80d3cSBarry Smith the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again 1551dbc8283eSBarry Smith */ 1552b770b1f6SSatish Balay { 1553dbc8283eSBarry Smith PetscCommCounter *counter; 1554dbc8283eSBarry Smith PetscMPIInt flg; 1555dbc8283eSBarry Smith MPI_Comm icomm; 1556265f3f35SJed Brown union {MPI_Comm comm; void *ptr;} ucomm; 155747435625SJed Brown ierr = MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr); 1558dbc8283eSBarry Smith if (flg) { 1559265f3f35SJed Brown icomm = ucomm.comm; 156047435625SJed Brown ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr); 1561dbc8283eSBarry 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"); 1562dbc8283eSBarry Smith 156347435625SJed Brown ierr = MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);CHKERRQ(ierr); 156447435625SJed Brown ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr); 1565efb80d3cSBarry Smith ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr); 1566dbc8283eSBarry Smith } 156747435625SJed Brown ierr = MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);CHKERRQ(ierr); 1568dbc8283eSBarry Smith if (flg) { 1569265f3f35SJed Brown icomm = ucomm.comm; 157047435625SJed Brown ierr = MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr); 1571dbc8283eSBarry 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"); 1572dbc8283eSBarry Smith 157347435625SJed Brown ierr = MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);CHKERRQ(ierr); 157447435625SJed Brown ierr = MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);CHKERRQ(ierr); 1575efb80d3cSBarry Smith ierr = MPI_Comm_free(&icomm);CHKERRQ(ierr); 1576dbc8283eSBarry Smith } 1577b770b1f6SSatish Balay } 1578dbc8283eSBarry Smith 157947435625SJed Brown ierr = MPI_Comm_free_keyval(&Petsc_Counter_keyval);CHKERRQ(ierr); 158047435625SJed Brown ierr = MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);CHKERRQ(ierr); 158147435625SJed Brown ierr = MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);CHKERRQ(ierr); 15825f7487a0SJunchao Zhang ierr = MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);CHKERRQ(ierr); 1583480cf27aSJed Brown 15845ad9ad5bSBarry Smith ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);CHKERRQ(ierr); 15855ad9ad5bSBarry Smith ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);CHKERRQ(ierr); 15865ad9ad5bSBarry Smith ierr = PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);CHKERRQ(ierr); 1587ef19f930SBarry Smith ierr = PetscSpinlockDestroy(&PetscCommSpinLock);CHKERRQ(ierr); 1588ef19f930SBarry Smith 1589e5c89e4eSSatish Balay if (PetscBeganMPI) { 159099608316SBarry Smith #if defined(PETSC_HAVE_MPI_FINALIZED) 159199b1327fSBarry Smith PetscMPIInt flag; 159299b1327fSBarry Smith ierr = MPI_Finalized(&flag);CHKERRQ(ierr); 1593e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()"); 159499608316SBarry Smith #endif 1595e5c89e4eSSatish Balay ierr = MPI_Finalize();CHKERRQ(ierr); 1596e5c89e4eSSatish Balay } 1597e5c89e4eSSatish Balay /* 1598e5c89e4eSSatish Balay 1599e5c89e4eSSatish Balay Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because 1600e5c89e4eSSatish Balay the communicator has some outstanding requests on it. Specifically if the 1601e5c89e4eSSatish Balay flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See 1602e5c89e4eSSatish Balay src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate() 1603e5c89e4eSSatish Balay is never freed as it should be. Thus one may obtain messages of the form 16040e5e90baSSatish Balay [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the 1605e5c89e4eSSatish Balay memory was not freed. 1606e5c89e4eSSatish Balay 1607e5c89e4eSSatish Balay */ 16081d1a0024SBarry Smith ierr = PetscMallocClear();CHKERRQ(ierr); 1609a297a907SKarl Rupp 1610e5c89e4eSSatish Balay PetscInitializeCalled = PETSC_FALSE; 1611e5c89e4eSSatish Balay PetscFinalizeCalled = PETSC_TRUE; 16123db9a53dSBarry Smith PetscFunctionReturn(0); 1613e5c89e4eSSatish Balay } 1614e5c89e4eSSatish Balay 161543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame_) 16168cc058d9SJed Brown PETSC_EXTERN int lsame_(char *a,char *b) 161743db4dbbSBarry Smith { 161843db4dbbSBarry Smith if (*a == *b) return 1; 161943db4dbbSBarry Smith if (*a + 32 == *b) return 1; 162043db4dbbSBarry Smith if (*a - 32 == *b) return 1; 162143db4dbbSBarry Smith return 0; 162243db4dbbSBarry Smith } 1623a70650f6SBarry Smith #endif 162443db4dbbSBarry Smith 162543db4dbbSBarry Smith #if defined(PETSC_MISSING_LAPACK_lsame) 16268cc058d9SJed Brown PETSC_EXTERN int lsame(char *a,char *b) 162743db4dbbSBarry Smith { 162843db4dbbSBarry Smith if (*a == *b) return 1; 162943db4dbbSBarry Smith if (*a + 32 == *b) return 1; 163043db4dbbSBarry Smith if (*a - 32 == *b) return 1; 163143db4dbbSBarry Smith return 0; 163243db4dbbSBarry Smith } 163343db4dbbSBarry Smith #endif 1634