Actual source code: pinit.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: /*
  3:    This file defines the initialization of PETSc, including PetscInitialize()
  4: */
  5: #define PETSC_DESIRE_COMPLEX
  6: #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
  7: #include <petscviewer.h>

  9: #if defined(PETSC_HAVE_CUDA)
 10: #include <cublas.h>
 11: #endif

 13: #include <petscthreadcomm.h>

 15: #if defined(PETSC_USE_LOG)
 16: extern PetscErrorCode PetscLogBegin_Private(void);
 17: #endif

 19: #if defined(PETSC_SERIALIZE_FUNCTIONS)
 20: PetscFPT PetscFPTData = 0;
 21: #endif

 23: #if defined(PETSC_HAVE_SAWS)
 24: #include <petscviewersaws.h>
 25: #endif
 26: /* -----------------------------------------------------------------------------------------*/

 28: extern FILE *petsc_history;

 30: extern PetscErrorCode PetscInitialize_DynamicLibraries(void);
 31: extern PetscErrorCode PetscFinalize_DynamicLibraries(void);
 32: extern PetscErrorCode PetscFunctionListPrintAll(void);
 33: extern PetscErrorCode PetscSequentialPhaseBegin_Private(MPI_Comm,int);
 34: extern PetscErrorCode PetscSequentialPhaseEnd_Private(MPI_Comm,int);
 35: extern PetscErrorCode PetscCloseHistoryFile(FILE**);

 37: /* user may set this BEFORE calling PetscInitialize() */
 38: MPI_Comm PETSC_COMM_WORLD = MPI_COMM_NULL;

 40: PetscMPIInt Petsc_Counter_keyval   = MPI_KEYVAL_INVALID;
 41: PetscMPIInt Petsc_InnerComm_keyval = MPI_KEYVAL_INVALID;
 42: PetscMPIInt Petsc_OuterComm_keyval = MPI_KEYVAL_INVALID;

 44: /*
 45:      Declare and set all the string names of the PETSc enums
 46: */
 47: const char *const PetscBools[]     = {"FALSE","TRUE","PetscBool","PETSC_",0};
 48: const char *const PetscCopyModes[] = {"COPY_VALUES","OWN_POINTER","USE_POINTER","PetscCopyMode","PETSC_",0};
 49: const char *const PetscDataTypes[] = {"INT","DOUBLE","COMPLEX","LONG","SHORT","FLOAT",
 50:                                       "CHAR","LOGICAL","ENUM","BOOL","LONGDOUBLE","OBJECT","FUNCTION","PetscDataType","PETSC_",0};

 52: PetscBool PetscPreLoadingUsed = PETSC_FALSE;
 53: PetscBool PetscPreLoadingOn   = PETSC_FALSE;

 55: PetscInt PetscHotRegionDepth;

 57: /*
 58:        Checks the options database for initializations related to the
 59:     PETSc components
 60: */
 63: PetscErrorCode  PetscOptionsCheckInitial_Components(void)
 64: {
 65:   PetscBool      flg1;

 69:   PetscOptionsHasName(NULL,"-help",&flg1);
 70:   if (flg1) {
 71: #if defined(PETSC_USE_LOG)
 72:     MPI_Comm comm = PETSC_COMM_WORLD;
 73:     (*PetscHelpPrintf)(comm,"------Additional PETSc component options--------\n");
 74:     (*PetscHelpPrintf)(comm," -log_summary_exclude: <vec,mat,pc.ksp,snes>\n");
 75:     (*PetscHelpPrintf)(comm," -info_exclude: <null,vec,mat,pc,ksp,snes,ts>\n");
 76:     (*PetscHelpPrintf)(comm,"-----------------------------------------------\n");
 77: #endif
 78:   }
 79:   return(0);
 80: }

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

 87:    Collective

 89:    Level: advanced

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

 95:      Turns off PETSc signal handling because that can interact with MATLAB's signal handling causing random crashes.

 97: .seealso: PetscInitialize(), PetscInitializeFortran(), PetscInitializeNoArguments()
 98: */
 99: PetscErrorCode  PetscInitializeNoPointers(int argc,char **args,const char *filename,const char *help)
100: {
102:   int            myargc   = argc;
103:   char           **myargs = args;

106:   PetscInitialize(&myargc,&myargs,filename,help);
107:   PetscPopSignalHandler();
108:   PetscBeganMPI = PETSC_FALSE;
109:   PetscFunctionReturn(ierr);
110: }

114: /*
115:       Used by MATLAB and Julia interface to get communicator
116: */
117: PetscErrorCode  PetscGetPETSC_COMM_SELF(MPI_Comm *comm)
118: {
120:   *comm = PETSC_COMM_SELF;
121:   return(0);
122: }

126: /*@C
127:       PetscInitializeNoArguments - Calls PetscInitialize() from C/C++ without
128:         the command line arguments.

130:    Collective

132:    Level: advanced

134: .seealso: PetscInitialize(), PetscInitializeFortran()
135: @*/
136: PetscErrorCode  PetscInitializeNoArguments(void)
137: {
139:   int            argc   = 0;
140:   char           **args = 0;

143:   PetscInitialize(&argc,&args,NULL,NULL);
144:   PetscFunctionReturn(ierr);
145: }

149: /*@
150:       PetscInitialized - Determine whether PETSc is initialized.

152:    Level: beginner

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

164: /*@
165:       PetscFinalized - Determine whether PetscFinalize() has been called yet

167:    Level: developer

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

177: extern PetscErrorCode PetscOptionsCheckInitial_Private(void);

179: /*
180:        This function is the MPI reduction operation used to compute the sum of the
181:    first half of the datatype and the max of the second half.
182: */
183: MPI_Op PetscMaxSum_Op = 0;

