Actual source code: petscimpl.h

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: /*
  3:     Defines the basic header of all PETSc objects.
  4: */

  6: #if !defined(PETSCIMPL_H)
  7: #define PETSCIMPL_H
  8:  #include <petscsys.h>

 10: /* These are used internally by PETSc ASCII IO routines*/
 11: #include <stdarg.h>
 12: PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);

 14: #if defined(PETSC_HAVE_CLOSURE)
 15: PETSC_EXTERN PetscErrorCode PetscVFPrintfSetClosure(int (^)(const char*));
 16: #endif

 18: /*
 19:    All major PETSc data structures have a common core; this is defined
 20:    below by PETSCHEADER.

 22:    PetscHeaderCreate() should be used whenever creating a PETSc structure.
 23: */

 25: /*
 26:    PetscOps: structure of core operations that all PETSc objects support.

 28:       getcomm()         - Gets the object's communicator.
 29:       view()            - Is the routine for viewing the entire PETSc object; for
 30:                           example, MatView() is the general matrix viewing routine.
 31:                           This is used by PetscObjectView((PetscObject)obj) to allow
 32:                           viewing any PETSc object.
 33:       destroy()         - Is the routine for destroying the entire PETSc object;
 34:                           for example,MatDestroy() is the general matrix
 35:                           destruction routine.
 36:                           This is used by PetscObjectDestroy((PetscObject*)&obj) to allow
 37:                           destroying any PETSc object.
 38:       compose()         - Associates a PETSc object with another PETSc object with a name
 39:       query()           - Returns a different PETSc object that has been associated
 40:                           with the first object using a name.
 41:       composefunction() - Attaches an a function to a PETSc object with a name.
 42:       queryfunction()   - Requests a registered function that has been attached to a PETSc object.
 43: */

 45: typedef struct {
 46:    PetscErrorCode (*getcomm)(PetscObject,MPI_Comm *);
 47:    PetscErrorCode (*view)(PetscObject,PetscViewer);
 48:    PetscErrorCode (*destroy)(PetscObject*);
 49:    PetscErrorCode (*compose)(PetscObject,const char[],PetscObject);
 50:    PetscErrorCode (*query)(PetscObject,const char[],PetscObject *);
 51:    PetscErrorCode (*composefunction)(PetscObject,const char[],void (*)(void));
 52:    PetscErrorCode (*queryfunction)(PetscObject,const char[],void (**)(void));
 53: } PetscOps;

 55: typedef enum {PETSC_FORTRAN_CALLBACK_CLASS,PETSC_FORTRAN_CALLBACK_SUBTYPE,PETSC_FORTRAN_CALLBACK_MAXTYPE} PetscFortranCallbackType;
 56: typedef int PetscFortranCallbackId;
 57: #define PETSC_SMALLEST_FORTRAN_CALLBACK ((PetscFortranCallbackId)1000)
 58: PETSC_EXTERN PetscErrorCode PetscFortranCallbackRegister(PetscClassId,const char*,PetscFortranCallbackId*);
 59: PETSC_EXTERN PetscErrorCode PetscFortranCallbackGetSizes(PetscClassId,PetscInt*,PetscInt*);

 61: typedef struct {
 62:   void (*func)(void);
 63:   void *ctx;
 64: } PetscFortranCallback;

 66: /*
 67:    All PETSc objects begin with the fields defined in PETSCHEADER.
 68:    The PetscObject is a way of examining these fields regardless of
 69:    the specific object. In C++ this could be a base abstract class
 70:    from which all objects are derived.
 71: */
 72: #define PETSC_MAX_OPTIONS_HANDLER 5
 73: typedef struct _p_PetscObject {
 74:   PetscClassId         classid;
 75:   PetscOps             bops[1];
 76:   MPI_Comm             comm;
 77:   PetscInt             type;
 78:   PetscLogDouble       flops,time,mem,memchildren;
 79:   PetscObjectId        id;
 80:   PetscInt             refct;
 81:   PetscMPIInt          tag;
 82:   PetscFunctionList    qlist;
 83:   PetscObjectList      olist;
 84:   char                 *class_name;    /*  for example, "Vec" */
 85:   char                 *description;
 86:   char                 *mansec;
 87:   char                 *type_name;     /*  this is the subclass, for example VECSEQ which equals "seq" */
 88:   PetscObject          parent;
 89:   PetscObjectId        parentid;
 90:   char*                name;
 91:   char                 *prefix;
 92:   PetscInt             tablevel;
 93:   void                 *cpp;
 94:   PetscObjectState     state;
 95:   PetscInt             int_idmax,        intstar_idmax;
 96:   PetscObjectState     *intcomposedstate,*intstarcomposedstate;
 97:   PetscInt             *intcomposeddata, **intstarcomposeddata;
 98:   PetscInt             real_idmax,        realstar_idmax;
 99:   PetscObjectState     *realcomposedstate,*realstarcomposedstate;
100:   PetscReal            *realcomposeddata, **realstarcomposeddata;
101:   PetscInt             scalar_idmax,        scalarstar_idmax;
102:   PetscObjectState     *scalarcomposedstate,*scalarstarcomposedstate;
103:   PetscScalar          *scalarcomposeddata, **scalarstarcomposeddata;
104:   void                 (**fortran_func_pointers)(void);                  /* used by Fortran interface functions to stash user provided Fortran functions */
105:   PetscInt             num_fortran_func_pointers;                        /* number of Fortran function pointers allocated */
106:   PetscFortranCallback *fortrancallback[PETSC_FORTRAN_CALLBACK_MAXTYPE];
107:   PetscInt             num_fortrancallback[PETSC_FORTRAN_CALLBACK_MAXTYPE];
108:   void                 *python_context;
109:   PetscErrorCode       (*python_destroy)(void*);

111:   PetscInt             noptionhandler;
112:   PetscErrorCode       (*optionhandler[PETSC_MAX_OPTIONS_HANDLER])(PetscOptionItems*,PetscObject,void*);
113:   PetscErrorCode       (*optiondestroy[PETSC_MAX_OPTIONS_HANDLER])(PetscObject,void*);
114:   void                 *optionctx[PETSC_MAX_OPTIONS_HANDLER];
115:   PetscBool            optionsprinted;
116: #if defined(PETSC_HAVE_SAWS)
117:   PetscBool            amsmem;          /* if PETSC_TRUE then this object is registered with SAWs and visible to clients */
118:   PetscBool            amspublishblock; /* if PETSC_TRUE and publishing objects then will block at PetscObjectSAWsBlock() */
119: #endif
120:   PetscOptions         options;         /* options database used, NULL means default */
121:   PetscBool            donotPetscObjectPrintClassNamePrefixType;
122: } _p_PetscObject;

