Actual source code: pinit.c

  1: #define PETSC_DESIRE_FEATURE_TEST_MACROS
  2: /*
  3:    This file defines the initialization of PETSc, including PetscInitialize()
  4: */
  5: #include <petsc/private/petscimpl.h>
  6: #include <petscviewer.h>

  8: #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
  9: #include <petsc/private/valgrind/valgrind.h>
 10: #endif

 12: #if defined(PETSC_HAVE_FORTRAN)
 13: #include <petsc/private/fortranimpl.h>
 14: #endif

 16: #if defined(PETSC_HAVE_CUDA)
 17: #include <petsc/private/deviceimpl.h>
 18: PETSC_EXTERN cudaEvent_t petsc_gputimer_begin;
 19: PETSC_EXTERN cudaEvent_t petsc_gputimer_end;
 20: #endif

 22: #if defined(PETSC_HAVE_HIP)
 23: #include <petsc/private/deviceimpl.h>
 24: PETSC_EXTERN hipEvent_t petsc_gputimer_begin;
 25: PETSC_EXTERN hipEvent_t petsc_gputimer_end;
 26: #endif

 28: #if defined(PETSC_USE_GCOV)
 29: EXTERN_C_BEGIN
 30: void  __gcov_flush(void);
 31: EXTERN_C_END
 32: #endif

 34: #if PetscDefined(USE_LOG)
 35: PETSC_INTERN PetscErrorCode PetscLogFinalize(void);
 36: #endif

 38: #if defined(PETSC_SERIALIZE_FUNCTIONS)
 39: PETSC_INTERN PetscFPT PetscFPTData;
 40: PetscFPT PetscFPTData = 0;
 41: #endif

 43: #if PetscDefined(HAVE_SAWS)
 44: #include <petscviewersaws.h>
 45: #endif

 47: /* -----------------------------------------------------------------------------------------*/

 49: PETSC_INTERN FILE *petsc_history;

 51: PETSC_INTERN PetscErrorCode PetscInitialize_DynamicLibraries(void);
 52: PETSC_INTERN PetscErrorCode PetscFinalize_DynamicLibraries(void);
 53: PETSC_INTERN PetscErrorCode PetscFunctionListPrintAll(void);
 54: PETSC_INTERN PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
 55: PETSC_INTERN PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
 56: PETSC_INTERN PetscErrorCode PetscCloseHistoryFile(FILE**);

 58: /* user may set these BEFORE calling PetscInitialize() */
 59: MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;
 60: #if PetscDefined(HAVE_MPI_INIT_THREAD)
 61: PetscMPIInt PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_FUNNELED;
 62: #else
 63: PetscMPIInt PETSC_MPI_THREAD_REQUIRED = 0;
 64: #endif

 66: PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
 67: PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
 68: PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;
 69: PetscMPIInt Petsc_ShmComm_keyval   = MPI_KEYVAL_INVALID;

 71: /*
 72:      Declare and set all the string names of the PETSc enums
 73: */
 74: const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",NULL};
 75: const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",NULL};

 77: PetscBool PetscPreLoadingUsed = PETSC_FALSE;
 78: PetscBool PetscPreLoadingOn   = PETSC_FALSE;

 80: PetscInt PetscHotRegionDepth;

 82: PetscBool PETSC_RUNNING_ON_VALGRIND = PETSC_FALSE;

 84: #if defined(PETSC_HAVE_THREADSAFETY)
 85: PetscSpinlock PetscViewerASCIISpinLockOpen;
 86: PetscSpinlock PetscViewerASCIISpinLockStdout;
 87: PetscSpinlock PetscViewerASCIISpinLockStderr;
 88: PetscSpinlock PetscCommSpinLock;
 89: #endif

 91: /*
 92:       PetscInitializeNoPointers - Calls PetscInitialize() from C/C++ without the pointers to argc and args

 94:    Collective

 96:    Level: advanced

 98:     Notes:
 99:     this is called only by the PETSc Julia interface. Even though it might start MPI it sets the flag to
100:      indicate that it did NOT start MPI so that the PetscFinalize() does not end MPI, thus allowing PetscInitialize() to
101:      be called multiple times from Julia without the problem of trying to initialize MPI more than once.

103:      Developer Note: Turns off PETSc signal handling to allow Julia to manage signals

105: .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
106: */
107: PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
108: {
110:   int            myargc   = argc;
111:   char           **myargs = args;

114:   PetscInitialize(&myargc,&myargs,filename,help);if (ierr) PetscFunctionReturn(ierr);
115:   PetscPopSignalHandler();
116:   PetscBeganMPI = PETSC_FALSE;
117:   return(0);
118: }

120: /*
121:       Used by Julia interface to get communicator
122: */
123: PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
124: {
127:   *comm = PETSC_COMM_SELF;
128:   return(0);
129: }

131: /*@C
132:       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
133:         the command line arguments.

135:    Collective

137:    Level: advanced

139: .seealso: PetscInitialize(), PetscInitializeFortran()
140: @*/
141: PetscErrorCode  PetscInitializeNoArguments(void)
142: {
144:   int            argc   = 0;
145:   char           **args = NULL;

148:   PetscInitialize(&argc,&args,NULL,NULL);
149:   PetscFunctionReturn(ierr);
150: }

152: /*@
153:       PetscInitialized - Determine whether PETSc is initialized.

155:    Level: beginner

157: .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
158: @*/
159: PetscErrorCode PetscInitialized(PetscBool *isInitialized)
160: {
163:   *isInitialized = PetscInitializeCalled;
164:   return(0);
165: }

167: /*@
168:       PetscFinalized - Determine whether PetscFinalize() has been called yet

170:    Level: developer

172: .seealso: PetscInitialize(), PetscInitializeNoArguments(), PetscInitializeFortran()
173: @*/
174: PetscErrorCode  PetscFinalized(PetscBool  *isFinalized)
175: {
178:   *isFinalized = PetscFinalizeCalled;
179:   return(0);
180: }

182: PETSC_INTERN PetscErrorCode PetscOptionsCheckInitial_Private(const char []);

184: /*
185:        This function is the MPI reduction operation used to compute the sum of the
186:    first half of the datatype and the max of the second half.
187: */
188: MPI_Op MPIU_MAXSUM_OP = 0;

190: PETSC_INTERN void MPIAPI MPIU_MaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
191: {
192:   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;

195:   if (*datatype != MPIU_2INT) {
196:     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
197:     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
198:   }

200:   for (i=0; i<count; i++) {
201:     xout[2*i]    = PetscMax(xout[2*i],xin[2*i]);
202:     xout[2*i+1] += xin[2*i+1];
203:   }
204:   PetscFunctionReturnVoid();
205: }

207: /*
208:     Returns the max of the first entry owned by this processor and the
209: sum of the second entry.

211:     The reason sizes[2*i] contains lengths sizes[2*i+1] contains flag of 1 if length is nonzero
212: is so that the MPIU_MAXSUM_OP() can set TWO values, if we passed in only sizes[i] with lengths
213: there would be no place to store the both needed results.
214: */
215: PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
216: {

220: #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
221:   {
222:     struct {PetscInt max,sum;} work;
223:     MPI_Reduce_scatter_block((void*)sizes,&work,1,MPIU_2INT,MPIU_MAXSUM_OP,comm);
224:     *max = work.max;
225:     *sum = work.sum;
226:   }
227: #else
228:   {
229:     PetscMPIInt    size,rank;
230:     struct {PetscInt max,sum;} *work;
231:     MPI_Comm_size(comm,&size);
232:     MPI_Comm_rank(comm,&rank);
233:     PetscMalloc1(size,&work);
234:     MPIU_Allreduce((void*)sizes,work,size,MPIU_2INT,MPIU_MAXSUM_OP,comm);
235:     *max = work[rank].max;
236:     *sum = work[rank].sum;
237:     PetscFree(work);
238:   }
239: #endif
240:   return(0);
241: }

243: /* ----------------------------------------------------------------------------*/

245: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
246: MPI_Op MPIU_SUM = 0;

248: PETSC_EXTERN void MPIAPI PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
249: {
250:   PetscInt i,count = *cnt;

253:   if (*datatype == MPIU_REAL) {
254:     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
255:     for (i=0; i<count; i++) xout[i] += xin[i];
256:   }
257: #if defined(PETSC_HAVE_COMPLEX)
258:   else if (*datatype == MPIU_COMPLEX) {
259:     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
260:     for (i=0; i<count; i++) xout[i] += xin[i];
261:   }
262: #endif
263:   else {
264:     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
265:     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
266:   }
267:   PetscFunctionReturnVoid();
268: }
269: #endif

271: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
272: MPI_Op MPIU_MAX = 0;
273: MPI_Op MPIU_MIN = 0;