187: PETSC_EXTERN void MPIAPI PetscMaxSum_Local(void *in,void *out,int *cnt,MPI_Datatype *datatype)
188: {
189:   PetscInt *xin = (PetscInt*)in,*xout = (PetscInt*)out,i,count = *cnt;

192:   if (*datatype != MPIU_2INT) {
193:     (*PetscErrorPrintf)("Can only handle MPIU_2INT data types");
194:     MPI_Abort(MPI_COMM_WORLD,1);
195:   }

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

204: /*
205:     Returns the max of the first entry owned by this processor and the
206: sum of the second entry.

208:     The reason nprocs[2*i] contains lengths nprocs[2*i+1] contains flag of 1 if length is nonzero
209: is so that the PetscMaxSum_Op() can set TWO values, if we passed in only nprocs[i] with lengths
210: there would be no place to store the both needed results.
211: */
214: PetscErrorCode  PetscMaxSum(MPI_Comm comm,const PetscInt sizes[],PetscInt *max,PetscInt *sum)
215: {
216:   PetscMPIInt    size,rank;
217:   struct {PetscInt max,sum;} *work;

221:   MPI_Comm_size(comm,&size);
222:   MPI_Comm_rank(comm,&rank);
223:   PetscMalloc1(size,&work);
224:   MPI_Allreduce((void*)sizes,work,size,MPIU_2INT,PetscMaxSum_Op,comm);
225:   *max = work[rank].max;
226:   *sum = work[rank].sum;
227:   PetscFree(work);
228:   return(0);
229: }

231: /* ----------------------------------------------------------------------------*/
232: MPI_Op  PetscADMax_Op = 0;

236: PETSC_EXTERN void MPIAPI PetscADMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
237: {
238:   PetscScalar *xin = (PetscScalar*)in,*xout = (PetscScalar*)out;
239:   PetscInt    i,count = *cnt;

242:   if (*datatype != MPIU_2SCALAR) {
243:     (*PetscErrorPrintf)("Can only handle MPIU_2SCALAR data (i.e. double or complex) types");
244:     MPI_Abort(MPI_COMM_WORLD,1);
245:   }

247:   for (i=0; i<count; i++) {
248:     if (PetscRealPart(xout[2*i]) < PetscRealPart(xin[2*i])) {
249:       xout[2*i]   = xin[2*i];
250:       xout[2*i+1] = xin[2*i+1];
251:     }
252:   }
253:   PetscFunctionReturnVoid();
254: }

256: MPI_Op PetscADMin_Op = 0;

260: PETSC_EXTERN void MPIAPI PetscADMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
261: {
262:   PetscScalar *xin = (PetscScalar*)in,*xout = (PetscScalar*)out;
263:   PetscInt    i,count = *cnt;

266:   if (*datatype != MPIU_2SCALAR) {
267:     (*PetscErrorPrintf)("Can only handle MPIU_2SCALAR data (i.e. double or complex) types");
268:     MPI_Abort(MPI_COMM_WORLD,1);
269:   }

271:   for (i=0; i<count; i++) {
272:     if (PetscRealPart(xout[2*i]) > PetscRealPart(xin[2*i])) {
273:       xout[2*i]   = xin[2*i];
274:       xout[2*i+1] = xin[2*i+1];
275:     }
276:   }
277:   PetscFunctionReturnVoid();
278: }
279: /* ---------------------------------------------------------------------------------------*/

281: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
282: MPI_Op MPIU_SUM = 0;

286: PETSC_EXTERN void PetscSum_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
287: {
288:   PetscInt i,count = *cnt;

291:   if (*datatype == MPIU_REAL) {
292:     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
293:     for (i=0; i<count; i++) xout[i] += xin[i];
294:   }
295: #if defined(PETSC_HAVE_COMPLEX)
296:   else if (*datatype == MPIU_COMPLEX) {
297:     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
298:     for (i=0; i<count; i++) xout[i] += xin[i];
299:   }
300: #endif
301:   else {
302:     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
303:     MPI_Abort(MPI_COMM_WORLD,1);
304:   }
305:   PetscFunctionReturnVoid();
306: }
307: #endif

309: #if defined(PETSC_USE_REAL___FLOAT128)
310: MPI_Op MPIU_MAX = 0;
311: MPI_Op MPIU_MIN = 0;

315: PETSC_EXTERN void PetscMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
316: {
317:   PetscInt i,count = *cnt;

320:   if (*datatype == MPIU_REAL) {
321:     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
322:     for (i=0; i<count; i++) xout[i] = PetscMax(xout[i],xin[i]);
323:   }
324: #if defined(PETSC_HAVE_COMPLEX)
325:   else if (*datatype == MPIU_COMPLEX) {
326:     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
327:     for (i=0; i<count; i++) {
328:       xout[i] = PetscRealPartComplex(xout[i])<PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
329:     }
330:   }
331: #endif
332:   else {
333:     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_COMPLEX data types");
334:     MPI_Abort(MPI_COMM_WORLD,1);
335:   }
336:   PetscFunctionReturnVoid();
337: }

341: PETSC_EXTERN void PetscMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
342: {
343:   PetscInt    i,count = *cnt;

346:   if (*datatype == MPIU_REAL) {
347:     PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
348:     for (i=0; i<count; i++) xout[i] = PetscMin(xout[i],xin[i]);
349:   }
350: #if defined(PETSC_HAVE_COMPLEX)
351:   else if (*datatype == MPIU_COMPLEX) {
352:     PetscComplex *xin = (PetscComplex*)in,*xout = (PetscComplex*)out;
353:     for (i=0; i<count; i++) {
354:       xout[i] = PetscRealPartComplex(xout[i])>PetscRealPartComplex(xin[i]) ? xin[i] : xout[i];
355:     }
356:   }
357: #endif
358:   else {
359:     (*PetscErrorPrintf)("Can only handle MPIU_REAL or MPIU_SCALAR data (i.e. double or complex) types");
360:     MPI_Abort(MPI_COMM_WORLD,1);
361:   }
362:   PetscFunctionReturnVoid();
363: }
364: #endif

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

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

373:    Note: this is declared extern "C" because it is passed to MPI_Keyval_create()

375: */
376: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelCounter(MPI_Comm comm,PetscMPIInt keyval,void *count_val,void *extra_state)
377: {

381:   PetscInfo1(0,"Deleting counter data in an MPI_Comm %ld\n",(long)comm);if (ierr) PetscFunctionReturn((PetscMPIInt)ierr);
382:   PetscFree(count_val);if (ierr) PetscFunctionReturn((PetscMPIInt)ierr);
383:   PetscFunctionReturn(MPI_SUCCESS);
384: }

388: /*
389:   This is invoked on the outer comm as a result of either PetscCommDestroy() (via MPI_Attr_delete) or when the user
390:   calls MPI_Comm_free().

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

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

396:   Note: this is declared extern "C" because it is passed to MPI_Keyval_create()

398: */
399: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Outer(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
400: {
402:   PetscMPIInt    flg;
403:   union {MPI_Comm comm; void *ptr;} icomm,ocomm;

406:   if (keyval != Petsc_InnerComm_keyval) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Unexpected keyval");
407:   icomm.ptr = attr_val;

409:   MPI_Attr_get(icomm.comm,Petsc_OuterComm_keyval,&ocomm,&flg);
410:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected reference to outer comm");
411:   if (ocomm.comm != comm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm has reference to non-matching outer comm");
412:   MPI_Attr_delete(icomm.comm,Petsc_OuterComm_keyval); /* Calls Petsc_DelComm_Inner */
413:   PetscInfo1(0,"User MPI_Comm %ld is being freed after removing reference from inner PETSc comm to this outer comm\n",(long)comm);if (ierr) PetscFunctionReturn((PetscMPIInt)ierr);
414:   PetscFunctionReturn(MPI_SUCCESS);
415: }

419: /*
420:  * This is invoked on the inner comm when Petsc_DelComm_Outer calls MPI_Attr_delete.  It should not be reached any other way.
421:  */
422: PETSC_EXTERN PetscMPIInt MPIAPI Petsc_DelComm_Inner(MPI_Comm comm,PetscMPIInt keyval,void *attr_val,void *extra_state)
423: {

427:   PetscInfo1(0,"Removing reference to PETSc communicator embedded in a user MPI_Comm %ld\n",(long)comm);if (ierr) PetscFunctionReturn((PetscMPIInt)ierr);
428:   PetscFunctionReturn(MPI_SUCCESS);
429: }

431: #if defined(PETSC_USE_PETSC_MPI_EXTERNAL32)
432: #if !defined(PETSC_WORDS_BIGENDIAN)
433: PETSC_EXTERN PetscMPIInt PetscDataRep_extent_fn(MPI_Datatype,MPI_Aint*,void*);
434: PETSC_EXTERN PetscMPIInt PetscDataRep_read_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
435: PETSC_EXTERN PetscMPIInt PetscDataRep_write_conv_fn(void*, MPI_Datatype,PetscMPIInt,void*,MPI_Offset,void*);
436: #endif
437: #endif

439: int  PetscGlobalArgc   = 0;
440: char **PetscGlobalArgs = 0;
441: PetscSegBuffer PetscCitationsList;

445: PetscErrorCode PetscCitationsInitialize()
446: {

450:   PetscSegBufferCreate(1,10000,&PetscCitationsList);
451:   PetscCitationsRegister("@TechReport{petsc-user-ref,\n  Author = {Satish Balay and Shrirang Abhyankar and Mark F. Adams and Jed Brown and Peter Brune\n            and Kris Buschelman and Victor Eijkhout and William D. Gropp\n            and Dinesh Kaushik and Matthew G. Knepley\n            and Lois Curfman McInnes and Karl Rupp and Barry F. Smith\n            and Hong Zhang},\n  Title = {{PETS}c Users Manual},\n  Number = {ANL-95/11 - Revision 3.5},\n  Institution = {Argonne National Laboratory},\n  Year = {2014}\n}\n",NULL);
452:   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);
453:   return(0);
454: }

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

462:    Not Collective

464:    Output Parameters:
465: +  argc - count of number of command line arguments
466: -  args - the command line arguments

468:    Level: intermediate

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

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

476:    Concepts: command line arguments

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

480: @*/
481: PetscErrorCode  PetscGetArgs(int *argc,char ***args)
482: {
484:   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
485:   *argc = PetscGlobalArgc;
486:   *args = PetscGlobalArgs;
487:   return(0);
488: }

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

496:    Not Collective

498:    Output Parameters:
499: .  args - the command line arguments

501:    Level: intermediate

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

506:    Concepts: command line arguments

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

510: @*/
511: PetscErrorCode  PetscGetArguments(char ***args)
512: {
513:   PetscInt       i,argc = PetscGlobalArgc;

517:   if (!PetscInitializeCalled && PetscFinalizeCalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call after PetscInitialize() but before PetscFinalize()");
518:   if (!argc) {*args = 0; return(0);}
519:   PetscMalloc1(argc,args);
520:   for (i=0; i<argc-1; i++) {
521:     PetscStrallocpy(PetscGlobalArgs[i+1],&(*args)[i]);
522:   }
523:   (*args)[argc-1] = 0;
524:   return(0);
525: }

529: /*@C
530:    PetscFreeArguments - Frees the memory obtained with PetscGetArguments()

532:    Not Collective

534:    Output Parameters:
535: .  args - the command line arguments

537:    Level: intermediate

539:    Concepts: command line arguments

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

543: @*/
544: PetscErrorCode  PetscFreeArguments(char **args)
545: {
546:   PetscInt       i = 0;

550:   if (!args) return(0);
551:   while (args[i]) {
552:     PetscFree(args[i]);
553:     i++;
554:   }
555:   PetscFree(args);
556:   return(0);
557: }

559: #if defined(PETSC_HAVE_SAWS)
560: #include <petscconfiginfo.h>

564: PetscErrorCode  PetscInitializeSAWs(const char help[])
565: {
566:   if (!PetscGlobalRank) {
567:     char           cert[PETSC_MAX_PATH_LEN],root[PETSC_MAX_PATH_LEN],*intro,programname[64],*appline,*options,version[64];
568:     int            port;
569:     PetscBool      flg,rootlocal = PETSC_FALSE,flg2;
570:     size_t         applinelen,introlen;

573:     PetscOptionsHasName(NULL,"-saws_log",&flg);
574:     if (flg) {
575:       char  sawslog[PETSC_MAX_PATH_LEN];

577:       PetscOptionsGetString(NULL,"-saws_log",sawslog,PETSC_MAX_PATH_LEN,NULL);
578:       if (sawslog[0]) {
579:         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(sawslog));
580:       } else {
581:         PetscStackCallSAWs(SAWs_Set_Use_Logfile,(NULL));
582:       }
583:     }
584:     PetscOptionsGetString(NULL,"-saws_https",cert,PETSC_MAX_PATH_LEN,&flg);
585:     if (flg) {
586:       PetscStackCallSAWs(SAWs_Set_Use_HTTPS,(cert));
587:     }
588:     PetscOptionsGetInt(NULL,"-saws_port",&port,&flg);
589:     if (flg) {
590:       PetscStackCallSAWs(SAWs_Set_Port,(port));
591:     }
592:     PetscOptionsGetString(NULL,"-saws_root",root,PETSC_MAX_PATH_LEN,&flg);
593:     if (flg) {
594:       PetscStackCallSAWs(SAWs_Set_Document_Root,(root));
595:       PetscStrcmp(root,".",&rootlocal);
596:     }
597:     PetscOptionsHasName(NULL,"-saws_local",&flg2);
598:     if (flg2) {
599:       char jsdir[PETSC_MAX_PATH_LEN];
600:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"-saws_local option requires -saws_root option");
601:       PetscSNPrintf(jsdir,PETSC_MAX_PATH_LEN,"%s/js",root);
602:       PetscTestDirectory(jsdir,'r',&flg);
603:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"-saws_local option requires js directory in root directory");
604:       PetscStackCallSAWs(SAWs_Push_Local_Header,());
605:     }
606:     PetscGetProgramName(programname,64);
607:     PetscStrlen(help,&applinelen);
608:     introlen   = 4096 + applinelen;
609:     applinelen += 256;
610:     PetscMalloc(applinelen,&appline);
611:     PetscMalloc(introlen,&intro);

613:     if (rootlocal) {
614:       PetscSNPrintf(appline,applinelen,"%s.c.html",programname);
615:       PetscTestFile(appline,'r',&rootlocal);
616:     }
617:     PetscOptionsGetAll(&options);
618:     if (rootlocal && help) {
619:       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);
620:     } else if (help) {
621:       PetscSNPrintf(appline,applinelen,"<center>Running %s %s</center><br><center><pre>%s</pre></center><br>\n",programname,options,help);
622:     } else {
623:       PetscSNPrintf(appline,applinelen,"<center> Running %s %s</center><br>\n",programname,options);
624:     }
625:     PetscFree(options);
626:     PetscGetVersion(version,sizeof(version));
627:     PetscSNPrintf(intro,introlen,"<body>\n"
628:                                     "<center><h2> <a href=\"http://www.mcs.anl.gov/petsc\">PETSc</a> Application Web server powered by <a href=\"https://bitbucket.org/saws/saws\">SAWs</a> </h2></center>\n"
629:                                     "<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"
630:                                     "%s",version,petscconfigureoptions,appline);
631:     PetscStackCallSAWs(SAWs_Push_Body,("index.html",0,intro));
632:     PetscFree(intro);
633:     PetscFree(appline);
634:     PetscStackCallSAWs(SAWs_Initialize,());
635:     PetscCitationsRegister("@TechReport{ saws,"
636:                                   "Author = {Matt Otten and Jed Brown and Barry Smith},"
637:                                   "Title  = {Scientific Application Web Server (SAWs) Users Manual},"
638:                                   "Institution = {Argonne National Laboratory},"
639:                                   "Year   = 2013}",NULL);
640:   }
641:   return(0);
642: }
643: #endif

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

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