124: #define PETSCHEADER(ObjectOps) \
125:   _p_PetscObject hdr;          \
126:   ObjectOps      ops[1]

128: #define  PETSCFREEDHEADER -1

130: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscObjectDestroyFunction)(PetscObject*); /* force cast in next macro to NEVER use extern "C" style */
131: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscObjectViewFunction)(PetscObject,PetscViewer);

133: /*@C
134:     PetscHeaderCreate - Creates a PETSc object of a particular class

136:     Input Parameters:
137: +   classid - the classid associated with this object (for example VEC_CLASSID)
138: .   class_name - string name of class; should be static (for example "Vec")
139: .   descr - string containing short description; should be static (for example "Vector")
140: .   mansec - string indicating section in manual pages; should be static (for example "Vec")
141: .   comm - the MPI Communicator
142: .   destroy - the destroy routine for this object (for example VecDestroy())
143: -   view - the view routine for this object (for example VecView())

145:     Output Parameter:
146: .   h - the newly created object

148:     Level: developer

150: .seealso: PetscHeaderDestroy(), PetscClassIdRegister()

152: @*/
153: #define PetscHeaderCreate(h,classid,class_name,descr,mansec,comm,destroy,view) \
154:   (PetscNew(&(h)) || \
155:    PetscHeaderCreate_Private((PetscObject)(h),classid,class_name,descr,mansec,comm,(PetscObjectDestroyFunction)(destroy),(PetscObjectViewFunction)(view)) || \
156:    PetscLogObjectCreate(h) || \
157:    PetscLogObjectMemory((PetscObject)(h),sizeof(*(h))))

159: PETSC_EXTERN PetscErrorCode PetscComposedQuantitiesDestroy(PetscObject obj);
160: PETSC_EXTERN PetscErrorCode PetscHeaderCreate_Private(PetscObject,PetscClassId,const char[],const char[],const char[],MPI_Comm,PetscObjectDestroyFunction,PetscObjectViewFunction);

162: /*@C
163:     PetscHeaderDestroy - Final step in destroying a PetscObject

165:     Input Parameters:
166: .   h - the header created with PetscHeaderCreate()

168:     Level: developer

170: .seealso: PetscHeaderCreate()
171: @*/
172: #define PetscHeaderDestroy(h) (PetscHeaderDestroy_Private((PetscObject)(*(h))) || PetscFree(*(h)))

174: PETSC_EXTERN PetscErrorCode PetscHeaderDestroy_Private(PetscObject);
175: PETSC_EXTERN PetscErrorCode PetscObjectCopyFortranFunctionPointers(PetscObject,PetscObject);
176: PETSC_EXTERN PetscErrorCode PetscObjectSetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId*,void(*)(void),void *ctx);
177: PETSC_EXTERN PetscErrorCode PetscObjectGetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId,void(**)(void),void **ctx);

179: PETSC_INTERN PetscErrorCode PetscCitationsInitialize(void);
180: PETSC_INTERN PetscErrorCode PetscFreeMPIResources(void);



185: #if defined(PETSC_HAVE_CUDA)
187: #endif
188: /*
189:     Macros to test if a PETSc object is valid and if pointers are valid
190: */
191: #if !defined(PETSC_USE_DEBUG)