275: PETSC_EXTERN void MPIAPI PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
276: {
277:   PetscInt i,count = *cnt;

280:   if (*datatype == MPIU_REAL) {
281:     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
282:     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
283:   }
284: #if defined(PETSC_HAVE_COMPLEX)
285:   else if (*datatype == MPIU_COMPLEX) {
286:     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
287:     for (i=0; i<count; i++) {
288:       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
289:     }
290:   }
291: #endif
292:   else {
293:     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
294:     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
295:   }
296:   PetscFunctionReturnVoid();
297: }

299: PETSC_EXTERN void MPIAPI PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
300: {
301:   PetscInt    i,count = *cnt;

304:   if (*datatype == MPIU_REAL) {
305:     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
306:     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
307:   }
308: #if defined(PETSC_HAVE_COMPLEX)
309:   else if (*datatype == MPIU_COMPLEX) {
310:     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
311:     for (i=0; i<count; i++) {
312:       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
313:     }
314:   }
315: #endif
316:   else {
317:     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
318:     PETSCABORT(MPI_COMM_SELF,PETSC_ERR_ARG_WRONG);
319:   }
320:   PetscFunctionReturnVoid();
321: }
322: #endif

324: /*
325:    Private routine to delete internal tag/name counter storage when a communicator is freed.

327:    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.

329:    Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()

331: */
332: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_Counter_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
333: {
334:   PetscErrorCode   ierr;
335:   PetscCommCounter *counter=(PetscCommCounter*)count_val;

338:   PetscInfo1(NULL,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);
339:   PetscFree(counter->iflags);
340:   PetscFree(counter);
341:   PetscFunctionReturn(MPI_SUCCESS);
342: }

344: /*
345:   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Comm_delete_attr) or when the user
346:   calls MPI_Comm_free().

348:   This is the only entry point for breaking the links between inner and outer comms.

350:   This is called by MPI, not by users. This is called when MPI_Comm_free() is called on the communicator.

352:   Note: this is declared extern "C" because it is passed to MPI_Comm_create_keyval()

354: */
355: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_InnerComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
356: {
357:   PetscErrorCode                    ierr;
358:   union {MPI_Comm comm; void *ptr;} icomm;

361:   if (keyval != Petsc_InnerComm_keyval) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
362:   icomm.ptr = attr_val;
363:   if (PetscDefined(USE_DEBUG)) {
364:     /* Error out if the inner/outer comms are not correctly linked through their Outer/InnterComm attributes */
365:     PetscMPIInt flg;
366:     union {MPI_Comm comm; void *ptr;} ocomm;
367:     MPI_Comm_get_attr(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);
368:     if (!flg) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm does not have OuterComm attribute");
369:     if (ocomm.comm != comm) SETERRMPI(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner comm's OuterComm attribute does not point to outer PETSc comm");
370:   }
371:   MPI_Comm_delete_attr(icomm.comm,Petsc_OuterComm_keyval);
372:   PetscInfo2(NULL,"User MPI_Comm %ld is being unlinked from inner PETSc comm %ld\n",(long)comm,(long)icomm.comm);
373:   PetscFunctionReturn(MPI_SUCCESS);
374: }

376: /*
377:  * This is invoked on the inner comm when Petsc_InnerComm_Attr_Delete_Fn calls MPI_Comm_delete_attr().  It should not be reached any other way.
378:  */
379: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_OuterComm_Attr_Delete_Fn(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
380: {

384:   PetscInfo1(NULL,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);
385:   PetscFunctionReturn(MPI_SUCCESS);
386: }

388: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_ShmComm_Attr_Delete_Fn(MPI_Comm,PetscMPIInt,void *,void *);

390: #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
391: PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
392: PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
393: PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
394: #endif

396: PetscMPIInt PETSC_MPI_ERROR_CLASS=MPI_ERR_LASTCODE,PETSC_MPI_ERROR_CODE;

398: PETSC_INTERN int  PetscGlobalArgc;
399: PETSC_INTERN char **PetscGlobalArgs;
400: int  PetscGlobalArgc   = 0;
401: char **PetscGlobalArgs = NULL;
402: PetscSegBuffer PetscCitationsList;

404: PetscErrorCode PetscCitationsInitialize(void)
405: {

409:   PetscSegBufferCreate(1,10000,&PetscCitationsList);
410:   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);
411:   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);
412:   return(0);
413: }

415: static char programname[PETSC_MAX_PATH_LEN] = ""; /* HP includes entire path in name */

417: PetscErrorCode  PetscSetProgramName(const char name[])
418: {

422:   PetscStrncpy(programname,name,sizeof(programname));
423:   return(0);
424: }

426: /*@C
427:     PetscGetProgramName - Gets the name of the running program.

429:     Not Collective

431:     Input Parameter:
432: .   len - length of the string name

434:     Output Parameter:
435: .   name - the name of the running program

437:    Level: advanced

439:     Notes:
440:     The name of the program is copied into the user-provided character
441:     array of length len.  On some machines the program name includes
442:     its entire path, so one should generally set len >= PETSC_MAX_PATH_LEN.
443: @*/
444: PetscErrorCode  PetscGetProgramName(char name[],size_t len)
445: {

449:   PetscStrncpy(name,programname,len);
450:   return(0);
451: }

453: /*@C
454:    PetscGetArgs - Allows you to access the raw command line arguments anywhere
455:      after PetscInitialize() is called but before PetscFinalize().

457:    Not Collective

459:    Output Parameters:
460: +  argc - count of number of command line arguments
461: -  args - the command line arguments

463:    Level: intermediate

465:    Notes:
466:       This is usually used to pass the command line arguments into other libraries
467:    that are called internally deep in PETSc or the application.

469:       The first argument contains the program name as is normal for C arguments.

471: .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArguments()

473: @*/
474: PetscErrorCode  PetscGetArgs(int *argc,char ***args)
475: {
477:   if (PetscUnlikely(!PetscInitializeCalled && PetscFinalizeCalled)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
478:   *argc = PetscGlobalArgc;
479:   *args = PetscGlobalArgs;
480:   return(0);
481: }

483: /*@C
484:    PetscGetArguments - Allows you to access the  command line arguments anywhere
485:      after PetscInitialize() is called but before PetscFinalize().

487:    Not Collective

489:    Output Parameters:
490: .  args - the command line arguments

492:    Level: intermediate

494:    Notes:
495:       This does NOT start with the program name and IS null terminated (final arg is void)

497: .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscFreeArguments()

499: @*/
500: PetscErrorCode  PetscGetArguments(char ***args)
501: {
502:   PetscInt       i,argc = PetscGlobalArgc;

506:   if (PetscUnlikely(!PetscInitializeCalled && PetscFinalizeCalled)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
507:   if (!argc) {*args = NULL; return(0);}
508:   PetscMalloc1(argc,args);
509:   for (i=0; i<argc-1; i++) {
510:     PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);
511:   }
512:   (*args)[argc-1] = NULL;
513:   return(0);
514: }

516: /*@C
517:    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()

519:    Not Collective

521:    Output Parameters:
522: .  args - the command line arguments

524:    Level: intermediate

526: .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscGetArguments()

528: @*/
529: PetscErrorCode  PetscFreeArguments(char **args)
530: {
531:   PetscInt       i = 0;

535:   if (!args) return(0);
536:   while (args[i]) {
537:     PetscFree(args[i]);
538:     i++;
539:   }
540:   PetscFree(args);
541:   return(0);
542: }

544: #if PetscDefined(HAVE_SAWS)
545: #include <petscconfiginfo.h>

547: PETSC_INTERN PetscErrorCode PetscInitializeSAWs(const char help[])
548: {
550:   if (!PetscGlobalRank) {
551:     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
552:     int            port;
553:     PetscBool      flg,rootlocal = PETSC_FALSE,flg2,selectport = PETSC_FALSE;
554:     size_t         applinelen,introlen;
556:     char           sawsurl[256];

558:     PetscOptionsHasName(NULL,NULL,"-saws_log",&flg);
559:     if (flg) {
560:       char  sawslog[PETSC_MAX_PATH_LEN];

562:       PetscOptionsGetString(NULL,NULL,"-saws_log",sawslog,sizeof(sawslog),NULL);
563:       if (sawslog[0]) {
564:         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
565:       } else {
566:         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
567:       }
568:     }
569:     PetscOptionsGetString(NULL,NULL,"-saws_https",cert,sizeof(cert),&flg);
570:     if (flg) {
571:       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
572:     }
573:     PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select",&selectport,NULL);
574:     if (selectport) {
575:         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
576:         PetscStackCallSAWs(SAWs_Set_Port,(port));
577:     } else {
578:       PetscOptionsGetInt(NULL,NULL,"-saws_port",&port,&flg);
579:       if (flg) {
580:         PetscStackCallSAWs(SAWs_Set_Port,(port));
581:       }
582:     }
583:     PetscOptionsGetString(NULL,NULL,"-saws_root",root,sizeof(root),&flg);
584:     if (flg) {
585:       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
586:       PetscStrcmp(root,".",&rootlocal);
587:     } else {
588:       PetscOptionsHasName(NULL,NULL,"-saws_options",&flg);
589:       if (flg) {
590:         PetscStrreplace(PETSC_COMM_WORLD,"${PETSC_DIR}/share/petsc/saws",root,sizeof(root));
591:         PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
592:       }
593:     }
594:     PetscOptionsHasName(NULL,NULL,"-saws_local",&flg2);
595:     if (flg2) {
596:       char jsdir[PETSC_MAX_PATH_LEN];
597:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
598:       PetscSNPrintf(jsdir,sizeof(jsdir),"%s/js",root);
599:       PetscTestDirectory(jsdir,'r',&flg);
600:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
601:       PetscStackCallSAWs(SAWs_Push_Local_Header,());
602:     }
603:     PetscGetProgramName(programname,sizeof(programname));
604:     PetscStrlen(help,&applinelen);
605:     introlen   = 4096 + applinelen;
606:     applinelen += 1024;
607:     PetscMalloc(applinelen,&appline);
608:     PetscMalloc(introlen,&intro);

610:     if (rootlocal) {
611:       PetscSNPrintf(appline,applinelen,"%s.c.html",programname);
612:       PetscTestFile(appline,'r',&rootlocal);
613:     }
614:     PetscOptionsGetAll(NULL,&options);
615:     if (rootlocal && help) {
616:       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);
617:     } else if (help) {
618:       PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>",programname,options,help);
619:     } else {
620:       PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);
621:     }
622:     PetscFree(options);
623:     PetscGetVersion(version,sizeof(version));
624:     PetscSNPrintf(intro,introlen,"<body>\n"
625:                                     "<center><h2> <a href=\"https://petsc.org/\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n"
626:                                     "<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"
627:                                     "%s",version,petscconfigureoptions,appline);
628:     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
629:     PetscFree(intro);
630:     PetscFree(appline);
631:     if (selectport) {
632:       PetscBool silent;

634:       SAWs_Initialize();
635:       /* another process may have grabbed the port so keep trying */
636:       while (ierr) {
637:         PetscStackCallSAWs(SAWs_Get_Available_Port,(&port));
638:         PetscStackCallSAWs(SAWs_Set_Port,(port));
639:         SAWs_Initialize();
640:       }

642:       PetscOptionsGetBool(NULL,NULL,"-saws_port_auto_select_silent",&silent,NULL);
643:       if (!silent) {
644:         PetscStackCallSAWs(SAWs_Get_FullURL,(sizeof(sawsurl),sawsurl));
645:         PetscPrintf(PETSC_COMM_WORLD,"Point your browser to %s for SAWs\n",sawsurl);
646:       }
647:     } else {
648:       PetscStackCallSAWs(SAWs_Initialize,());
649:     }
650:     PetscCitationsRegister("@TechReport{ saws,\n"
651:                                   "  Author = {Matt Otten and Jed Brown and Barry Smith},\n"
652:                                   "  Title  = {Scientific Application Web Server (SAWs) Users Manual},\n"
653:                                   "  Institution = {Argonne National Laboratory},\n"
654:                                   "  Year   = 2013\n}\n",NULL);
655:   }
656:   return(0);
657: }
658: #endif