655:    Input Parameters:
656: +  argc - count of number of command line arguments
657: .  args - the command line arguments
658: .  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use NULL to not check for
659:           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files
660: -  help - [optional] Help message to print, use NULL for no message

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

668:    Options Database Keys:
669: +  -start_in_debugger [noxterm,dbx,xdb,gdb,...] - Starts program in debugger
670: .  -on_error_attach_debugger [noxterm,dbx,xdb,gdb,...] - Starts debugger when error detected
671: .  -on_error_emacs <machinename> causes emacsclient to jump to error file
672: .  -on_error_abort calls abort() when error detected (no traceback)
673: .  -on_error_mpiabort calls MPI_abort() when error detected
674: .  -error_output_stderr prints error messages to stderr instead of the default stdout
675: .  -error_output_none does not print the error messages (but handles errors in the same way as if this was not called)
676: .  -debugger_nodes [node1,node2,...] - Indicates nodes to start in debugger
677: .  -debugger_pause [sleeptime] (in seconds) - Pauses debugger
678: .  -stop_for_debugger - Print message on how to attach debugger manually to
679:                         process and wait (-debugger_pause) seconds for attachment
680: .  -malloc - Indicates use of PETSc error-checking malloc (on by default for debug version of libraries)
681: .  -malloc no - Indicates not to use error-checking malloc
682: .  -malloc_debug - check for memory corruption at EVERY malloc or free
683: .  -malloc_dump - prints a list of all unfreed memory at the end of the run
684: .  -malloc_test - like -malloc_dump -malloc_debug, but only active for debugging builds
685: .  -fp_trap - Stops on floating point exceptions (Note that on the
686:               IBM RS6000 this slows code by at least a factor of 10.)
687: .  -no_signal_handler - Indicates not to trap error signals
688: .  -shared_tmp - indicates /tmp directory is shared by all processors
689: .  -not_shared_tmp - each processor has own /tmp
690: .  -tmp - alternative name of /tmp directory
691: .  -get_total_flops - returns total flops done by all processors
692: -  -memory_info - Print memory usage at end of run