204: #else

206: /*  This check is for subtype methods such as DMDAGetCorners() that do not use the PetscTryMethod() or PetscUseMethod() paradigm */
208:   do {   \
209:     PetscErrorCode _7_ierr; \
210:     PetscBool      _7_same; \
212:     _7_PetscObjectTypeCompare((PetscObject)(h),t,&_7_same);CHKERRQ(_7_ierr); \
213:     if (!_7_same) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong subtype object:Parameter # %d must have implementation %s it is %s",arg,t,((PetscObject)(h))->type_name); \
214:   } while (0)

217:   do {                                                                  \
218:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object: Parameter # %d",arg); \
220:     if (((PetscObject)(h))->classid != ck) {                            \
221:       if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free: Parameter # %d",arg); \
222:       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong type of object: Parameter # %d",arg); \
223:     }                                                                   \
224:   } while (0)

227:   do {                                                                  \
228:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object: Parameter # %d",arg); \
230:     if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free: Parameter # %d",arg); \
231:     else if (((PetscObject)(h))->classid < PETSC_SMALLEST_CLASSID || ((PetscObject)(h))->classid > PETSC_LARGEST_CLASSID) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Invalid type of object: Parameter # %d",arg); \
232:   } while (0)

235:   do {                                                                  \
236:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer: Parameter # %d",arg); \
238:   } while (0)

241:   do {                                                                  \
242:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer: Parameter # %d",arg);\
244:   } while (0)

247:   do {                                                                  \
248:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Null Pointer: Parameter # %d",arg); \
250:   } while (0)

253:   do {                                                                  \
254:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Null Pointer: Parameter # %d",arg); \
256:   } while (0)

259:   do {                                                                  \
260:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer: Parameter # %d",arg); \
262:   } while (0)

265:   do {                                                                  \
266:     if (!(h)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer: Parameter # %d",arg); \
268:   } while (0)

271:   do {                                                                  \
272:     if (!(f)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Function Pointer: Parameter # %d",arg); \
273:   } while (0)

275: #endif

277: #define PetscSorted(n,idx,sorted)           \
278:   do {                                      \
279:     PetscInt _i_;                           \
280:     (sorted) = PETSC_TRUE;                  \
281:     for (_i_ = 1; _i_ < (n); _i_++)         \
282:       if ((idx)[_i_] < (idx)[_i_ - 1])      \
283:         { (sorted) = PETSC_FALSE; break; }  \
284:   } while(0)

286: #if !defined(PETSC_USE_DEBUG)


302: #else

304: /*
305:     For example, in the dot product between two vectors,
306:   both vectors must be either Seq or MPI, not one of each
307: */
309:   do {                                                                  \
310:     if (((PetscObject)(a))->type != ((PetscObject)(b))->type) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Objects not of same type: Argument # %d and %d",arga,argb); \
311:   } while (0)
312: /*
313:     Check type_name
314: */
316:   do {                                                                  \
317:     PetscBool      _7_match;                                            \
318:     PetscErrorCode _7_ierr;                                             \
319:     _7_PetscObjectTypeCompare(((PetscObject)(a)),(type),&_7_match);CHKERRQ(_7_ierr); \
320:     if (!_7_match) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Object (%s) is not %s",(char*)(((PetscObject)(a))->type_name),type); \
321:   } while (0)

324:   do {                                                                  \
325:     PetscBool      _7_match;                                            \
326:     PetscErrorCode _7_ierr;                                             \
327:     _7_PetscObjectTypeCompareAny(((PetscObject)(a)),&_7_match,(type1),(type2),"");CHKERRQ(_7_ierr); \
328:     if (!_7_match) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Object (%s) is not %s or %s",(char*)(((PetscObject)(a))->type_name),type1,type2); \
329:   } while (0)
330: /*
331:    Use this macro to check if the type is set
332: */
334:   do {                                                                  \
335:     if (!((PetscObject)(a))->type_name) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"%s object's type is not set: Argument # %d",((PetscObject)(a))->class_name,arg); \
336:   } while (0)
337: /*
338:    Sometimes object must live on same communicator to inter-operate
339: */
341:   do {                                                                  \
342:     PetscErrorCode _7_ierr;                                             \
343:     PetscMPIInt    _7_flag;                                             \
344:     _7_MPI_Comm_compare(PetscObjectComm((PetscObject)(a)),PetscObjectComm((PetscObject)(b)),&_7_flag);CHKERRQ(_7_ierr); \
345:     if (_7_flag != MPI_CONGRUENT && _7_flag != MPI_IDENT) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"Different communicators in the two objects: Argument # %d and %d flag %d",arga,argb,_7_flag); \
346:   } while (0)

349:   do {                                                  \
352:   } while (0)

355:   do {                                                                  \
356:     PetscErrorCode _7_ierr;                                             \
357:     PetscScalar b0=(b);                                                 \
358:     PetscReal b1[5],b2[5];                                              \
359:     if (PetscIsNanScalar(b0)) {b1[4] = 1;} else {b1[4] = 0;};           \
360:     b1[0] = -PetscRealPart(b0); b1[1] = PetscRealPart(b0); b1[2] = -PetscImaginaryPart(b0); b1[3] = PetscImaginaryPart(b0); \
361:     _7_MPI_Allreduce(b1,b2,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)(a)));CHKERRQ(_7_ierr); \
362:     if (!(b2[4] > 0) && !(PetscEqualReal(-b2[0],b2[1]) && PetscEqualReal(-b2[2],b2[3]))) SETERRQ1(PetscObjectComm((PetscObject)(a)),PETSC_ERR_ARG_WRONG,"Scalar value must be same on all processes, argument # %d",arg); \
363:   } while (0)