660: /* Things must be done before MPI_Init() when MPI is not yet initialized, and can be shared between C init and Fortran init */
661: PETSC_INTERN PetscErrorCode PetscPreMPIInit_Private(void)
662: {
664: #if defined(PETSC_HAVE_HWLOC_SOLARIS_BUG)
665:     /* see MPI.py for details on this bug */
666:     (void) setenv("HWLOC_COMPONENTS","-x86",1);
667: #endif
668:   return(0);
669: }

671: #if defined(PETSC_HAVE_ADIOS)
672: #include <adios.h>
673: #include <adios_read.h>
674: int64_t Petsc_adios_group;
675: #endif
676: #if defined(PETSC_HAVE_OPENMP)
677: #include <omp.h>
678: PetscInt PetscNumOMPThreads;
679: #endif

681: #if PetscDefined(HAVE_DLFCN_H)
682: #include <dlfcn.h>
683: #endif

685: /*
686:   PetscInitialize_Common  - shared code between C and Fortran initialization

688:   prog:     program name
689:   file:     optional PETSc database file name. Might be in Fortran string format when 'ftn' is true
690:   help:     program help message
691:   ftn:      is it called from Fortran initilization (petscinitializef_)?
692:   readarguments,len: used when fortran is true
693: */
694: PETSC_INTERN PetscErrorCode PetscInitialize_Common(const char* prog,const char* file,const char *help,PetscBool ftn,PetscBool readarguments,PetscInt len)
695: {
697:   PetscMPIInt    size;
698:   PetscBool      flg = PETSC_TRUE;
699:   char           hostname[256];

702:   if (PetscInitializeCalled) return(0);
703:   /*
704:       The checking over compatible runtime libraries is complicated by the MPI ABI initiative
705:       https://wiki.mpich.org/mpich/index.php/ABI_Compatibility_Initiative which started with
706:         MPICH v3.1 (Released February 2014)
707:         IBM MPI v2.1 (December 2014)
708:         Intel MPI Library v5.0 (2014)
709:         Cray MPT v7.0.0 (June 2014)
710:       As of July 31, 2017 the ABI number still appears to be 12, that is all of the versions
711:       listed above and since that time are compatible.

713:       Unfortunately the MPI ABI initiative has not defined a way to determine the ABI number
714:       at compile time or runtime. Thus we will need to systematically track the allowed versions
715:       and how they are represented in the mpi.h and MPI_Get_library_version() output in order
716:       to perform the checking.

718:       Currently we only check for pre MPI ABI versions (and packages that do not follow the MPI ABI).

720:       Questions:

722:         Should the checks for ABI incompatibility be only on the major version number below?
723:         Presumably the output to stderr will be removed before a release.
724:   */

726: #if defined(PETSC_HAVE_MPI_GET_LIBRARY_VERSION)
727:   {
728:     char        mpilibraryversion[MPI_MAX_LIBRARY_VERSION_STRING];
729:     PetscMPIInt mpilibraryversionlength;
730:     MPI_Get_library_version(mpilibraryversion,&mpilibraryversionlength);
731:     if (ierr) PetscFunctionReturn(ierr);
732:     /* check for MPICH versions before MPI ABI initiative */
733: #if defined(MPICH_VERSION)
734: #if MPICH_NUMVERSION < 30100000
735:     {
736:       char      *ver,*lf;
737:       PetscBool flg = PETSC_FALSE;
738:       PetscStrstr(mpilibraryversion,"MPICH Version:",&ver);
739:       if (ierr) PetscFunctionReturn(ierr);
740:       else if (ver) {
741:         PetscStrchr(ver,'\n',&lf);
742:         if (ierr) PetscFunctionReturn(ierr);
743:         else if (lf) {
744:           *lf = 0;
745:           PetscStrendswith(ver,MPICH_VERSION,&flg);if (ierr) PetscFunctionReturn(ierr);
746:         }
747:       }
748:       if (!flg) {
749:         PetscInfo2(NULL,"PETSc warning --- MPICH library version \n%s does not match what PETSc was compiled with %s.\n",mpilibraryversion,MPICH_VESION);
750:         flg = PETSC_TRUE;
751:       }
752:     }
753: #endif
754:     /* check for OpenMPI version, it is not part of the MPI ABI initiative (is it part of another initiative that needs to be handled?) */
755: #elif defined(OMPI_MAJOR_VERSION)
756:     {
757:       char *ver,bs[MPI_MAX_LIBRARY_VERSION_STRING],*bsf;
758:       PetscBool flg = PETSC_FALSE;
759: #define PSTRSZ 2
760:       char ompistr1[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"Open MPI","FUJITSU MPI"};
761:       char ompistr2[PSTRSZ][MPI_MAX_LIBRARY_VERSION_STRING] = {"v","Library "};
762:       int i;
763:       for (i=0; i<PSTRSZ; i++) {
764:         PetscStrstr(mpilibraryversion,ompistr1[i],&ver);
765:         if (ierr) PetscFunctionReturn(ierr);
766:         else if (ver) {
767:           PetscSNPrintf(bs,MPI_MAX_LIBRARY_VERSION_STRING,"%s%d.%d",ompistr2[i],OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
768:           PetscStrstr(ver,bs,&bsf);
769:           if (ierr) PetscFunctionReturn(ierr);
770:           else if (bsf) flg = PETSC_TRUE;
771:           break;
772:         }
773:       }
774:       if (!flg) {
775:         PetscInfo3(NULL,"PETSc warning --- Open MPI library version \n%s does not match what PETSc was compiled with %d.%d.\n",mpilibraryversion,OMPI_MAJOR_VERSION,OMPI_MINOR_VERSION);
776:         flg = PETSC_TRUE;
777:       }
778:     }
779: #endif
780:   }
781: #endif

783: #if defined(PETSC_HAVE_DLSYM)
784:   /* These symbols are currently in the OpenMPI and MPICH libraries; they may not always be, in that case the test will simply not detect the problem */
785:   if (PetscUnlikely(dlsym(RTLD_DEFAULT,"ompi_mpi_init") && dlsym(RTLD_DEFAULT,"MPID_Abort"))) {
786:     fprintf(stderr,"PETSc Error --- Application was linked against both OpenMPI and MPICH based MPI libraries and will not run correctly\n");
787:     PetscStackView(stderr);
788:     PetscFunctionReturn(PETSC_ERR_MPI_LIB_INCOMP);
789:   }
790: #endif

792:   /* these must be initialized in a routine, not as a constant declaration*/
793:   PETSC_STDOUT = stdout;
794:   PETSC_STDERR = stderr;

796:   /*CHKERRQ can be used from now */
797:   PetscErrorHandlingInitialized = PETSC_TRUE;

799:   /* on Windows - set printf to default to printing 2 digit exponents */
800: #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
801:   _set_output_format(_TWO_DIGIT_EXPONENT);
802: #endif

804:   PetscOptionsCreateDefault();

806:   PetscFinalizeCalled = PETSC_FALSE;

808:   PetscSetProgramName(prog);
809:   PetscSpinlockCreate(&PetscViewerASCIISpinLockOpen);
810:   PetscSpinlockCreate(&PetscViewerASCIISpinLockStdout);
811:   PetscSpinlockCreate(&PetscViewerASCIISpinLockStderr);
812:   PetscSpinlockCreate(&PetscCommSpinLock);

814:   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
815:   MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);

817:   if (PETSC_MPI_ERROR_CLASS == MPI_ERR_LASTCODE) {
818:     MPI_Add_error_class(&PETSC_MPI_ERROR_CLASS);
819:     MPI_Add_error_code(PETSC_MPI_ERROR_CLASS,&PETSC_MPI_ERROR_CODE);
820:   }

822:   /* Done after init due to a bug in MPICH-GM? */
823:   PetscErrorPrintfInitialize();

825:   MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);
826:   MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);