694:    Options Database Keys for Profiling:
695:    See Users-Manual: Chapter 12 Profiling for details.
696: +  -info <optional filename> - Prints verbose information to the screen
697: .  -info_exclude <null,vec,mat,pc,ksp,snes,ts> - Excludes some of the verbose messages
698: .  -log_sync - Log the synchronization in scatters, inner products and norms
699: .  -log_trace [filename] - Print traces of all PETSc calls to the screen (useful to determine where a program
700:         hangs without running in the debugger).  See PetscLogTraceBegin().
701: .  -log_summary [filename] - Prints summary of flop and timing information to screen. If the filename is specified the
702:         summary is written to the file.  See PetscLogView().
703: .  -log_all [filename] - Logs extensive profiling information  See PetscLogDump().
704: .  -log [filename] - Logs basic profiline information  See PetscLogDump().
705: -  -log_mpe [filename] - Creates a logfile viewable by the utility Jumpshot (in MPICH distribution)

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

709:    Environmental Variables:
710: +   PETSC_TMP - alternative tmp directory
711: .   PETSC_SHARED_TMP - tmp is shared by all processes
712: .   PETSC_NOT_SHARED_TMP - each process has its own private tmp
713: .   PETSC_VIEWER_SOCKET_PORT - socket number to use for socket viewer
714: -   PETSC_VIEWER_SOCKET_MACHINE - machine to use for socket viewer to connect to


717:    Level: beginner

719:    Notes:
720:    If for some reason you must call MPI_Init() separately, call
721:    it before PetscInitialize().

723:    Fortran Version:
724:    In Fortran this routine has the format
725: $       call PetscInitialize(file,ierr)

727: +   ierr - error return code
728: -  file - [optional] PETSc database file, also checks ~username/.petscrc and .petscrc use PETSC_NULL_CHARACTER to not check for
729:           code specific file. Use -skip_petscrc in the code specific file to skip the .petscrc files

731:    Important Fortran Note:
732:    In Fortran, you MUST use PETSC_NULL_CHARACTER to indicate a
733:    null character string; you CANNOT just use NULL as
734:    in the C version. See Users-Manual: Chapter 11 PETSc for Fortran Users for details.

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

739:    Concepts: initializing PETSc

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