366:   do {                                                                  \
367:     PetscErrorCode _7_ierr;                                             \
368:     PetscReal b0=(b),b1[3],b2[3];                                       \
369:     if (PetscIsNanReal(b0)) {b1[2] = 1;} else {b1[2] = 0;};             \
370:     b1[0] = -b0; b1[1] = b0;                                            \
371:     _7_MPI_Allreduce(b1,b2,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)(a)));CHKERRQ(_7_ierr); \
372:     if (!(b2[2] > 0) && !PetscEqualReal(-b2[0],b2[1])) SETERRQ1(PetscObjectComm((PetscObject)(a)),PETSC_ERR_ARG_WRONG,"Real value must be same on all processes, argument # %d",arg); \
373:   } while (0)

376:   do {                                                                  \
377:     PetscErrorCode _7_ierr;                                             \
378:     PetscInt b0=(b),b1[2],b2[2];                                        \
379:     b1[0] = -b0; b1[1] = b0;                                            \
380:     _7_MPIU_Allreduce(b1,b2,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)(a)));CHKERRQ(_7_ierr); \
381:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)(a)),PETSC_ERR_ARG_WRONG,"Int value must be same on all processes, argument # %d",arg); \
382:   } while (0)

385:   do {                                                                  \
386:     PetscErrorCode _7_ierr;                                             \
387:     PetscMPIInt b0=(b),b1[2],b2[2];                                     \
388:     b1[0] = -b0; b1[1] = b0;                                            \
389:     _7_MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)(a)));CHKERRQ(_7_ierr); \
390:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)(a)),PETSC_ERR_ARG_WRONG,"PetscMPIInt value must be same on all processes, argument # %d",arg); \
391:   } while (0)

394:   do {                                                                  \
395:     PetscErrorCode _7_ierr;                                             \
396:     PetscMPIInt b0=(PetscMPIInt)(b),b1[2],b2[2];                        \
397:     b1[0] = -b0; b1[1] = b0;                                            \
398:     _7_MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)(a)));CHKERRQ(_7_ierr); \
399:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)(a)),PETSC_ERR_ARG_WRONG,"Bool value must be same on all processes, argument # %d",arg); \
400:   } while (0)

403:   do {                                                                  \
404:     PetscErrorCode _7_ierr;                                             \
405:     PetscMPIInt b0=(PetscMPIInt)(b),b1[2],b2[2];                        \
406:     b1[0] = -b0; b1[1] = b0;                                            \
407:     _7_MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)(a)));CHKERRQ(_7_ierr); \
408:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)(a)),PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes, argument # %d",arg); \
409:   } while (0)