828:   MPIU_BOOL = MPI_INT;
829:   MPIU_ENUM = MPI_INT;
830:   MPIU_FORTRANADDR = (sizeof(void*) == sizeof(int)) ? MPI_INT : MPIU_INT64;
831:   if (sizeof(size_t) == sizeof(unsigned)) MPIU_SIZE_T = MPI_UNSIGNED;
832:   else if (sizeof(size_t) == sizeof(unsigned long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG;
833: #if defined(PETSC_SIZEOF_LONG_LONG)
834:   else if (sizeof(size_t) == sizeof(unsigned long long)) MPIU_SIZE_T = MPI_UNSIGNED_LONG_LONG;
835: #endif
836:   else {
837:     (*PetscErrorPrintf)("PetscInitialize_Common: Could not find MPI type for size_t\n");
838:     PetscFunctionReturn(PETSC_ERR_SUP_SYS);
839:   }

841:   /*
842:      Initialized the global complex variable; this is because with
843:      shared libraries the constructors for global variables
844:      are not called; at least on IRIX.
845:   */
846: #if defined(PETSC_HAVE_COMPLEX)
847:   {
848: #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_REAL___FLOAT128)
849:     PetscComplex ic(0.0,1.0);
850:     PETSC_i = ic;
851: #else
852:     PETSC_i = _Complex_I;
853: #endif
854:   }
855: #endif /* PETSC_HAVE_COMPLEX */

857:   /*
858:      Create the PETSc MPI reduction operator that sums of the first
859:      half of the entries and maxes the second half.
860:   */
861:   MPI_Op_create(MPIU_MaxSum_Local,1,&MPIU_MAXSUM_OP);

863: #if defined(PETSC_USE_REAL___FLOAT128)
864:   MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);
865:   MPI_Type_commit(&MPIU___FLOAT128);
866: #if defined(PETSC_HAVE_COMPLEX)
867:   MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);
868:   MPI_Type_commit(&MPIU___COMPLEX128);
869: #endif
870:   MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);
871:   MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);
872: #elif defined(PETSC_USE_REAL___FP16)
873:   MPI_Type_contiguous(2,MPI_CHAR,&MPIU___FP16);
874:   MPI_Type_commit(&MPIU___FP16);
875:   MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);
876:   MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);
877: #endif

879: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
880:   MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);
881: #endif

883:   MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);
884:   MPI_Type_commit(&MPIU_2SCALAR);

886:   /* create datatypes used by MPIU_MAXLOC, MPIU_MINLOC and PetscSplitReduction_Op */
887: #if !defined(PETSC_HAVE_MPIUNI)
888:   {
889:     struct PetscRealInt { PetscReal v; PetscInt i; };
890:     PetscMPIInt  blockSizes[2] = {1,1};
891:     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscRealInt,v),offsetof(struct PetscRealInt,i)};
892:     MPI_Datatype blockTypes[2] = {MPIU_REAL,MPIU_INT}, tmpStruct;

894:     MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct);
895:     MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscRealInt),&MPIU_REAL_INT);
896:     MPI_Type_free(&tmpStruct);
897:     MPI_Type_commit(&MPIU_REAL_INT);
898:   }
899:   {
900:     struct PetscScalarInt { PetscScalar v; PetscInt i; };
901:     PetscMPIInt  blockSizes[2] = {1,1};
902:     MPI_Aint     blockOffsets[2] = {offsetof(struct PetscScalarInt,v),offsetof(struct PetscScalarInt,i)};
903:     MPI_Datatype blockTypes[2] = {MPIU_SCALAR,MPIU_INT}, tmpStruct;

905:     MPI_Type_create_struct(2,blockSizes,blockOffsets,blockTypes,&tmpStruct);
906:     MPI_Type_create_resized(tmpStruct,0,sizeof(struct PetscScalarInt),&MPIU_SCALAR_INT);
907:     MPI_Type_free(&tmpStruct);
908:     MPI_Type_commit(&MPIU_SCALAR_INT);
909:   }
910: #endif

912: #if defined(PETSC_USE_64BIT_INDICES)
913:   MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);
914:   MPI_Type_commit(&MPIU_2INT);
915: #endif
916:   MPI_Type_contiguous(4,MPI_INT,&MPI_4INT);
917:   MPI_Type_commit(&MPI_4INT);
918:   MPI_Type_contiguous(4,MPIU_INT,&MPIU_4INT);
919:   MPI_Type_commit(&MPIU_4INT);

921:   /*
922:      Attributes to be set on PETSc communicators
923:   */
924:   MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_Counter_Attr_Delete_Fn,&Petsc_Counter_keyval,(void*)0);
925:   MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_InnerComm_Attr_Delete_Fn,&Petsc_InnerComm_keyval,(void*)0);
926:   MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_OuterComm_Attr_Delete_Fn,&Petsc_OuterComm_keyval,(void*)0);
927:   MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,Petsc_ShmComm_Attr_Delete_Fn,&Petsc_ShmComm_keyval,(void*)0);

929: #if defined(PETSC_HAVE_FORTRAN)
930:   if (ftn) {PetscInitFortran_Private(readarguments,file,len);}
931:   else
932: #endif
933:   {PetscOptionsInsert(NULL,&PetscGlobalArgc,&PetscGlobalArgs,file);}

935:   /* call a second time so it can look in the options database */
936:   PetscErrorPrintfInitialize();

938:   /*
939:      Check system options and print help
940:   */
941:   PetscOptionsCheckInitial_Private(help);

943:   PetscCitationsInitialize();

945: #if defined(PETSC_HAVE_SAWS)
946:   PetscInitializeSAWs(ftn ? NULL : help);
947:   flg = PETSC_FALSE;
948:   PetscOptionsHasName(NULL,NULL,"-stack_view",&flg);
949:   if (flg) PetscStackViewSAWs();
950: #endif

952:   /*
953:      Load the dynamic libraries (on machines that support them), this registers all
954:      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
955:   */
956:   PetscInitialize_DynamicLibraries();

958:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
959:   PetscInfo1(NULL,"PETSc successfully started: number of processors = %d\n",size);
960:   PetscGetHostName(hostname,256);
961:   PetscInfo1(NULL,"Running on machine: %s\n",hostname);
962: #if defined(PETSC_HAVE_OPENMP)
963:   {
964:     PetscBool omp_view_flag;
965:     char      *threads = getenv("OMP_NUM_THREADS");

967:     if (threads) {
968:       PetscInfo1(NULL,"Number of OpenMP threads %s (as given by OMP_NUM_THREADS)\n",threads);
969:       (void) sscanf(threads, "%" PetscInt_FMT,&PetscNumOMPThreads);
970:     } else {
971:       PetscNumOMPThreads = (PetscInt) omp_get_max_threads();
972:       PetscInfo1(NULL,"Number of OpenMP threads %D (as given by omp_get_max_threads())\n",PetscNumOMPThreads);
973:     }
974:     PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"OpenMP options","Sys");
975:     PetscOptionsInt("-omp_num_threads","Number of OpenMP threads to use (can also use environmental variable OMP_NUM_THREADS","None",PetscNumOMPThreads,&PetscNumOMPThreads,&flg);
976:     PetscOptionsName("-omp_view","Display OpenMP number of threads",NULL,&omp_view_flag);
977:     PetscOptionsEnd();
978:     if (flg) {
979:       PetscInfo1(NULL,"Number of OpenMP theads %D (given by -omp_num_threads)\n",PetscNumOMPThreads);
980:       omp_set_num_threads((int)PetscNumOMPThreads);
981:     }
982:     if (omp_view_flag) {
983:       PetscPrintf(PETSC_COMM_WORLD,"OpenMP: number of threads %D\n",PetscNumOMPThreads);
984:     }
985:   }
986: #endif

988: #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
989:   /*
990:       Tell MPI about our own data representation converter, this would/should be used if extern32 is not supported by the MPI

992:       Currently not used because it is not supported by MPICH.
993:   */
994:   if (!PetscBinaryBigEndian()) {
995:     MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);
996:   }
997: #endif

999: #if defined(PETSC_SERIALIZE_FUNCTIONS)
1000:   PetscFPTCreate(10000);
1001: #endif

1003: #if defined(PETSC_HAVE_HWLOC)
1004:   {
1005:     PetscViewer viewer;
1006:     PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-process_view",&viewer,NULL,&flg);
1007:     if (flg) {
1008:       PetscProcessPlacementView(viewer);
1009:       PetscViewerDestroy(&viewer);
1010:     }
1011:   }
1012: #endif

1014:   flg  = PETSC_TRUE;
1015:   PetscOptionsGetBool(NULL,NULL,"-viewfromoptions",&flg,NULL);
1016:   if (!flg) {PetscOptionsPushGetViewerOff(PETSC_TRUE);}

1018: #if defined(PETSC_HAVE_ADIOS)
1019:   adios_init_noxml(PETSC_COMM_WORLD);
1020:   adios_declare_group(&Petsc_adios_group,"PETSc","",adios_stat_default);
1021:   adios_select_method(Petsc_adios_group,"MPI","","");
1022:   adios_read_init_method(ADIOS_READ_METHOD_BP,PETSC_COMM_WORLD,"");
1023: #endif

1025: #if defined(__VALGRIND_H)
1026:   PETSC_RUNNING_ON_VALGRIND = RUNNING_ON_VALGRIND? PETSC_TRUE: PETSC_FALSE;
1027: #if defined(PETSC_USING_DARWIN) && defined(PETSC_BLASLAPACK_SDOT_RETURNS_DOUBLE)
1028:   if (PETSC_RUNNING_ON_VALGRIND) {
1029:     PetscPrintf(PETSC_COMM_WORLD,"WARNING: Running valgrind with the MacOS native BLAS and LAPACK can fail. If it fails suggest configuring with --download-fblaslapack or --download-f2cblaslapack");
1030:     }
1031: #endif
1032: #endif

1034: #if (defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)) && defined(PETSC_EXPERIMENTAL)
1035:   PetscDeviceInitializeDefaultDevices_Internal();
1036:   PetscDeviceContextInitializeRootContext_Internal(PETSC_COMM_WORLD,NULL);
1037: #endif

1039:   /*
1040:       Set flag that we are completely initialized
1041:   */
1042:   PetscInitializeCalled = PETSC_TRUE;

1044:   PetscOptionsHasName(NULL,NULL,"-python",&flg);
1045:   if (flg) {PetscPythonInitialize(NULL,NULL);}
1046:   return(0);
1047: }

1049: /*@C
1050:    PetscInitialize - Initializes the PETSc database and MPI.
1051:    PetscInitialize() calls MPI_Init() if that has yet to be called,
1052:    so this routine should always be called near the beginning of
1053:    your program -- usually the very first line!

1055:    Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set

1057:    Input Parameters:
1058: +  argc - count of number of command line arguments
1059: .  args - the command line arguments
1060: .  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
1061:           Use NULL or empty string to not check for code specific file.
1062:           Also checks ~/.petscrc, .petscrc and petscrc.
1063:           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
1064: -  help - [optional] Help message to print, use NULL for no message

1066:    If you wish PETSc code to run ONLY on a subcommunicator of MPI_COMM_WORLD, create that
1067:    communicator first and assign it to PETSC_COMM_WORLD BEFORE calling PetscInitialize(). Thus if you are running a
1068:    four process job and two processes will run PETSc and have PetscInitialize() and PetscFinalize() and two process will not,
1069:    then do this. If ALL processes in the job are using PetscInitialize() and PetscFinalize() then you don't need to do this, even
1070:    if different subcommunicators of the job are doing different things with PETSc.

1072:    Options Database Keys:
1073: +  -help [intro] - prints help method for each option; if intro is given the program stops after printing the introductory help message
1074: .  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
1075: .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
1076: .  -on_error_emacs <machinename> - causes emacsclient to jump to error file
1077: .  -on_error_abort - calls abort() when error detected (no traceback)
1078: .  -on_error_mpiabort - calls MPI_abort() when error detected
1079: .  -error_output_stderr - prints error messages to stderr instead of the default stdout
1080: .  -error_output_none - does not print the error messages (but handles errors in the same way as if this was not called)
1081: .  -debugger_ranks [rank1,rank2,...] - Indicates ranks to start in debugger
1082: .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
1083: .  -stop_for_debugger - Print message on how to attach debugger manually to
1084:                         process and wait (-debugger_pause) seconds for attachment
1085: .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries) (deprecated, use -malloc_debug)
1086: .  -malloc no - Indicates not to use error-checking malloc (deprecated, use -malloc_debug no)
1087: .  -malloc_debug - check for memory corruption at EVERY malloc or free, see PetscMallocSetDebug()
1088: .  -malloc_dump - prints a list of all unfreed memory at the end of the run
1089: .  -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
1090: .  -malloc_view - show a list of all allocated memory during PetscFinalize()
1091: .  -malloc_view_threshold <t> - only list memory allocations of size greater than t with -malloc_view
1092: .  -malloc_requested_size - malloc logging will record the requested size rather than size after alignment
1093: .  -fp_trap - Stops on floating point exceptions
1094: .  -no_signal_handler - Indicates not to trap error signals
1095: .  -shared_tmp - indicates /tmp directory is shared by all processors
1096: .  -not_shared_tmp - each processor has own /tmp
1097: .  -tmp - alternative name of /tmp directory
1098: .  -get_total_flops - returns total flops done by all processors
1099: -  -memory_view - Print memory usage at end of run

1101:    Options Database Keys for Option Database:
1102: +  -skip_petscrc - skip the default option files ~/.petscrc, .petscrc, petscrc
1103: .  -options_monitor - monitor all set options to standard output for the whole program run
1104: -  -options_monitor_cancel - cancel options monitoring hard-wired using PetscOptionsMonitorSet()

1106:    Options -options_monitor_{all,cancel} are
1107:    position-independent and apply to all options set since the PETSc start.
1108:    They can be used also in option files.

1110:    See PetscOptionsMonitorSet() to do monitoring programmatically.

1112:    Options Database Keys for Profiling:
1113:    See Users-Manual: ch_profiling for details.
1114: +  -info [filename][:[~]<list,of,classnames>[:[~]self]] - Prints verbose information. See PetscInfo().
1115: .  -log_sync - Enable barrier synchronization for all events. This option is useful to debug imbalance within each event,
1116:         however it slows things down and gives a distorted view of the overall runtime.
1117: .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
1118:         hangs without running in the debugger).  See PetscLogTraceBegin().
1119: .  -log_view [:filename:format] - Prints summary of flop and timing information to screen or file, see PetscLogView().
1120: .  -log_view_memory - Includes in the summary from -log_view the memory used in each method, see PetscLogView().
1121: .  -log_summary [filename] - (Deprecated, use -log_view) Prints summary of flop and timing information to screen. If the filename is specified the
1122:         summary is written to the file.  See PetscLogView().
1123: .  -log_exclude: <vec,mat,pc,ksp,snes> - excludes subset of object classes from logging
1124: .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
1125: .  -log [filename] - Logs basic profiline information  See PetscLogDump().
1126: .  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)
1127: .  -viewfromoptions on,off - Enable or disable XXXSetFromOptions() calls, for applications with many small solves turn this off
1128: -  -check_pointer_intensity 0,1,2 - if pointers are checked for validity (debug version only), using 0 will result in faster code

1130:     Only one of -log_trace, -log_view, -log_all, -log, or -log_mpe may be used at a time

1132:    Options Database Keys for SAWs:
1133: +  -saws_port <portnumber> - port number to publish SAWs data, default is 8080
1134: .  -saws_port_auto_select - have SAWs select a new unique port number where it publishes the data, the URL is printed to the screen
1135:                             this is useful when you are running many jobs that utilize SAWs at the same time
1136: .  -saws_log <filename> - save a log of all SAWs communication
1137: .  -saws_https <certificate file> - have SAWs use HTTPS instead of HTTP
1138: -  -saws_root <directory> - allow SAWs to have access to the given directory to search for requested resources and files

1140:    Environmental Variables:
1141: +   PETSC_TMP - alternative tmp directory
1142: .   PETSC_SHARED_TMP - tmp is shared by all processes
1143: .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
1144: .   PETSC_OPTIONS - a string containing additional options for petsc in the form of command line "-key value" pairs
1145: .   PETSC_OPTIONS_YAML - (requires configuring PETSc to use libyaml) a string containing additional options for petsc in the form of a YAML document
1146: .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
1147: -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to

1149:    Level: beginner

1151:    Notes:
1152:    If for some reason you must call MPI_Init() separately, call
1153:    it before PetscInitialize().

1155:    Fortran Version:
1156:    In Fortran this routine has the format
1157: $       call PetscInitialize(file,ierr)

1159: +  ierr - error return code
1160: -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
1161:           Use PETSC_NULL_CHARACTER to not check for code specific file.
1162:           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.

1164:    Important Fortran Note:
1165:    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
1166:    null character string; you CANNOT just use NULL as
1167:    in the C version. See Users-Manual: ch_fortran for details.

1169:    If your main program is C but you call Fortran code that also uses PETSc you need to call PetscInitializeFortran() soon after
1170:    calling PetscInitialize().

1172: .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()

1174: @*/
1175: PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
1176: {
1178:   PetscMPIInt    flag;
1179:   const char     *prog = "Unknown Name";

1182:   if (PetscInitializeCalled) return(0);
1183:   MPI_Initialized(&flag);
1184:   if (!flag) {
1185:     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");
1186:     PetscPreMPIInit_Private();
1187: #if defined(PETSC_HAVE_MPI_INIT_THREAD)
1188:     {
1189:       PetscMPIInt provided;
1190:       MPI_Init_thread(argc,args,PETSC_MPI_THREAD_REQUIRED,&provided);
1191:     }
1192: #else
1193:     MPI_Init(argc,args);
1194: #endif
1195:     PetscBeganMPI = PETSC_TRUE;
1196:   }

1198:   if (argc && *argc) prog = **args;
1199:   if (argc && args) {
1200:     PetscGlobalArgc = *argc;
1201:     PetscGlobalArgs = *args;
1202:   }
1203:   PetscInitialize_Common(prog,file,help,PETSC_FALSE/*C*/,PETSC_FALSE,0);
1204:   return(0);
1205: }