743: @*/
744: PetscErrorCode  PetscInitialize(int *argc,char ***args,const char file[],const char help[])
745: {
747:   PetscMPIInt    flag, size;
748:   PetscBool      flg;
749:   char           hostname[256];

752:   if (PetscInitializeCalled) return(0);

754:   /* these must be initialized in a routine, not as a constant declaration*/
755:   PETSC_STDOUT = stdout;
756:   PETSC_STDERR = stderr;

758:   /* on Windows - set printf to default to printing 2 digit exponents */
759: #if defined(PETSC_HAVE__SET_OUTPUT_FORMAT)
760:   _set_output_format(_TWO_DIGIT_EXPONENT);
761: #endif

763:   PetscOptionsCreate();

765:   /*
766:      We initialize the program name here (before MPI_Init()) because MPICH has a bug in
767:      it that it sets args[0] on all processors to be args[0] on the first processor.
768:   */
769:   if (argc && *argc) {
770:     PetscSetProgramName(**args);
771:   } else {
772:     PetscSetProgramName("Unknown Name");
773:   }

775:   MPI_Initialized(&flag);
776:   if (!flag) {
777:     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");
778: #if defined(PETSC_HAVE_MPI_INIT_THREAD)
779:     {
780:       PetscMPIInt provided;
781:       MPI_Init_thread(argc,args,MPI_THREAD_FUNNELED,&provided);
782:     }
783: #else
784:     MPI_Init(argc,args);
785: #endif
786:     PetscBeganMPI = PETSC_TRUE;
787:   }
788:   if (argc && args) {
789:     PetscGlobalArgc = *argc;
790:     PetscGlobalArgs = *args;
791:   }
792:   PetscFinalizeCalled = PETSC_FALSE;

794:   if (PETSC_COMM_WORLD == MPI_COMM_NULL) PETSC_COMM_WORLD = MPI_COMM_WORLD;
795:   MPI_Comm_set_errhandler(PETSC_COMM_WORLD,MPI_ERRORS_RETURN);

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

800:   MPI_Comm_rank(MPI_COMM_WORLD,&PetscGlobalRank);
801:   MPI_Comm_size(MPI_COMM_WORLD,&PetscGlobalSize);

803:   MPIU_BOOL = MPI_INT;
804:   MPIU_ENUM = MPI_INT;

806:   /*
807:      Initialized the global complex variable; this is because with
808:      shared libraries the constructors for global variables
809:      are not called; at least on IRIX.
810:   */
811: #if defined(PETSC_HAVE_COMPLEX)
812:   {
813: #if defined(PETSC_CLANGUAGE_CXX)
814:     PetscComplex ic(0.0,1.0);
815:     PETSC_i = ic;
816: #elif defined(PETSC_CLANGUAGE_C)
817:     PETSC_i = _Complex_I;
818: #endif
819:   }

821: #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
822:   MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU_C_DOUBLE_COMPLEX);
823:   MPI_Type_commit(&MPIU_C_DOUBLE_COMPLEX);
824:   MPI_Type_contiguous(2,MPI_FLOAT,&MPIU_C_COMPLEX);
825:   MPI_Type_commit(&MPIU_C_COMPLEX);
826: #endif
827: #endif /* PETSC_HAVE_COMPLEX */

829:   /*
830:      Create the PETSc MPI reduction operator that sums of the first
831:      half of the entries and maxes the second half.
832:   */
833:   MPI_Op_create(PetscMaxSum_Local,1,&PetscMaxSum_Op);

835: #if defined(PETSC_USE_REAL___FLOAT128)
836:   MPI_Type_contiguous(2,MPI_DOUBLE,&MPIU___FLOAT128);
837:   MPI_Type_commit(&MPIU___FLOAT128);
838: #if defined(PETSC_HAVE_COMPLEX)
839:   MPI_Type_contiguous(4,MPI_DOUBLE,&MPIU___COMPLEX128);
840:   MPI_Type_commit(&MPIU___COMPLEX128);
841: #endif
842:   MPI_Op_create(PetscMax_Local,1,&MPIU_MAX);
843:   MPI_Op_create(PetscMin_Local,1,&MPIU_MIN);
844: #endif

846: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
847:   MPI_Op_create(PetscSum_Local,1,&MPIU_SUM);
848: #endif

850:   MPI_Type_contiguous(2,MPIU_SCALAR,&MPIU_2SCALAR);
851:   MPI_Type_commit(&MPIU_2SCALAR);
852:   MPI_Op_create(PetscADMax_Local,1,&PetscADMax_Op);
853:   MPI_Op_create(PetscADMin_Local,1,&PetscADMin_Op);

855: #if defined(PETSC_USE_64BIT_INDICES) || !defined(MPI_2INT)
856:   MPI_Type_contiguous(2,MPIU_INT,&MPIU_2INT);
857:   MPI_Type_commit(&MPIU_2INT);
858: #endif


861:   /*
862:      Attributes to be set on PETSc communicators
863:   */
864:   MPI_Keyval_create(MPI_NULL_COPY_FN,Petsc_DelCounter,&Petsc_Counter_keyval,(void*)0);
865:   MPI_Keyval_create(MPI_NULL_COPY_FN,Petsc_DelComm_Outer,&Petsc_InnerComm_keyval,(void*)0);
866:   MPI_Keyval_create(MPI_NULL_COPY_FN,Petsc_DelComm_Inner,&Petsc_OuterComm_keyval,(void*)0);

868:   /*
869:      Build the options database
870:   */
871:   PetscOptionsInsert(argc,args,file);


874:   /*
875:      Print main application help message
876:   */
877:   PetscOptionsHasName(NULL,"-help",&flg);
878:   if (help && flg) {
879:     PetscPrintf(PETSC_COMM_WORLD,help);
880:   }
881:   PetscOptionsCheckInitial_Private();

883:   PetscCitationsInitialize();

885: #if defined(PETSC_HAVE_SAWS)
886:   PetscInitializeSAWs(help);
887: #endif

889:   /* SHOULD PUT IN GUARDS: Make sure logging is initialized, even if we do not print it out */
890: #if defined(PETSC_USE_LOG)
891:   PetscLogBegin_Private();
892: #endif

894:   /*
895:      Load the dynamic libraries (on machines that support them), this registers all
896:      the solvers etc. (On non-dynamic machines this initializes the PetscDraw and PetscViewer classes)
897:   */
898:   PetscInitialize_DynamicLibraries();

900:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
901:   PetscInfo1(0,"PETSc successfully started: number of processors = %d\n",size);
902:   PetscGetHostName(hostname,256);
903:   PetscInfo1(0,"Running on machine: %s\n",hostname);

905:   /* Ensure that threadcomm-related keyval exists, so that PetscOptionsSetFromOptions can use PetscCommDuplicate. */
906:   PetscThreadCommInitializePackage();

908:   PetscOptionsCheckInitial_Components();
909:   /* Check the options database for options related to the options database itself */
910:   PetscOptionsSetFromOptions();

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

916:       Currently not used because it is not supported by MPICH.
917:   */
918: #if !defined(PETSC_WORDS_BIGENDIAN)
919:   MPI_Register_datarep((char*)"petsc",PetscDataRep_read_conv_fn,PetscDataRep_write_conv_fn,PetscDataRep_extent_fn,NULL);
920: #endif
921: #endif

923: #if defined(PETSC_HAVE_CUDA)
924:   flg  = PETSC_TRUE;
925:   PetscOptionsGetBool(NULL,"-cublas",&flg,NULL);
926:   if (flg) {
927:     PetscMPIInt p;
928:     for (p = 0; p < PetscGlobalSize; ++p) {
929:       if (p == PetscGlobalRank) cublasInit();
930:       MPI_Barrier(PETSC_COMM_WORLD);
931:     }
932:   }
933: #endif

935:   PetscOptionsHasName(NULL,"-python",&flg);
936:   if (flg) {
937:     PetscInitializeCalled = PETSC_TRUE;
938:     PetscPythonInitialize(NULL,NULL);
939:   }

941:   /*
942:       Setup building of stack frames for all function calls
943:   */
944:   PetscThreadLocalRegister((PetscThreadKey*)&petscstack); /* Creates pthread_key */
945: #if defined(PETSC_USE_DEBUG)
946:   PetscStackCreate();
947: #endif

949: #if defined(PETSC_SERIALIZE_FUNCTIONS)
950:   PetscFPTCreate(10000);
951: #endif