412:   do {                                                                                            \
413:     PetscBool _1_flg;                                                                             \
414:     PetscSorted(n,idx,_1_flg);                                                                    \
415:     if (!_1_flg) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input array needs to be sorted"); \
416:   } while (0)

418: #endif

420: /*
421:    PetscTryMethod - Queries an object for a method, if it exists then calls it.
422:               These are intended to be used only inside PETSc functions.

424:    Level: developer

426: .seealso: PetscUseMethod()
427: */
428: #define  PetscTryMethod(obj,A,B,C) \
429:   0; do { PetscErrorCode (*_7_f)B, _7_ierr; \
430:     _7_PetscObjectQueryFunction((PetscObject)(obj),A,&_7_f);CHKERRQ(_7_ierr); \
431:     if (_7_f) {_7_(*_7_f)C;CHKERRQ(_7_ierr);} \
432:   } while(0)

434: /*
435:    PetscUseMethod - Queries an object for a method, if it exists then calls it, otherwise generates an error.
436:               These are intended to be used only inside PETSc functions.

438:    Level: developer

440: .seealso: PetscTryMethod()
441: */
442: #define  PetscUseMethod(obj,A,B,C) \
443:   0; do { PetscErrorCode (*_7_f)B, _7_ierr; \
444:     _7_PetscObjectQueryFunction((PetscObject)(obj),A,&_7_f);CHKERRQ(_7_ierr); \
445:     if (_7_f) {_7_(*_7_f)C;CHKERRQ(_7_ierr);} \
446:     else SETERRQ1(PetscObjectComm((PetscObject)(obj)),PETSC_ERR_SUP,"Cannot locate function %s in object",A); \
447:   } while(0)

449: /*MC
450:    PetscObjectStateIncrease - Increases the state of any PetscObject

452:    Synopsis:
453:    #include "petsc/private/petscimpl.h"
454:    PetscErrorCode PetscObjectStateIncrease(PetscObject obj)

456:    Logically Collective

458:    Input Parameter:
459: .  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
460:          cast with a (PetscObject), for example,
461:          PetscObjectStateIncrease((PetscObject)mat);

463:    Notes:
464:     object state is an integer which gets increased every time
465:    the object is changed internally. By saving and later querying the object state
466:    one can determine whether information about the object is still current.
467:    Currently, state is maintained for Vec and Mat objects.

469:    This routine is mostly for internal use by PETSc; a developer need only
470:    call it after explicit access to an object's internals. Routines such
471:    as VecSet() or MatScale() already call this routine. It is also called, as a
472:    precaution, in VecRestoreArray(), MatRestoreRow(), MatDenseRestoreArray().

474:    This routine is logically collective because state equality comparison needs to be possible without communication.

476:    Level: developer

478:    seealso: PetscObjectStateGet()

480: M*/
481: #define PetscObjectStateIncrease(obj) ((obj)->state++,0)

483: PETSC_EXTERN PetscErrorCode PetscObjectStateGet(PetscObject,PetscObjectState*);
484: PETSC_EXTERN PetscErrorCode PetscObjectStateSet(PetscObject,PetscObjectState);
485: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataRegister(PetscInt*);
486: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseInt(PetscObject);
487: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseIntstar(PetscObject);
488: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseReal(PetscObject);
489: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseRealstar(PetscObject);
490: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalar(PetscObject);
491: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalarstar(PetscObject);
492: PETSC_EXTERN PetscInt       PetscObjectComposedDataMax;
493: /*MC
494:    PetscObjectComposedDataSetInt - attach integer data to a PetscObject

496:    Synopsis:
497:    #include "petsc/private/petscimpl.h"
498:    PetscErrorCode PetscObjectComposedDataSetInt(PetscObject obj,int id,int data)

500:    Not collective

502:    Input parameters:
503: +  obj - the object to which data is to be attached
504: .  id - the identifier for the data
505: -  data - the data to  be attached

507:    Notes
508:    The data identifier can best be created through a call to  PetscObjectComposedDataRegister()

510:    Level: developer
511: M*/
512: #define PetscObjectComposedDataSetInt(obj,id,data)                                      \
513:   ((((obj)->int_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseInt(obj)) ||  \
514:    ((obj)->intcomposeddata[id] = data,(obj)->intcomposedstate[id] = (obj)->state, 0))

516: /*MC
517:    PetscObjectComposedDataGetInt - retrieve integer data attached to an object

519:    Synopsis:
520:    #include "petsc/private/petscimpl.h"
521:    PetscErrorCode PetscObjectComposedDataGetInt(PetscObject obj,int id,int data,PetscBool  flag)

523:    Not collective

525:    Input parameters:
526: +  obj - the object from which data is to be retrieved
527: -  id - the identifier for the data

529:    Output parameters:
530: +  data - the data to be retrieved
531: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

533:    The 'data' and 'flag' variables are inlined, so they are not pointers.

535:    Level: developer
536: M*/
537: #define PetscObjectComposedDataGetInt(obj,id,data,flag)                            \
538:   ((((obj)->intcomposedstate && ((obj)->intcomposedstate[id] == (obj)->state)) ?   \
539:    (data = (obj)->intcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

541: /*MC
542:    PetscObjectComposedDataSetIntstar - attach integer array data to a PetscObject

544:    Synopsis:
545:    #include "petsc/private/petscimpl.h"
546:    PetscErrorCode PetscObjectComposedDataSetIntstar(PetscObject obj,int id,int *data)

548:    Not collective

550:    Input parameters:
551: +  obj - the object to which data is to be attached
552: .  id - the identifier for the data
553: -  data - the data to  be attached

555:    Notes
556:    The data identifier can best be determined through a call to
557:    PetscObjectComposedDataRegister()

559:    Level: developer
560: M*/
561: #define PetscObjectComposedDataSetIntstar(obj,id,data)                                          \
562:   ((((obj)->intstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseIntstar(obj)) ||  \
563:    ((obj)->intstarcomposeddata[id] = data,(obj)->intstarcomposedstate[id] = (obj)->state, 0))

565: /*MC
566:    PetscObjectComposedDataGetIntstar - retrieve integer array data
567:    attached to an object

569:    Synopsis:
570:    #include "petsc/private/petscimpl.h"
571:    PetscErrorCode PetscObjectComposedDataGetIntstar(PetscObject obj,int id,int *data,PetscBool  flag)

573:    Not collective

575:    Input parameters:
576: +  obj - the object from which data is to be retrieved
577: -  id - the identifier for the data

579:    Output parameters:
580: +  data - the data to be retrieved
581: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

583:    The 'data' and 'flag' variables are inlined, so they are not pointers.

585:    Level: developer
586: M*/
587: #define PetscObjectComposedDataGetIntstar(obj,id,data,flag)                               \
588:   ((((obj)->intstarcomposedstate && ((obj)->intstarcomposedstate[id] == (obj)->state)) ?  \
589:    (data = (obj)->intstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

591: /*MC
592:    PetscObjectComposedDataSetReal - attach real data to a PetscObject

594:    Synopsis:
595:    #include "petsc/private/petscimpl.h"
596:    PetscErrorCode PetscObjectComposedDataSetReal(PetscObject obj,int id,PetscReal data)

598:    Not collective

600:    Input parameters:
601: +  obj - the object to which data is to be attached
602: .  id - the identifier for the data
603: -  data - the data to  be attached

605:    Notes
606:    The data identifier can best be determined through a call to
607:    PetscObjectComposedDataRegister()

609:    Level: developer
610: M*/
611: #define PetscObjectComposedDataSetReal(obj,id,data)                                       \
612:   ((((obj)->real_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseReal(obj)) ||  \
613:    ((obj)->realcomposeddata[id] = data,(obj)->realcomposedstate[id] = (obj)->state, 0))

615: /*MC
616:    PetscObjectComposedDataGetReal - retrieve real data attached to an object

618:    Synopsis:
619:    #include "petsc/private/petscimpl.h"
620:    PetscErrorCode PetscObjectComposedDataGetReal(PetscObject obj,int id,PetscReal data,PetscBool  flag)

622:    Not collective

624:    Input parameters:
625: +  obj - the object from which data is to be retrieved
626: -  id - the identifier for the data

628:    Output parameters:
629: +  data - the data to be retrieved
630: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

632:    The 'data' and 'flag' variables are inlined, so they are not pointers.

634:    Level: developer
635: M*/
636: #define PetscObjectComposedDataGetReal(obj,id,data,flag)                            \
637:   ((((obj)->realcomposedstate && ((obj)->realcomposedstate[id] == (obj)->state)) ?  \
638:    (data = (obj)->realcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

640: /*MC
641:    PetscObjectComposedDataSetRealstar - attach real array data to a PetscObject

643:    Synopsis:
644:    #include "petsc/private/petscimpl.h"
645:    PetscErrorCode PetscObjectComposedDataSetRealstar(PetscObject obj,int id,PetscReal *data)

647:    Not collective

649:    Input parameters:
650: +  obj - the object to which data is to be attached
651: .  id - the identifier for the data
652: -  data - the data to  be attached

654:    Notes
655:    The data identifier can best be determined through a call to
656:    PetscObjectComposedDataRegister()

658:    Level: developer
659: M*/
660: #define PetscObjectComposedDataSetRealstar(obj,id,data)                                           \
661:   ((((obj)->realstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseRealstar(obj)) ||  \
662:    ((obj)->realstarcomposeddata[id] = data, (obj)->realstarcomposedstate[id] = (obj)->state, 0))

664: /*MC
665:    PetscObjectComposedDataGetRealstar - retrieve real array data
666:    attached to an object

668:    Synopsis:
669:    #include "petsc/private/petscimpl.h"
670:    PetscErrorCode PetscObjectComposedDataGetRealstar(PetscObject obj,int id,PetscReal *data,PetscBool  flag)

672:    Not collective

674:    Input parameters:
675: +  obj - the object from which data is to be retrieved
676: -  id - the identifier for the data

678:    Output parameters:
679: +  data - the data to be retrieved
680: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

682:    The 'data' and 'flag' variables are inlined, so they are not pointers.

684:    Level: developer
685: M*/
686: #define PetscObjectComposedDataGetRealstar(obj,id,data,flag)                                \
687:   ((((obj)->realstarcomposedstate && ((obj)->realstarcomposedstate[id] == (obj)->state)) ?  \
688:    (data = (obj)->realstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

690: /*MC
691:    PetscObjectComposedDataSetScalar - attach scalar data to a PetscObject

693:    Synopsis:
694:    #include "petsc/private/petscimpl.h"
695:    PetscErrorCode PetscObjectComposedDataSetScalar(PetscObject obj,int id,PetscScalar data)

697:    Not collective

699:    Input parameters:
700: +  obj - the object to which data is to be attached
701: .  id - the identifier for the data
702: -  data - the data to  be attached

704:    Notes
705:    The data identifier can best be determined through a call to
706:    PetscObjectComposedDataRegister()

708:    Level: developer
709: M*/
710: #if defined(PETSC_USE_COMPLEX)
711: #define PetscObjectComposedDataSetScalar(obj,id,data)                                        \
712:   ((((obj)->scalar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalar(obj)) || \
713:    ((obj)->scalarcomposeddata[id] = data,(obj)->scalarcomposedstate[id] = (obj)->state, 0))
714: #else
715: #define PetscObjectComposedDataSetScalar(obj,id,data) \
716:         PetscObjectComposedDataSetReal(obj,id,data)
717: #endif
718: /*MC
719:    PetscObjectComposedDataGetScalar - retrieve scalar data attached to an object

721:    Synopsis:
722:    #include "petsc/private/petscimpl.h"
723:    PetscErrorCode PetscObjectComposedDataGetScalar(PetscObject obj,int id,PetscScalar data,PetscBool  flag)

725:    Not collective

727:    Input parameters:
728: +  obj - the object from which data is to be retrieved
729: -  id - the identifier for the data

731:    Output parameters:
732: +  data - the data to be retrieved
733: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

735:    The 'data' and 'flag' variables are inlined, so they are not pointers.

737:    Level: developer
738: M*/
739: #if defined(PETSC_USE_COMPLEX)
740: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                              \
741:   ((((obj)->scalarcomposedstate && ((obj)->scalarcomposedstate[id] == (obj)->state) ) ? \
742:    (data = (obj)->scalarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
743: #else
744: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                             \
745:         PetscObjectComposedDataGetReal(obj,id,data,flag)
746: #endif

748: /*MC
749:    PetscObjectComposedDataSetScalarstar - attach scalar array data to a PetscObject

751:    Synopsis:
752:    #include "petsc/private/petscimpl.h"
753:    PetscErrorCode PetscObjectComposedDataSetScalarstar(PetscObject obj,int id,PetscScalar *data)

755:    Not collective

757:    Input parameters:
758: +  obj - the object to which data is to be attached
759: .  id - the identifier for the data
760: -  data - the data to  be attached

762:    Notes
763:    The data identifier can best be determined through a call to
764:    PetscObjectComposedDataRegister()

766:    Level: developer
767: M*/
768: #if defined(PETSC_USE_COMPLEX)
769: #define PetscObjectComposedDataSetScalarstar(obj,id,data)                                             \
770:   ((((obj)->scalarstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalarstar(obj)) ||  \
771:    ((obj)->scalarstarcomposeddata[id] = data,(obj)->scalarstarcomposedstate[id] = (obj)->state, 0))
772: #else
773: #define PetscObjectComposedDataSetScalarstar(obj,id,data) \
774:         PetscObjectComposedDataSetRealstar(obj,id,data)
775: #endif
776: /*MC
777:    PetscObjectComposedDataGetScalarstar - retrieve scalar array data
778:    attached to an object

780:    Synopsis:
781:    #include "petsc/private/petscimpl.h"
782:    PetscErrorCode PetscObjectComposedDataGetScalarstar(PetscObject obj,int id,PetscScalar *data,PetscBool  flag)

784:    Not collective

786:    Input parameters:
787: +  obj - the object from which data is to be retrieved
788: -  id - the identifier for the data

790:    Output parameters:
791: +  data - the data to be retrieved
792: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

794:    The 'data' and 'flag' variables are inlined, so they are not pointers.

796:    Level: developer
797: M*/
798: #if defined(PETSC_USE_COMPLEX)
799: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)                                 \
800:   ((((obj)->scalarstarcomposedstate && ((obj)->scalarstarcomposedstate[id] == (obj)->state)) ? \
801:        (data = (obj)->scalarstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
802: #else
803: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)         \
804:         PetscObjectComposedDataGetRealstar(obj,id,data,flag)
805: #endif

807: PETSC_EXTERN PetscMPIInt Petsc_Counter_keyval;
808: PETSC_EXTERN PetscMPIInt Petsc_InnerComm_keyval;
809: PETSC_EXTERN PetscMPIInt Petsc_OuterComm_keyval;
810: PETSC_EXTERN PetscMPIInt Petsc_Seq_keyval;
811: PETSC_EXTERN PetscMPIInt Petsc_ShmComm_keyval;

813: /*
814:   PETSc communicators have this attribute, see
815:   PetscCommDuplicate(), PetscCommDestroy(), PetscCommGetNewTag(), PetscObjectGetName()
816: */
817: typedef struct {
818:   PetscMPIInt tag;              /* next free tag value */
819:   PetscInt    refcount;         /* number of references, communicator can be freed when this reaches 0 */
820:   PetscInt    namecount;        /* used to generate the next name, as in Vec_0, Mat_1, ... */
821:   PetscMPIInt *iflags;          /* length of comm size, shared by all calls to PetscCommBuildTwoSided_Allreduce/RedScatter on this comm */
822: } PetscCommCounter;

824: /*E
825:     PetscOffloadMask - indicates which memory (CPU, GPU, or none) contains valid data

827:    PETSC_OFFLOAD_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
828:    PETSC_OFFLOAD_GPU - GPU has valid vector/matrix entries
829:    PETSC_OFFLOAD_CPU - CPU has valid vector/matrix entries
830:    PETSC_OFFLOAD_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

832:    Level: developer
833: E*/
834: typedef enum {PETSC_OFFLOAD_UNALLOCATED=0x0,PETSC_OFFLOAD_CPU=0x1,PETSC_OFFLOAD_GPU=0x2,PETSC_OFFLOAD_BOTH=0x3} PetscOffloadMask;

836: typedef enum {STATE_BEGIN, STATE_PENDING, STATE_END} SRState;

838: typedef enum {PETSC_SR_REDUCE_SUM=0,PETSC_SR_REDUCE_MAX=1,PETSC_SR_REDUCE_MIN=2} PetscSRReductionType;

840: typedef struct {
841:   MPI_Comm    comm;
842:   MPI_Request request;
843:   PetscBool   async;
844:   PetscScalar *lvalues;     /* this are the reduced values before call to MPI_Allreduce() */
845:   PetscScalar *gvalues;     /* values after call to MPI_Allreduce() */
846:   void        **invecs;     /* for debugging only, vector/memory used with each op */
847:   PetscInt    *reducetype;  /* is particular value to be summed or maxed? */
848:   SRState     state;        /* are we calling xxxBegin() or xxxEnd()? */
849:   PetscInt    maxops;       /* total amount of space we have for requests */
850:   PetscInt    numopsbegin;  /* number of requests that have been queued in */
851:   PetscInt    numopsend;    /* number of requests that have been gotten by user */
852: } PetscSplitReduction;

854: PETSC_EXTERN PetscErrorCode PetscSplitReductionGet(MPI_Comm,PetscSplitReduction**);
855: PETSC_EXTERN PetscErrorCode PetscSplitReductionEnd(PetscSplitReduction*);
856: PETSC_EXTERN PetscErrorCode PetscSplitReductionExtend(PetscSplitReduction*);

858: #if !defined(PETSC_SKIP_SPINLOCK)
859: #if defined(PETSC_HAVE_THREADSAFETY)
860: #  if defined(PETSC_HAVE_CONCURRENCYKIT)
861: #if defined(__cplusplus)
862: /*  CK does not have extern "C" protection in their include files */
863: extern "C" {
864: #endif
865: #include <ck_spinlock.h>
866: #if defined(__cplusplus)
867: }
868: #endif
869: typedef ck_spinlock_t PetscSpinlock;
870: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockCreate(PetscSpinlock *ck_spinlock)
871: {
872:   ck_spinlock_init(ck_spinlock);
873:   return 0;
874: }
875: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockLock(PetscSpinlock *ck_spinlock)
876: {
877:   ck_spinlock_lock(ck_spinlock);
878:   return 0;
879: }
880: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockUnlock(PetscSpinlock *ck_spinlock)
881: {
882:   ck_spinlock_unlock(ck_spinlock);
883:   return 0;
884: }
885: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockDestroy(PetscSpinlock *ck_spinlock)
886: {
887:   return 0;
888: }
889: #  elif defined(PETSC_HAVE_OPENMP)

891: #include <omp.h>
892: typedef omp_lock_t PetscSpinlock;
893: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockCreate(PetscSpinlock *omp_lock)
894: {
895:   omp_init_lock(omp_lock);
896:   return 0;
897: }
898: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockLock(PetscSpinlock *omp_lock)
899: {
900:   omp_set_lock(omp_lock);
901:   return 0;
902: }
903: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockUnlock(PetscSpinlock *omp_lock)
904: {
905:   omp_unset_lock(omp_lock);
906:   return 0;
907: }
908: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockDestroy(PetscSpinlock *omp_lock)
909: {
910:   omp_destroy_lock(omp_lock);
911:   return 0;
912: }
913: #else
914: Thread safety requires either --with-openmp or --download-concurrencykit
915: #endif

917: #else
918: typedef int PetscSpinlock;
919: #define PetscSpinlockCreate(a)  0
920: #define PetscSpinlockLock(a)    0
921: #define PetscSpinlockUnlock(a)  0
922: #define PetscSpinlockDestroy(a) 0
923: #endif

925: #if defined(PETSC_HAVE_THREADSAFETY)
926: PETSC_INTERN PetscSpinlock PetscViewerASCIISpinLockOpen;
927: PETSC_INTERN PetscSpinlock PetscViewerASCIISpinLockStdout;
928: PETSC_INTERN PetscSpinlock PetscViewerASCIISpinLockStderr;
929: PETSC_INTERN PetscSpinlock PetscCommSpinLock;
930: #endif
931: #endif

933: PETSC_EXTERN PetscLogEvent PETSC_Barrier;
934: PETSC_EXTERN PetscLogEvent PETSC_BuildTwoSided;
935: PETSC_EXTERN PetscLogEvent PETSC_BuildTwoSidedF;
936: PETSC_EXTERN PetscBool     use_gpu_aware_mpi;

938: #if defined(PETSC_HAVE_ADIOS)
939: PETSC_EXTERN int64_t Petsc_adios_group;
940: #endif

942: #endif /* PETSCIMPL_H */