1207: #if defined(PETSC_USE_LOG)
1208: PETSC_INTERN PetscObject *PetscObjects;
1209: PETSC_INTERN PetscInt    PetscObjectsCounts;
1210: PETSC_INTERN PetscInt    PetscObjectsMaxCounts;
1211: PETSC_INTERN PetscBool   PetscObjectsLog;
1212: #endif

1214: /*
1215:     Frees all the MPI types and operations that PETSc may have created
1216: */
1217: PetscErrorCode  PetscFreeMPIResources(void)
1218: {

1222: #if defined(PETSC_USE_REAL___FLOAT128)
1223:   MPI_Type_free(&MPIU___FLOAT128);
1224: #if defined(PETSC_HAVE_COMPLEX)
1225:   MPI_Type_free(&MPIU___COMPLEX128);
1226: #endif
1227:   MPI_Op_free(&MPIU_MAX);
1228:   MPI_Op_free(&MPIU_MIN);
1229: #elif defined(PETSC_USE_REAL___FP16)
1230:   MPI_Type_free(&MPIU___FP16);
1231:   MPI_Op_free(&MPIU_MAX);
1232:   MPI_Op_free(&MPIU_MIN);
1233: #endif

1235: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1236:   MPI_Op_free(&MPIU_SUM);
1237: #endif

1239:   MPI_Type_free(&MPIU_2SCALAR);
1240:   MPI_Type_free(&MPIU_REAL_INT);
1241:   MPI_Type_free(&MPIU_SCALAR_INT);
1242: #if defined(PETSC_USE_64BIT_INDICES)
1243:   MPI_Type_free(&MPIU_2INT);
1244: #endif
1245:   MPI_Type_free(&MPI_4INT);
1246:   MPI_Type_free(&MPIU_4INT);
1247:   MPI_Op_free(&MPIU_MAXSUM_OP);
1248:   return(0);
1249: }

1251: /*@C
1252:    PetscFinalize - Checks for options to be called at the conclusion
1253:    of the program. MPI_Finalize() is called only if the user had not
1254:    called MPI_Init() before calling PetscInitialize().

1256:    Collective on PETSC_COMM_WORLD

1258:    Options Database Keys:
1259: +  -options_view - Calls PetscOptionsView()
1260: .  -options_left - Prints unused options that remain in the database
1261: .  -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
1262: .  -mpidump - Calls PetscMPIDump()
1263: .  -malloc_dump <optional filename> - Calls PetscMallocDump(), displays all memory allocated that has not been freed
1264: .  -malloc_info - Prints total memory usage
1265: -  -malloc_view <optional filename> - Prints list of all memory allocated and where

1267:    Level: beginner

1269:    Note:
1270:    See PetscInitialize() for more general runtime options.

1272: .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
1273: @*/
1274: PetscErrorCode  PetscFinalize(void)
1275: {
1277:   PetscMPIInt    rank;
1278:   PetscInt       nopt;
1279:   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
1280:   PetscBool      flg;
1281: #if defined(PETSC_USE_LOG)
1282:   char           mname[PETSC_MAX_PATH_LEN];
1283: #endif

1286:   if (PetscUnlikely(!PetscInitializeCalled)) {
1287:     fprintf(PETSC_STDOUT,"PetscInitialize() must be called before PetscFinalize()\n");
1288:     PetscStackView(PETSC_STDOUT);
1289:     PetscStackClearTop;
1290:     return PETSC_ERR_ARG_WRONGSTATE;
1291:   }
1292:   PetscInfo(NULL,"PetscFinalize() called\n");

1294:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
1295: #if defined(PETSC_HAVE_ADIOS)
1296:   adios_read_finalize_method(ADIOS_READ_METHOD_BP_AGGREGATE);
1297:   adios_finalize(rank);
1298: #endif
1299:   PetscOptionsHasName(NULL,NULL,"-citations",&flg);
1300:   if (flg) {
1301:     char  *cits, filename[PETSC_MAX_PATH_LEN];
1302:     FILE  *fd = PETSC_STDOUT;

1304:     PetscOptionsGetString(NULL,NULL,"-citations",filename,sizeof(filename),NULL);
1305:     if (filename[0]) {
1306:       PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);
1307:     }
1308:     PetscSegBufferGet(PetscCitationsList,1,&cits);
1309:     cits[0] = 0;
1310:     PetscSegBufferExtractAlloc(PetscCitationsList,&cits);
1311:     PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");
1312:     PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");
1313:     PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);
1314:     PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");
1315:     PetscFClose(PETSC_COMM_WORLD,fd);
1316:     PetscFree(cits);
1317:   }
1318:   PetscSegBufferDestroy(&PetscCitationsList);

1320: #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
1321:   /* TextBelt is run for testing purposes only, please do not use this feature often */
1322:   {
1323:     PetscInt nmax = 2;
1324:     char     **buffs;
1325:     PetscMalloc1(2,&buffs);
1326:     PetscOptionsGetStringArray(NULL,NULL,"-textbelt",buffs,&nmax,&flg1);
1327:     if (flg1) {
1328:       if (PetscUnlikely(!nmax)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
1329:       if (nmax == 1) {
1330:         PetscMalloc1(128,&buffs[1]);
1331:         PetscGetProgramName(buffs[1],32);
1332:         PetscStrcat(buffs[1]," has completed");
1333:       }
1334:       PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);
1335:       PetscFree(buffs[0]);
1336:       PetscFree(buffs[1]);
1337:     }
1338:     PetscFree(buffs);
1339:   }
1340:   {
1341:     PetscInt nmax = 2;
1342:     char     **buffs;
1343:     PetscMalloc1(2,&buffs);
1344:     PetscOptionsGetStringArray(NULL,NULL,"-tellmycell",buffs,&nmax,&flg1);
1345:     if (flg1) {
1346:       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-tellmycell requires either the phone number or number,\"message\"");
1347:       if (nmax == 1) {
1348:         PetscMalloc1(128,&buffs[1]);
1349:         PetscGetProgramName(buffs[1],32);
1350:         PetscStrcat(buffs[1]," has completed");
1351:       }
1352:       PetscTellMyCell(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);
1353:       PetscFree(buffs[0]);
1354:       PetscFree(buffs[1]);
1355:     }
1356:     PetscFree(buffs);
1357:   }
1358: #endif