954:   /*
955:       Once we are completedly initialized then we can set this variables
956:   */
957:   PetscInitializeCalled = PETSC_TRUE;
958:   return(0);
959: }

961: #if defined(PETSC_USE_LOG)
962: extern PetscObject *PetscObjects;
963: extern PetscInt    PetscObjectsCounts, PetscObjectsMaxCounts;
964: extern PetscBool   PetscObjectsLog;
965: #endif

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

974:    Collective on PETSC_COMM_WORLD

976:    Options Database Keys:
977: +  -options_table - Calls PetscOptionsView()
978: .  -options_left - Prints unused options that remain in the database
979: .  -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
980: .  -mpidump - Calls PetscMPIDump()
981: .  -malloc_dump - Calls PetscMallocDump()
982: .  -malloc_info - Prints total memory usage
983: -  -malloc_log - Prints summary of memory usage

985:    Level: beginner

987:    Note:
988:    See PetscInitialize() for more general runtime options.

990: .seealso: PetscInitialize(), PetscOptionsView(), PetscMallocDump(), PetscMPIDump(), PetscEnd()
991: @*/
992: PetscErrorCode  PetscFinalize(void)
993: {
995:   PetscMPIInt    rank;
996:   PetscInt       nopt;
997:   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE;
998:   PetscBool      flg;
999: #if defined(PETSC_USE_LOG)
1000:   char           mname[PETSC_MAX_PATH_LEN];
1001: #endif

1004:   if (!PetscInitializeCalled) {
1005:     printf("PetscInitialize() must be called before PetscFinalize()\n");
1006:     PetscFunctionReturn(PETSC_ERR_ARG_WRONGSTATE);
1007:   }
1008:   PetscInfo(NULL,"PetscFinalize() called\n");

1010:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

1012:   PetscOptionsHasName(NULL,"-citations",&flg);
1013:   if (flg) {
1014:     char  *cits, filename[PETSC_MAX_PATH_LEN];
1015:     FILE  *fd = PETSC_STDOUT;

1017:     PetscOptionsGetString(NULL,"-citations",filename,PETSC_MAX_PATH_LEN,NULL);
1018:     if (filename[0]) {
1019:       PetscFOpen(PETSC_COMM_WORLD,filename,"w",&fd);
1020:     }
1021:     PetscSegBufferGet(PetscCitationsList,1,&cits);
1022:     cits[0] = 0;
1023:     PetscSegBufferExtractAlloc(PetscCitationsList,&cits);
1024:     PetscFPrintf(PETSC_COMM_WORLD,fd,"If you publish results based on this computation please cite the following:\n");
1025:     PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");
1026:     PetscFPrintf(PETSC_COMM_WORLD,fd,"%s",cits);
1027:     PetscFPrintf(PETSC_COMM_WORLD,fd,"===========================================================================\n");
1028:     PetscFClose(PETSC_COMM_WORLD,fd);
1029:     PetscFree(cits);
1030:   }
1031:   PetscSegBufferDestroy(&PetscCitationsList);

1033: #if defined(PETSC_HAVE_SSL) && defined(PETSC_USE_SOCKET_VIEWER)
1034:   /* TextBelt is run for testing purposes only, please do not use this feature often */
1035:   {
1036:     PetscInt nmax = 2;
1037:     char     **buffs;
1038:     PetscMalloc1(2,&buffs);
1039:     PetscOptionsGetStringArray(NULL,"-textbelt",buffs,&nmax,&flg1);
1040:     if (flg1) {
1041:       if (!nmax) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-textbelt requires either the phone number or number,\"message\"");
1042:       if (nmax == 1) {
1043:         PetscMalloc1(128,&buffs[1]);
1044:         PetscGetProgramName(buffs[1],32);
1045:         PetscStrcat(buffs[1]," has completed");
1046:       }
1047:       PetscTextBelt(PETSC_COMM_WORLD,buffs[0],buffs[1],NULL);
1048:       PetscFree(buffs[0]);
1049:       PetscFree(buffs[1]);
1050:     }
1051:     PetscFree(buffs);
1052:   }
1053: #endif
1054:   /*
1055:     It should be safe to cancel the options monitors, since we don't expect to be setting options
1056:     here (at least that are worth monitoring).  Monitors ought to be released so that they release
1057:     whatever memory was allocated there before -malloc_dump reports unfreed memory.
1058:   */
1059:   PetscOptionsMonitorCancel();

1061: #if defined(PETSC_SERIALIZE_FUNCTIONS)
1062:   PetscFPTDestroy();
1063: #endif


1066: #if defined(PETSC_HAVE_SAWS)
1067:   flg = PETSC_FALSE;
1068:   PetscOptionsGetBool(NULL,"-saw_options",&flg,NULL);
1069:   if (flg) {
1070:     PetscOptionsSAWsDestroy();
1071:   }
1072: #endif

1074: #if defined(PETSC_HAVE_X)
1075:   flg1 = PETSC_FALSE;
1076:   PetscOptionsGetBool(NULL,"-x_virtual",&flg1,NULL);
1077:   if (flg1) {
1078:     /*  this is a crude hack, but better than nothing */
1079:     PetscPOpen(PETSC_COMM_WORLD,NULL,"pkill -9 Xvfb","r",NULL);
1080:   }
1081: #endif

1083:   PetscOptionsGetBool(NULL,"-malloc_info",&flg2,NULL);
1084:   if (!flg2) {
1085:     flg2 = PETSC_FALSE;
1086:     PetscOptionsGetBool(NULL,"-memory_info",&flg2,NULL);
1087:   }
1088:   if (flg2) {
1089:     PetscMemoryShowUsage(PETSC_VIEWER_STDOUT_WORLD,"Summary of Memory Usage in PETSc\n");
1090:   }

1092: #if defined(PETSC_USE_LOG)
1093:   flg1 = PETSC_FALSE;
1094:   PetscOptionsGetBool(NULL,"-get_total_flops",&flg1,NULL);
1095:   if (flg1) {
1096:     PetscLogDouble flops = 0;
1097:     MPI_Reduce(&petsc_TotalFlops,&flops,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);
1098:     PetscPrintf(PETSC_COMM_WORLD,"Total flops over all processors %g\n",flops);
1099:   }
1100: #endif


1103: #if defined(PETSC_USE_LOG)
1104: #if defined(PETSC_HAVE_MPE)
1105:   mname[0] = 0;

1107:   PetscOptionsGetString(NULL,"-log_mpe",mname,PETSC_MAX_PATH_LEN,&flg1);
1108:   if (flg1) {
1109:     if (mname[0]) {PetscLogMPEDump(mname);}
1110:     else          {PetscLogMPEDump(0);}
1111:   }
1112: #endif
1113:   mname[0] = 0;