1360: #if defined(PETSC_SERIALIZE_FUNCTIONS)
1361:   PetscFPTDestroy();
1362: #endif

1364: #if defined(PETSC_HAVE_SAWS)
1365:   flg = PETSC_FALSE;
1366:   PetscOptionsGetBool(NULL,NULL,"-saw_options",&flg,NULL);
1367:   if (flg) {
1368:     PetscOptionsSAWsDestroy();
1369:   }
1370: #endif

1372: #if defined(PETSC_HAVE_X)
1373:   flg1 = PETSC_FALSE;
1374:   PetscOptionsGetBool(NULL,NULL,"-x_virtual",&flg1,NULL);
1375:   if (flg1) {
1376:     /*  this is a crude hack, but better than nothing */
1377:     PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);
1378:   }
1379: #endif

1381: #if !defined(PETSC_HAVE_THREADSAFETY)
1382:   PetscOptionsGetBool(NULL,NULL,"-malloc_info",&flg2,NULL);
1383:   if (!flg2) {
1384:     flg2 = PETSC_FALSE;
1385:     PetscOptionsGetBool(NULL,NULL,"-memory_view",&flg2,NULL);
1386:   }
1387:   if (flg2) {
1388:     PetscMemoryView(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");
1389:   }
1390: #endif

1392: #if defined(PETSC_USE_LOG)
1393:   flg1 = PETSC_FALSE;
1394:   PetscOptionsGetBool(NULL,NULL,"-get_total_flops",&flg1,NULL);
1395:   if (flg1) {
1396:     PetscLogDouble flops = 0;
1397:     MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);
1398:     PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);
1399:   }
1400: #endif

1402: #if defined(PETSC_USE_LOG)
1403: #if defined(PETSC_HAVE_MPE)
1404:   mname[0] = 0;
1405:   PetscOptionsGetString(NULL,NULL,"-log_mpe",mname,sizeof(mname),&flg1);
1406:   if (flg1) {
1407:     if (mname[0]) {PetscLogMPEDump(mname);}
1408:     else          {PetscLogMPEDump(0);}
1409:   }
1410: #endif
1411: #endif

1413:   /*
1414:      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1415:   */
1416:   PetscObjectRegisterDestroyAll();

1418: #if defined(PETSC_USE_LOG)
1419:   PetscOptionsPushGetViewerOff(PETSC_FALSE);
1420:   PetscLogViewFromOptions();
1421:   PetscOptionsPopGetViewerOff();

1423:   mname[0] = 0;
1424:   PetscOptionsGetString(NULL,NULL,"-log_summary",mname,sizeof(mname),&flg1);
1425:   if (flg1) {
1426:     PetscViewer viewer;
1427:     (*PetscHelpPrintf)(PETSC_COMM_WORLD,"\n\n WARNING:   -log_summary is being deprecated; switch to -log_view\n\n\n");
1428:     if (mname[0]) {
1429:       PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);
1430:       PetscLogView(viewer);
1431:       PetscViewerDestroy(&viewer);
1432:     } else {
1433:       viewer = PETSC_VIEWER_STDOUT_WORLD;
1434:       PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);
1435:       PetscLogView(viewer);
1436:       PetscViewerPopFormat(viewer);
1437:     }
1438:   }

1440:   /*
1441:      Free any objects created by the last block of code.
1442:   */
1443:   PetscObjectRegisterDestroyAll();

1445:   mname[0] = 0;
1446:   PetscOptionsGetString(NULL,NULL,"-log_all",mname,sizeof(mname),&flg1);
1447:   PetscOptionsGetString(NULL,NULL,"-log",mname,sizeof(mname),&flg2);
1448:   if (flg1 || flg2) {PetscLogDump(mname);}
1449: #endif

1451:   flg1 = PETSC_FALSE;
1452:   PetscOptionsGetBool(NULL,NULL,"-no_signal_handler",&flg1,NULL);
1453:   if (!flg1) { PetscPopSignalHandler();}
1454:   flg1 = PETSC_FALSE;
1455:   PetscOptionsGetBool(NULL,NULL,"-mpidump",&flg1,NULL);
1456:   if (flg1) {
1457:     PetscMPIDump(stdout);
1458:   }
1459:   flg1 = PETSC_FALSE;
1460:   flg2 = PETSC_FALSE;
1461:   /* preemptive call to avoid listing this option in options table as unused */
1462:   PetscOptionsHasName(NULL,NULL,"-malloc_dump",&flg1);
1463:   PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);
1464:   PetscOptionsGetBool(NULL,NULL,"-options_view",&flg2,NULL);

1466:   if (flg2) {
1467:     PetscViewer viewer;
1468:     PetscViewerCreate(PETSC_COMM_WORLD,&viewer);
1469:     PetscViewerSetType(viewer,PETSCVIEWERASCII);
1470:     PetscOptionsView(NULL,viewer);
1471:     PetscViewerDestroy(&viewer);
1472:   }

1474:   /* to prevent PETSc -options_left from warning */
1475:   PetscOptionsHasName(NULL,NULL,"-nox",&flg1);
1476:   PetscOptionsHasName(NULL,NULL,"-nox_warning",&flg1);

1478:   flg3 = PETSC_FALSE; /* default value is required */
1479:   PetscOptionsGetBool(NULL,NULL,"-options_left",&flg3,&flg1);
1480:   if (PetscUnlikelyDebug(!flg1)) flg3 = PETSC_TRUE;
1481:   if (flg3) {
1482:     if (!flg2 && flg1) { /* have not yet printed the options */
1483:       PetscViewer viewer;
1484:       PetscViewerCreate(PETSC_COMM_WORLD,&viewer);
1485:       PetscViewerSetType(viewer,PETSCVIEWERASCII);
1486:       PetscOptionsView(NULL,viewer);
1487:       PetscViewerDestroy(&viewer);
1488:     }
1489:     PetscOptionsAllUsed(NULL,&nopt);
1490:     if (nopt) {
1491:       PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");
1492:       PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");
1493:       if (nopt == 1) {
1494:         PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");
1495:       } else {
1496:         PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);
1497:       }
1498:     } else if (flg3 && flg1) {
1499:       PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");
1500:     }
1501:     PetscOptionsLeft(NULL);
1502:   }

1504: #if defined(PETSC_HAVE_SAWS)
1505:   if (!PetscGlobalRank) {
1506:     PetscStackSAWsViewOff();
1507:     PetscStackCallSAWs(SAWs_Finalize,());
1508:   }
1509: #endif

1511: #if defined(PETSC_USE_LOG)
1512:   /*
1513:        List all objects the user may have forgot to free
1514:   */
1515:   if (PetscObjectsLog) {
1516:     PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);
1517:     if (flg1) {
1518:       MPI_Comm local_comm;
1519:       char     string[64];

1521:       PetscOptionsGetString(NULL,NULL,"-objects_dump",string,sizeof(string),NULL);
1522:       MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);
1523:       PetscSequentialPhaseBegin_Private(local_comm,1);
1524:       PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);
1525:       PetscSequentialPhaseEnd_Private(local_comm,1);
1526:       MPI_Comm_free(&local_comm);
1527:     }
1528:   }
1529: #endif

1531: #if defined(PETSC_USE_LOG)
1532:   PetscObjectsCounts    = 0;
1533:   PetscObjectsMaxCounts = 0;
1534:   PetscFree(PetscObjects);
1535: #endif

1537:   /*
1538:      Destroy any packages that registered a finalize
1539:   */
1540:   PetscRegisterFinalizeAll();

1542: #if defined(PETSC_USE_LOG)
1543:   PetscLogFinalize();
1544: #endif

1546:   /*
1547:      Print PetscFunctionLists that have not been properly freed

1549:   PetscFunctionListPrintAll();
1550:   */

1552:   if (petsc_history) {
1553:     PetscCloseHistoryFile(&petsc_history);
1554:     petsc_history = NULL;
1555:   }
1556:   PetscOptionsHelpPrintedDestroy(&PetscOptionsHelpPrintedSingleton);
1557:   PetscInfoDestroy();

1559: #if !defined(PETSC_HAVE_THREADSAFETY)
1560:   if (!(PETSC_RUNNING_ON_VALGRIND)) {
1561:     char fname[PETSC_MAX_PATH_LEN];
1562:     char sname[PETSC_MAX_PATH_LEN];
1563:     FILE *fd;
1564:     int  err;

1566:     flg2 = PETSC_FALSE;
1567:     flg3 = PETSC_FALSE;
1568:     if (PetscDefined(USE_DEBUG)) {PetscOptionsGetBool(NULL,NULL,"-malloc_test",&flg2,NULL);}
1569:     PetscOptionsGetBool(NULL,NULL,"-malloc_debug",&flg3,NULL);
1570:     fname[0] = 0;
1571:     PetscOptionsGetString(NULL,NULL,"-malloc_dump",fname,sizeof(fname),&flg1);
1572:     if (flg1 && fname[0]) {

1574:       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1575:       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1576:       PetscMallocDump(fd);
1577:       err  = fclose(fd);
1578:       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1579:     } else if (flg1 || flg2 || flg3) {
1580:       MPI_Comm local_comm;

1582:       MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);
1583:       PetscSequentialPhaseBegin_Private(local_comm,1);
1584:       PetscMallocDump(stdout);
1585:       PetscSequentialPhaseEnd_Private(local_comm,1);
1586:       MPI_Comm_free(&local_comm);
1587:     }
1588:     fname[0] = 0;
1589:     PetscOptionsGetString(NULL,NULL,"-malloc_view",fname,sizeof(fname),&flg1);
1590:     if (flg1 && fname[0]) {

1592:       PetscSNPrintf(sname,sizeof(sname),"%s_%d",fname,rank);
1593:       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1594:       PetscMallocView(fd);
1595:       err  = fclose(fd);
1596:       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1597:     } else if (flg1) {
1598:       MPI_Comm local_comm;

1600:       MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);
1601:       PetscSequentialPhaseBegin_Private(local_comm,1);
1602:       PetscMallocView(stdout);
1603:       PetscSequentialPhaseEnd_Private(local_comm,1);
1604:       MPI_Comm_free(&local_comm);
1605:     }
1606:   }
1607: #endif

1609:   /*
1610:      Close any open dynamic libraries
1611:   */
1612:   PetscFinalize_DynamicLibraries();

1614:   /* Can be destroyed only after all the options are used */
1615:   PetscOptionsDestroyDefault();

1617:   PetscGlobalArgc = 0;
1618:   PetscGlobalArgs = NULL;

1620: #if defined(PETSC_HAVE_KOKKOS)
1621:   if (PetscBeganKokkos) {
1622:     PetscKokkosFinalize_Private();
1623:     PetscBeganKokkos = PETSC_FALSE;
1624:     PetscKokkosInitialized = PETSC_FALSE;
1625:   }
1626: #endif

1628: #if defined(PETSC_HAVE_NVSHMEM)
1629:   if (PetscBeganNvshmem) {
1630:     PetscNvshmemFinalize();
1631:     PetscBeganNvshmem = PETSC_FALSE;
1632:   }
1633: #endif

1635: #if defined(PETSC_HAVE_CUDA)
1636:   if (PetscDefaultCudaStream) {cudaError_t cerr = cudaStreamDestroy(PetscDefaultCudaStream);CHKERRCUDA(cerr);}
1637:   if (petsc_gputimer_begin) {
1638:     cudaError_t cerr = cudaEventDestroy(petsc_gputimer_begin);CHKERRCUDA(cerr);
1639:   }
1640:   if (petsc_gputimer_end) {
1641:     cudaError_t cerr = cudaEventDestroy(petsc_gputimer_end);CHKERRCUDA(cerr);
1642:   }
1643: #endif

1645: #if defined(PETSC_HAVE_HIP)
1646:   if (PetscDefaultHipStream)  {hipError_t cerr  = hipStreamDestroy(PetscDefaultHipStream);CHKERRHIP(cerr);}
1647:   if (petsc_gputimer_begin) {
1648:     hipError_t cerr = hipEventDestroy(petsc_gputimer_begin);CHKERRHIP(cerr);
1649:   }
1650:   if (petsc_gputimer_end) {
1651:     hipError_t cerr = hipEventDestroy(petsc_gputimer_end);CHKERRHIP(cerr);
1652:   }
1653: #endif

1655:   PetscFreeMPIResources();

1657:   /*
1658:      Destroy any known inner MPI_Comm's and attributes pointing to them
1659:      Note this will not destroy any new communicators the user has created.

1661:      If all PETSc objects were not destroyed those left over objects will have hanging references to
1662:      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1663:  */
1664:   {
1665:     PetscCommCounter *counter;
1666:     PetscMPIInt      flg;
1667:     MPI_Comm         icomm;
1668:     union {MPI_Comm comm; void *ptr;} ucomm;
1669:     MPI_Comm_get_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);
1670:     if (flg) {
1671:       icomm = ucomm.comm;
1672:       MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);
1673:       if (PetscUnlikely(!flg)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");

1675:       MPI_Comm_delete_attr(PETSC_COMM_SELF,Petsc_InnerComm_keyval);
1676:       MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);
1677:       MPI_Comm_free(&icomm);
1678:     }
1679:     MPI_Comm_get_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);
1680:     if (flg) {
1681:       icomm = ucomm.comm;
1682:       MPI_Comm_get_attr(icomm,Petsc_Counter_keyval,&counter,&flg);
1683:       if (PetscUnlikely(!flg)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");

1685:       MPI_Comm_delete_attr(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);
1686:       MPI_Comm_delete_attr(icomm,Petsc_Counter_keyval);
1687:       MPI_Comm_free(&icomm);
1688:     }
1689:   }

1691:   MPI_Comm_free_keyval(&Petsc_Counter_keyval);
1692:   MPI_Comm_free_keyval(&Petsc_InnerComm_keyval);
1693:   MPI_Comm_free_keyval(&Petsc_OuterComm_keyval);
1694:   MPI_Comm_free_keyval(&Petsc_ShmComm_keyval);

1696:   PetscSpinlockDestroy(&PetscViewerASCIISpinLockOpen);
1697:   PetscSpinlockDestroy(&PetscViewerASCIISpinLockStdout);
1698:   PetscSpinlockDestroy(&PetscViewerASCIISpinLockStderr);
1699:   PetscSpinlockDestroy(&PetscCommSpinLock);

1701:   if (PetscBeganMPI) {
1702: #if defined(PETSC_HAVE_MPI_FINALIZED)
1703:     PetscMPIInt flag;
1704:     MPI_Finalized(&flag);
1705:     if (PetscUnlikely(flag)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
1706: #endif
1707:     MPI_Finalize();
1708:   }
1709: /*

1711:      Note: In certain cases PETSC_COMM_WORLD is never MPI_Comm_free()ed because
1712:    the communicator has some outstanding requests on it. Specifically if the
1713:    flag PETSC_HAVE_BROKEN_REQUEST_FREE is set (for IBM MPI implementation). See
1714:    src/vec/utils/vpscat.c. Due to this the memory allocated in PetscCommDuplicate()
1715:    is never freed as it should be. Thus one may obtain messages of the form
1716:    [ 1] 8 bytes PetscCommDuplicate() line 645 in src/sys/mpiu.c indicating the
1717:    memory was not freed.

1719: */
1720:   PetscMallocClear();
1721:   PetscStackReset();

1723:   PetscErrorHandlingInitialized = PETSC_FALSE;
1724:   PetscInitializeCalled = PETSC_FALSE;
1725:   PetscFinalizeCalled   = PETSC_TRUE;
1726: #if defined(PETSC_USE_GCOV)
1727:   /*
1728:      flush gcov, otherwise during CI the flushing continues into the next pipeline resulting in git not being able to delete directories since the
1729:      gcov files are still being added to the directories as git tries to remove the directories.
1730:    */
1731:   __gcov_flush();
1732: #endif
1734:   PetscStackClearTop;
1735:   return 0;
1736: }

1738: #if defined(PETSC_MISSING_LAPACK_lsame_)
1739: PETSC_EXTERN int lsame_(char *a,char *b)
1740: {
1741:   if (*a == *b) return 1;
1742:   if (*a + 32 == *b) return 1;
1743:   if (*a - 32 == *b) return 1;
1744:   return 0;
1745: }
1746: #endif

1748: #if defined(PETSC_MISSING_LAPACK_lsame)
1749: PETSC_EXTERN int lsame(char *a,char *b)
1750: {
1751:   if (*a == *b) return 1;
1752:   if (*a + 32 == *b) return 1;
1753:   if (*a - 32 == *b) return 1;
1754:   return 0;
1755: }
1756: #endif