1115:   PetscLogViewFromOptions();
1116:   PetscOptionsGetString(NULL,"-log_summary",mname,PETSC_MAX_PATH_LEN,&flg1);
1117:   if (flg1) {
1118:     PetscViewer viewer;
1119:     if (mname[0]) {
1120:       PetscViewerASCIIOpen(PETSC_COMM_WORLD,mname,&viewer);
1121:       PetscLogView(viewer);
1122:       PetscViewerDestroy(&viewer);
1123:     } else {
1124:       viewer = PETSC_VIEWER_STDOUT_WORLD;
1125:       PetscLogView(viewer);
1126:     }
1127:   }
1128:   mname[0] = 0;

1130:   PetscOptionsGetString(NULL,"-log_all",mname,PETSC_MAX_PATH_LEN,&flg1);
1131:   PetscOptionsGetString(NULL,"-log",mname,PETSC_MAX_PATH_LEN,&flg2);
1132:   if (flg1 || flg2) {
1133:     if (mname[0]) PetscLogDump(mname);
1134:     else          PetscLogDump(0);
1135:   }
1136: #endif

1138:   /*
1139:      Free all objects registered with PetscObjectRegisterDestroy() such as PETSC_VIEWER_XXX_().
1140:   */
1141:   PetscObjectRegisterDestroyAll();

1143:   PetscStackDestroy();
1144:   PetscThreadLocalDestroy((PetscThreadKey)petscstack); /* Deletes pthread_key */

1146:   flg1 = PETSC_FALSE;
1147:   PetscOptionsGetBool(NULL,"-no_signal_handler",&flg1,NULL);
1148:   if (!flg1) { PetscPopSignalHandler();}
1149:   flg1 = PETSC_FALSE;
1150:   PetscOptionsGetBool(NULL,"-mpidump",&flg1,NULL);
1151:   if (flg1) {
1152:     PetscMPIDump(stdout);
1153:   }
1154:   flg1 = PETSC_FALSE;
1155:   flg2 = PETSC_FALSE;
1156:   /* preemptive call to avoid listing this option in options table as unused */
1157:   PetscOptionsHasName(NULL,"-malloc_dump",&flg1);
1158:   PetscOptionsHasName(NULL,"-objects_dump",&flg1);
1159:   PetscOptionsGetBool(NULL,"-options_view",&flg2,NULL);

1161:   if (flg2) {
1162:     PetscViewer viewer;
1163:     PetscViewerCreate(PETSC_COMM_WORLD,&viewer);
1164:     PetscViewerSetType(viewer,PETSCVIEWERASCII);
1165:     PetscOptionsView(viewer);
1166:     PetscViewerDestroy(&viewer);
1167:   }

1169:   /* to prevent PETSc -options_left from warning */
1170:   PetscOptionsHasName(NULL,"-nox",&flg1);
1171:   PetscOptionsHasName(NULL,"-nox_warning",&flg1);

1173:   flg3 = PETSC_FALSE; /* default value is required */
1174:   PetscOptionsGetBool(NULL,"-options_left",&flg3,&flg1);
1175:   PetscOptionsAllUsed(&nopt);
1176:   if (flg3) {
1177:     if (!flg2) { /* have not yet printed the options */
1178:       PetscViewer viewer;
1179:       PetscViewerCreate(PETSC_COMM_WORLD,&viewer);
1180:       PetscViewerSetType(viewer,PETSCVIEWERASCII);
1181:       PetscOptionsView(viewer);
1182:       PetscViewerDestroy(&viewer);
1183:     }
1184:     if (!nopt) {
1185:       PetscPrintf(PETSC_COMM_WORLD,"There are no unused options.\n");
1186:     } else if (nopt == 1) {
1187:       PetscPrintf(PETSC_COMM_WORLD,"There is one unused database option. It is:\n");
1188:     } else {
1189:       PetscPrintf(PETSC_COMM_WORLD,"There are %D unused database options. They are:\n",nopt);
1190:     }
1191:   }
1192: #if defined(PETSC_USE_DEBUG)
1193:   if (nopt && !flg3 && !flg1) {
1194:     PetscPrintf(PETSC_COMM_WORLD,"WARNING! There are options you set that were not used!\n");
1195:     PetscPrintf(PETSC_COMM_WORLD,"WARNING! could be spelling mistake, etc!\n");
1196:     PetscOptionsLeft();
1197:   } else if (nopt && flg3) {
1198: #else
1199:   if (nopt && flg3) {
1200: #endif
1201:     PetscOptionsLeft();
1202:   }

1204: #if defined(PETSC_HAVE_SAWS)
1205:   if (!PetscGlobalRank) {
1206:     PetscStackSAWsViewOff();
1207:     PetscStackCallSAWs(SAWs_Finalize,());
1208:   }
1209: #endif

1211:   {
1212:     PetscThreadComm tcomm_world;
1213:     PetscGetThreadCommWorld(&tcomm_world);
1214:     /* Free global thread communicator */
1215:     PetscThreadCommDestroy(&tcomm_world);
1216:   }

1218: #if defined(PETSC_USE_LOG)
1219:   /*
1220:        List all objects the user may have forgot to free
1221:   */
1222:   if (PetscObjectsLog) {
1223:     PetscOptionsHasName(NULL,"-objects_dump",&flg1);
1224:     if (flg1) {
1225:       MPI_Comm local_comm;
1226:       char     string[64];

1228:       PetscOptionsGetString(NULL,"-objects_dump",string,64,NULL);
1229:       MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);
1230:       PetscSequentialPhaseBegin_Private(local_comm,1);
1231:       PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);
1232:       PetscSequentialPhaseEnd_Private(local_comm,1);
1233:       MPI_Comm_free(&local_comm);
1234:     }
1235:   }
1236: #endif

1238: #if defined(PETSC_USE_LOG)
1239:   PetscObjectsCounts    = 0;
1240:   PetscObjectsMaxCounts = 0;
1241:   PetscFree(PetscObjects);
1242: #endif

1244: #if defined(PETSC_USE_LOG)
1245:   PetscLogDestroy();
1246: #endif

1248:   /*
1249:      Destroy any packages that registered a finalize
1250:   */
1251:   PetscRegisterFinalizeAll();

1253:   /*
1254:      Destroy all the function registration lists created
1255:   */
1256:   PetscFinalize_DynamicLibraries();

1258:   /*
1259:      Print PetscFunctionLists that have not been properly freed

1261:   PetscFunctionListPrintAll();
1262:   */

1264:   if (petsc_history) {
1265:     PetscCloseHistoryFile(&petsc_history);
1266:     petsc_history = 0;
1267:   }

1269:   PetscInfoAllow(PETSC_FALSE,NULL);

1271:   {
1272:     char fname[PETSC_MAX_PATH_LEN];
1273:     FILE *fd;
1274:     int  err;

1276:     fname[0] = 0;

1278:     PetscOptionsGetString(NULL,"-malloc_dump",fname,250,&flg1);
1279:     flg2 = PETSC_FALSE;
1280:     PetscOptionsGetBool(NULL,"-malloc_test",&flg2,NULL);
1281: #if defined(PETSC_USE_DEBUG)
1282:     if (PETSC_RUNNING_ON_VALGRIND) flg2 = PETSC_FALSE;
1283: #else
1284:     flg2 = PETSC_FALSE;         /* Skip reporting for optimized builds regardless of -malloc_test */
1285: #endif
1286:     if (flg1 && fname[0]) {
1287:       char sname[PETSC_MAX_PATH_LEN];

1289:       sprintf(sname,"%s_%d",fname,rank);
1290:       fd   = fopen(sname,"w"); if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",sname);
1291:       PetscMallocDump(fd);
1292:       err  = fclose(fd);
1293:       if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1294:     } else if (flg1 || flg2) {
1295:       MPI_Comm local_comm;

1297:       MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);
1298:       PetscSequentialPhaseBegin_Private(local_comm,1);
1299:       PetscMallocDump(stdout);
1300:       PetscSequentialPhaseEnd_Private(local_comm,1);
1301:       MPI_Comm_free(&local_comm);
1302:     }
1303:   }

1305:   {
1306:     char fname[PETSC_MAX_PATH_LEN];
1307:     FILE *fd = NULL;

1309:     fname[0] = 0;

1311:     PetscOptionsGetString(NULL,"-malloc_log",fname,250,&flg1);
1312:     PetscOptionsHasName(NULL,"-malloc_log_threshold",&flg2);
1313:     if (flg1 && fname[0]) {
1314:       int err;

1316:       if (!rank) {
1317:         fd = fopen(fname,"w");
1318:         if (!fd) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open log file: %s",fname);
1319:       }
1320:       PetscMallocDumpLog(fd);
1321:       if (fd) {
1322:         err = fclose(fd);
1323:         if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
1324:       }
1325:     } else if (flg1 || flg2) {
1326:       PetscMallocDumpLog(stdout);
1327:     }
1328:   }

1330: #if defined(PETSC_HAVE_CUDA)
1331:   flg  = PETSC_TRUE;
1332:   PetscOptionsGetBool(NULL,"-cublas",&flg,NULL);
1333:   if (flg) {
1334:     PetscInt p;
1335:     for (p = 0; p < PetscGlobalSize; ++p) {
1336:       if (p == PetscGlobalRank) cublasShutdown();
1337:       MPI_Barrier(PETSC_COMM_WORLD);
1338:     }
1339:   }
1340: #endif

1342:   /* Can be destroyed only after all the options are used */
1343:   PetscOptionsDestroy();

1345:   PetscGlobalArgc = 0;
1346:   PetscGlobalArgs = 0;

1348: #if defined(PETSC_USE_REAL___FLOAT128)
1349:   MPI_Type_free(&MPIU___FLOAT128);
1350: #if defined(PETSC_HAVE_COMPLEX)
1351:   MPI_Type_free(&MPIU___COMPLEX128);
1352: #endif
1353:   MPI_Op_free(&MPIU_MAX);
1354:   MPI_Op_free(&MPIU_MIN);
1355: #endif

1357: #if defined(PETSC_HAVE_COMPLEX)
1358: #if !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
1359:   MPI_Type_free(&MPIU_C_DOUBLE_COMPLEX);
1360:   MPI_Type_free(&MPIU_C_COMPLEX);
1361: #endif
1362: #endif

1364: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1365:   MPI_Op_free(&MPIU_SUM);
1366: #endif

1368:   MPI_Type_free(&MPIU_2SCALAR);
1369: #if defined(PETSC_USE_64BIT_INDICES) || !defined(MPI_2INT)
1370:   MPI_Type_free(&MPIU_2INT);
1371: #endif
1372:   MPI_Op_free(&PetscMaxSum_Op);
1373:   MPI_Op_free(&PetscADMax_Op);
1374:   MPI_Op_free(&PetscADMin_Op);

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

1380:      If all PETSc objects were not destroyed those left over objects will have hanging references to
1381:      the MPI_Comms that were freed; but that is ok because those PETSc objects will never be used again
1382:  */
1383:   {
1384:     PetscCommCounter *counter;
1385:     PetscMPIInt      flg;
1386:     MPI_Comm         icomm;
1387:     union {MPI_Comm comm; void *ptr;} ucomm;
1388:     MPI_Attr_get(PETSC_COMM_SELF,Petsc_InnerComm_keyval,&ucomm,&flg);
1389:     if (flg) {
1390:       icomm = ucomm.comm;
1391:       MPI_Attr_get(icomm,Petsc_Counter_keyval,&counter,&flg);
1392:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");

1394:       MPI_Attr_delete(PETSC_COMM_SELF,Petsc_InnerComm_keyval);
1395:       MPI_Attr_delete(icomm,Petsc_Counter_keyval);
1396:       MPI_Comm_free(&icomm);
1397:     }
1398:     MPI_Attr_get(PETSC_COMM_WORLD,Petsc_InnerComm_keyval,&ucomm,&flg);
1399:     if (flg) {
1400:       icomm = ucomm.comm;
1401:       MPI_Attr_get(icomm,Petsc_Counter_keyval,&counter,&flg);
1402:       if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_CORRUPT,"Inner MPI_Comm does not have expected tag/name counter, problem with corrupted memory");

1404:       MPI_Attr_delete(PETSC_COMM_WORLD,Petsc_InnerComm_keyval);
1405:       MPI_Attr_delete(icomm,Petsc_Counter_keyval);
1406:       MPI_Comm_free(&icomm);
1407:     }
1408:   }

1410:   MPI_Keyval_free(&Petsc_Counter_keyval);
1411:   MPI_Keyval_free(&Petsc_InnerComm_keyval);
1412:   MPI_Keyval_free(&Petsc_OuterComm_keyval);

1414:   if (PetscBeganMPI) {
1415: #if defined(PETSC_HAVE_MPI_FINALIZED)
1416:     PetscMPIInt flag;
1417:     MPI_Finalized(&flag);
1418:     if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI_Finalize() has already been called, even though MPI_Init() was called by PetscInitialize()");
1419: #endif
1420:     MPI_Finalize();
1421:   }
1422: /*

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

1432: */
1433:   PetscMallocClear();

1435:   PetscInitializeCalled = PETSC_FALSE;
1436:   PetscFinalizeCalled   = PETSC_TRUE;
1437:   PetscFunctionReturn(ierr);
1438: }

1440: #if defined(PETSC_MISSING_LAPACK_lsame_)
1441: PETSC_EXTERN int lsame_(char *a,char *b)
1442: {
1443:   if (*a == *b) return 1;
1444:   if (*a + 32 == *b) return 1;
1445:   if (*a - 32 == *b) return 1;
1446:   return 0;
1447: }
1448: #endif

1450: #if defined(PETSC_MISSING_LAPACK_lsame)
1451: PETSC_EXTERN int lsame(char *a,char *b)
1452: {
1453:   if (*a == *b) return 1;
1454:   if (*a + 32 == *b) return 1;
1455:   if (*a - 32 == *b) return 1;
1456:   return 0;
1457: }
1458: #endif