Actual source code: petscimpl.h

petsc-3.9.4 2018-09-11
Report Typos and Errors

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

  6: #if !defined(_PETSCHEAD_H)
  7: #define _PETSCHEAD_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_MATLAB_ENGINE)
 15: PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
 16: #endif

 18: #if defined(PETSC_HAVE_CLOSURES)
 19: PETSC_EXTERN PetscErrorCode PetscVFPrintfSetClosure(int (^)(const char*));
 20: #endif


 23: #if defined(PETSC_HAVE_CUDA)
 24: #include <cuda.h>
 25: #include <cublas_v2.h>
 26: #endif

 28: /*
 29:    All major PETSc data structures have a common core; this is defined
 30:    below by PETSCHEADER.

 32:    PetscHeaderCreate() should be used whenever creating a PETSc structure.
 33: */

 35: /*
 36:    PetscOps: structure of core operations that all PETSc objects support.

 38:       getcomm()         - Gets the object's communicator.
 39:       view()            - Is the routine for viewing the entire PETSc object; for
 40:                           example, MatView() is the general matrix viewing routine.
 41:                           This is used by PetscObjectView((PetscObject)obj) to allow
 42:                           viewing any PETSc object.
 43:       destroy()         - Is the routine for destroying the entire PETSc object;
 44:                           for example,MatDestroy() is the general matrix
 45:                           destruction routine.
 46:                           This is used by PetscObjectDestroy((PetscObject*)&obj) to allow
 47:                           destroying any PETSc object.
 48:       compose()         - Associates a PETSc object with another PETSc object with a name
 49:       query()           - Returns a different PETSc object that has been associated
 50:                           with the first object using a name.
 51:       composefunction() - Attaches an a function to a PETSc object with a name.
 52:       queryfunction()   - Requests a registered function that has been attached to a PETSc object.
 53: */

 55: typedef struct {
 56:    PetscErrorCode (*getcomm)(PetscObject,MPI_Comm *);
 57:    PetscErrorCode (*view)(PetscObject,PetscViewer);
 58:    PetscErrorCode (*destroy)(PetscObject*);
 59:    PetscErrorCode (*compose)(PetscObject,const char[],PetscObject);
 60:    PetscErrorCode (*query)(PetscObject,const char[],PetscObject *);
 61:    PetscErrorCode (*composefunction)(PetscObject,const char[],void (*)(void));
 62:    PetscErrorCode (*queryfunction)(PetscObject,const char[],void (**)(void));
 63: } PetscOps;

 65: typedef enum {PETSC_FORTRAN_CALLBACK_CLASS,PETSC_FORTRAN_CALLBACK_SUBTYPE,PETSC_FORTRAN_CALLBACK_MAXTYPE} PetscFortranCallbackType;
 66: typedef int PetscFortranCallbackId;
 67: #define PETSC_SMALLEST_FORTRAN_CALLBACK ((PetscFortranCallbackId)1000)
 68: PETSC_EXTERN PetscErrorCode PetscFortranCallbackRegister(PetscClassId,const char*,PetscFortranCallbackId*);
 69: PETSC_EXTERN PetscErrorCode PetscFortranCallbackGetSizes(PetscClassId,PetscInt*,PetscInt*);

 71: typedef struct {
 72:   void (*func)(void);
 73:   void *ctx;
 74: } PetscFortranCallback;

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

121:   PetscInt             noptionhandler;
122:   PetscErrorCode       (*optionhandler[PETSC_MAX_OPTIONS_HANDLER])(PetscOptionItems*,PetscObject,void*);
123:   PetscErrorCode       (*optiondestroy[PETSC_MAX_OPTIONS_HANDLER])(PetscObject,void*);
124:   void                 *optionctx[PETSC_MAX_OPTIONS_HANDLER];
125:   PetscBool            optionsprinted;
126: #if defined(PETSC_HAVE_SAWS)
127:   PetscBool            amsmem;          /* if PETSC_TRUE then this object is registered with SAWs and visible to clients */
128:   PetscBool            amspublishblock; /* if PETSC_TRUE and publishing objects then will block at PetscObjectSAWsBlock() */
129: #endif
130:   PetscOptions         options;         /* options database used, NULL means default */
131:   PetscBool            donotPetscObjectPrintClassNamePrefixType;
132: } _p_PetscObject;

134: #define PETSCHEADER(ObjectOps) \
135:   _p_PetscObject hdr;          \
136:   ObjectOps      ops[1]

138: #define  PETSCFREEDHEADER -1

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

143: /*@C
144:     PetscHeaderCreate - Creates a PETSc object of a particular class

146:     Input Parameters:
147: +   classid - the classid associated with this object (for example VEC_CLASSID)
148: .   class_name - string name of class; should be static (for example "Vec")
149: .   descr - string containing short description; should be static (for example "Vector")
150: .   mansec - string indicating section in manual pages; should be static (for example "Vec")
151: .   comm - the MPI Communicator
152: .   destroy - the destroy routine for this object (for example VecDestroy())
153: -   view - the view routine for this object (for example VecView())

155:     Output Parameter:
156: .   h - the newly created object

158:     Level: developer

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

162: @*/
163: #define PetscHeaderCreate(h,classid,class_name,descr,mansec,comm,destroy,view) \
164:   (PetscNew(&(h)) || \
165:    PetscHeaderCreate_Private((PetscObject)h,classid,class_name,descr,mansec,comm,(PetscObjectDestroyFunction)destroy,(PetscObjectViewFunction)view) || \
166:    PetscLogObjectCreate(h) || \
167:    PetscLogObjectMemory((PetscObject)h,sizeof(*(h))))

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

172: /*@C
173:     PetscHeaderDestroy - Final step in destroying a PetscObject

175:     Input Parameters:
176: .   h - the header created with PetscHeaderCreate()

178:     Level: developer

180: .seealso: PetscHeaderCreate()
181: @*/
182: #define PetscHeaderDestroy(h) (PetscHeaderDestroy_Private((PetscObject)(*h)) || PetscFree(*h))

184: PETSC_EXTERN PetscErrorCode PetscHeaderDestroy_Private(PetscObject);
185: PETSC_EXTERN PetscErrorCode PetscObjectCopyFortranFunctionPointers(PetscObject,PetscObject);
186: PETSC_EXTERN PetscErrorCode PetscObjectSetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId*,void(*)(void),void *ctx);
187: PETSC_EXTERN PetscErrorCode PetscObjectGetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId,void(**)(void),void **ctx);

189: PETSC_INTERN PetscErrorCode PetscCitationsInitialize(void);
190: PETSC_INTERN PetscErrorCode PetscOptionsFindPair_Private(PetscOptions,const char[],const char[],char**,PetscBool*);
191: PETSC_INTERN PetscErrorCode PetscFreeMPIResources(void);




197: /*
198:     Macros to test if a PETSc object is valid and if pointers are valid
199: */
200: #if !defined(PETSC_USE_DEBUG)


211: #else

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

224:   do {                                                                  \
225:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object: Parameter # %d",arg); \
227:     if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free: Parameter # %d",arg); \
228:     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); \
229:   } while (0)

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

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

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

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

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

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

266: #endif

268: #if !defined(PETSC_USE_DEBUG)


283: #else

285: /*
286:     For example, in the dot product between two vectors,
287:   both vectors must be either Seq or MPI, not one of each
288: */
290:   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);
291: /*
292:     Check type_name
293: */
295:   do {                                                                  \
296:     PetscBool      __match;                                             \
297:     PetscErrorCode _7_ierr;                                             \
298:     _7_PetscObjectTypeCompare(((PetscObject)a),(type),&__match);CHKERRQ(_7_ierr); \
299:     if (!__match) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Object (%s) is not %s",(char*)(((PetscObject)a)->type_name),type); \
300:   } while (0)

303:   do {                                                                  \
304:     PetscBool       __match;                                            \
305:     PetscErrorCode _7_ierr;                                             \
306:     _7_PetscObjectTypeCompareAny(((PetscObject)a),&__match,(type1),(type2),"");CHKERRQ(_7_ierr); \
307:     if (!__match) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Object (%s) is not %s or %s",(char*)(((PetscObject)a)->type_name),type1,type2); \
308:   } while (0)
309: /*
310:    Use this macro to check if the type is set
311: */
313:   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);
314: /*
315:    Sometimes object must live on same communicator to inter-operate
316: */
318:   do {                                                                  \
319:     PetscErrorCode _6_ierr,__flag;                                      \
320:     _6_MPI_Comm_compare(PetscObjectComm((PetscObject)a),PetscObjectComm((PetscObject)b),&__flag);CHKERRQ(_6_ierr);                                                   \
321:     if (__flag != MPI_CONGRUENT && __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,__flag); \
322:   } while (0)

325:   do {                                                  \
328:   } while (0)

331:   do {                                                                  \
332:     PetscErrorCode _7_ierr;                                             \
333:     PetscReal b1[5],b2[5];                                              \
334:     if (PetscIsNanScalar(b)) {b1[4] = 1;} else {b1[4] = 0;};            \
335:     b1[0] = -PetscRealPart(b); b1[1] = PetscRealPart(b);b1[2] = -PetscImaginaryPart(b); b1[3] = PetscImaginaryPart(b);         \
336:     _7_MPI_Allreduce(b1,b2,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
337:     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",c); \
338:   } while (0)

341:   do {                                                                  \
342:     PetscErrorCode _7_ierr;                                             \
343:     PetscReal b1[3],b2[3];                                              \
344:     if (PetscIsNanReal(b)) {b1[2] = 1;} else {b1[2] = 0;};              \
345:     b1[0] = -b; b1[1] = b;                                              \
346:     _7_MPI_Allreduce(b1,b2,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
347:     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",c); \
348:   } while (0)

351:   do {                                                                  \
352:     PetscErrorCode _7_ierr;                                             \
353:     PetscInt b1[2],b2[2];                                               \
354:     b1[0] = -b; b1[1] = b;                                              \
355:     _7_MPIU_Allreduce(b1,b2,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
356:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Int value must be same on all processes, argument # %d",c); \
357:   } while (0)


362:   do {                                                                  \
363:     PetscErrorCode _7_ierr;                                             \
364:     PetscMPIInt b1[2],b2[2];                                            \
365:     b1[0] = -(PetscMPIInt)b; b1[1] = (PetscMPIInt)b;                    \
366:     _7_MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
367:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Bool value must be same on all processes, argument # %d",c); \
368:   } while (0)

371:   do {                                                                  \
372:     PetscErrorCode _7_ierr;                                             \
373:     PetscMPIInt b1[2],b2[2];                                            \
374:     b1[0] = -(PetscMPIInt)b; b1[1] = (PetscMPIInt)b;                    \
375:     _7_MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
376:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes, argument # %d",c); \
377:   } while (0)

379: #endif

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

385:    Level: developer

387: .seealso: PetscUseMethod()
388: */
389: #define  PetscTryMethod(obj,A,B,C) \
390:   0;{ PetscErrorCode (*f)B, __ierr; \
391:     __PetscObjectQueryFunction((PetscObject)obj,A,&f);CHKERRQ(__ierr); \
392:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
393:   }

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

399:    Level: developer

401: .seealso: PetscTryMethod()
402: */
403: #define  PetscUseMethod(obj,A,B,C) \
404:   0;{ PetscErrorCode (*f)B, __ierr; \
405:     __PetscObjectQueryFunction((PetscObject)obj,A,&f);CHKERRQ(__ierr); \
406:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
407:     else SETERRQ1(PetscObjectComm((PetscObject)obj),PETSC_ERR_SUP,"Cannot locate function %s in object",A); \
408:   }

410: /*MC
411:    PetscObjectStateIncrease - Increases the state of any PetscObject

413:    Synopsis:
414:    #include "petsc/private/petscimpl.h"
415:    PetscErrorCode PetscObjectStateIncrease(PetscObject obj)

417:    Logically Collective

419:    Input Parameter:
420: .  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
421:          cast with a (PetscObject), for example,
422:          PetscObjectStateIncrease((PetscObject)mat);

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

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

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

436:    Level: developer

438:    seealso: PetscObjectStateGet()

440:    Concepts: state

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

445: PETSC_EXTERN PetscErrorCode PetscObjectStateGet(PetscObject,PetscObjectState*);
446: PETSC_EXTERN PetscErrorCode PetscObjectStateSet(PetscObject,PetscObjectState);
447: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataRegister(PetscInt*);
448: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseInt(PetscObject);
449: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseIntstar(PetscObject);
450: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseReal(PetscObject);
451: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseRealstar(PetscObject);
452: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalar(PetscObject);
453: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalarstar(PetscObject);
454: PETSC_EXTERN PetscInt         PetscObjectComposedDataMax;
455: /*MC
456:    PetscObjectComposedDataSetInt - attach integer data to a PetscObject

458:    Synopsis:
459:    #include "petsc/private/petscimpl.h"
460:    PetscErrorCode PetscObjectComposedDataSetInt(PetscObject obj,int id,int data)

462:    Not collective

464:    Input parameters:
465: +  obj - the object to which data is to be attached
466: .  id - the identifier for the data
467: -  data - the data to  be attached

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

472:    Level: developer
473: M*/
474: #define PetscObjectComposedDataSetInt(obj,id,data)                                      \
475:   ((((obj)->int_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseInt(obj)) ||  \
476:    ((obj)->intcomposeddata[id] = data,(obj)->intcomposedstate[id] = (obj)->state, 0))

478: /*MC
479:    PetscObjectComposedDataGetInt - retrieve integer data attached to an object

481:    Synopsis:
482:    #include "petsc/private/petscimpl.h"
483:    PetscErrorCode PetscObjectComposedDataGetInt(PetscObject obj,int id,int data,PetscBool  flag)

485:    Not collective

487:    Input parameters:
488: +  obj - the object from which data is to be retrieved
489: -  id - the identifier for the data

491:    Output parameters
492: +  data - the data to be retrieved
493: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

497:    Level: developer
498: M*/
499: #define PetscObjectComposedDataGetInt(obj,id,data,flag)                            \
500:   ((((obj)->intcomposedstate && ((obj)->intcomposedstate[id] == (obj)->state)) ?   \
501:    (data = (obj)->intcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

503: /*MC
504:    PetscObjectComposedDataSetIntstar - attach integer array data to a PetscObject

506:    Synopsis:
507:    #include "petsc/private/petscimpl.h"
508:    PetscErrorCode PetscObjectComposedDataSetIntstar(PetscObject obj,int id,int *data)

510:    Not collective

512:    Input parameters:
513: +  obj - the object to which data is to be attached
514: .  id - the identifier for the data
515: -  data - the data to  be attached

517:    Notes
518:    The data identifier can best be determined through a call to
519:    PetscObjectComposedDataRegister()

521:    Level: developer
522: M*/
523: #define PetscObjectComposedDataSetIntstar(obj,id,data)                                          \
524:   ((((obj)->intstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseIntstar(obj)) ||  \
525:    ((obj)->intstarcomposeddata[id] = data,(obj)->intstarcomposedstate[id] = (obj)->state, 0))

527: /*MC
528:    PetscObjectComposedDataGetIntstar - retrieve integer array data
529:    attached to an object

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

535:    Not collective

537:    Input parameters:
538: +  obj - the object from which data is to be retrieved
539: -  id - the identifier for the data

541:    Output parameters
542: +  data - the data to be retrieved
543: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

547:    Level: developer
548: M*/
549: #define PetscObjectComposedDataGetIntstar(obj,id,data,flag)                               \
550:   ((((obj)->intstarcomposedstate && ((obj)->intstarcomposedstate[id] == (obj)->state)) ?  \
551:    (data = (obj)->intstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

553: /*MC
554:    PetscObjectComposedDataSetReal - attach real data to a PetscObject

556:    Synopsis:
557:    #include "petsc/private/petscimpl.h"
558:    PetscErrorCode PetscObjectComposedDataSetReal(PetscObject obj,int id,PetscReal data)

560:    Not collective

562:    Input parameters:
563: +  obj - the object to which data is to be attached
564: .  id - the identifier for the data
565: -  data - the data to  be attached

567:    Notes
568:    The data identifier can best be determined through a call to
569:    PetscObjectComposedDataRegister()

571:    Level: developer
572: M*/
573: #define PetscObjectComposedDataSetReal(obj,id,data)                                       \
574:   ((((obj)->real_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseReal(obj)) ||  \
575:    ((obj)->realcomposeddata[id] = data,(obj)->realcomposedstate[id] = (obj)->state, 0))

577: /*MC
578:    PetscObjectComposedDataGetReal - retrieve real data attached to an object

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

584:    Not collective

586:    Input parameters:
587: +  obj - the object from which data is to be retrieved
588: -  id - the identifier for the data

590:    Output parameters
591: +  data - the data to be retrieved
592: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

596:    Level: developer
597: M*/
598: #define PetscObjectComposedDataGetReal(obj,id,data,flag)                            \
599:   ((((obj)->realcomposedstate && ((obj)->realcomposedstate[id] == (obj)->state)) ?  \
600:    (data = (obj)->realcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

602: /*MC
603:    PetscObjectComposedDataSetRealstar - attach real array data to a PetscObject

605:    Synopsis:
606:    #include "petsc/private/petscimpl.h"
607:    PetscErrorCode PetscObjectComposedDataSetRealstar(PetscObject obj,int id,PetscReal *data)

609:    Not collective

611:    Input parameters:
612: +  obj - the object to which data is to be attached
613: .  id - the identifier for the data
614: -  data - the data to  be attached

616:    Notes
617:    The data identifier can best be determined through a call to
618:    PetscObjectComposedDataRegister()

620:    Level: developer
621: M*/
622: #define PetscObjectComposedDataSetRealstar(obj,id,data)                                           \
623:   ((((obj)->realstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseRealstar(obj)) ||  \
624:    ((obj)->realstarcomposeddata[id] = data, (obj)->realstarcomposedstate[id] = (obj)->state, 0))

626: /*MC
627:    PetscObjectComposedDataGetRealstar - retrieve real array data
628:    attached to an object

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

634:    Not collective

636:    Input parameters:
637: +  obj - the object from which data is to be retrieved
638: -  id - the identifier for the data

640:    Output parameters
641: +  data - the data to be retrieved
642: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

646:    Level: developer
647: M*/
648: #define PetscObjectComposedDataGetRealstar(obj,id,data,flag)                                \
649:   ((((obj)->realstarcomposedstate && ((obj)->realstarcomposedstate[id] == (obj)->state)) ?  \
650:    (data = (obj)->realstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

652: /*MC
653:    PetscObjectComposedDataSetScalar - attach scalar data to a PetscObject

655:    Synopsis:
656:    #include "petsc/private/petscimpl.h"
657:    PetscErrorCode PetscObjectComposedDataSetScalar(PetscObject obj,int id,PetscScalar data)

659:    Not collective

661:    Input parameters:
662: +  obj - the object to which data is to be attached
663: .  id - the identifier for the data
664: -  data - the data to  be attached

666:    Notes
667:    The data identifier can best be determined through a call to
668:    PetscObjectComposedDataRegister()

670:    Level: developer
671: M*/
672: #if defined(PETSC_USE_COMPLEX)
673: #define PetscObjectComposedDataSetScalar(obj,id,data)                                        \
674:   ((((obj)->scalar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalar(obj)) || \
675:    ((obj)->scalarcomposeddata[id] = data,(obj)->scalarcomposedstate[id] = (obj)->state, 0))
676: #else
677: #define PetscObjectComposedDataSetScalar(obj,id,data) \
678:         PetscObjectComposedDataSetReal(obj,id,data)
679: #endif
680: /*MC
681:    PetscObjectComposedDataGetScalar - retrieve scalar data attached to an object

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

687:    Not collective

689:    Input parameters:
690: +  obj - the object from which data is to be retrieved
691: -  id - the identifier for the data

693:    Output parameters
694: +  data - the data to be retrieved
695: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

699:    Level: developer
700: M*/
701: #if defined(PETSC_USE_COMPLEX)
702: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                              \
703:   ((((obj)->scalarcomposedstate && ((obj)->scalarcomposedstate[id] == (obj)->state) ) ? \
704:    (data = (obj)->scalarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
705: #else
706: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                             \
707:         PetscObjectComposedDataGetReal(obj,id,data,flag)
708: #endif

710: /*MC
711:    PetscObjectComposedDataSetScalarstar - attach scalar array data to a PetscObject

713:    Synopsis:
714:    #include "petsc/private/petscimpl.h"
715:    PetscErrorCode PetscObjectComposedDataSetScalarstar(PetscObject obj,int id,PetscScalar *data)

717:    Not collective

719:    Input parameters:
720: +  obj - the object to which data is to be attached
721: .  id - the identifier for the data
722: -  data - the data to  be attached

724:    Notes
725:    The data identifier can best be determined through a call to
726:    PetscObjectComposedDataRegister()

728:    Level: developer
729: M*/
730: #if defined(PETSC_USE_COMPLEX)
731: #define PetscObjectComposedDataSetScalarstar(obj,id,data)                                             \
732:   ((((obj)->scalarstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalarstar(obj)) ||  \
733:    ((obj)->scalarstarcomposeddata[id] = data,(obj)->scalarstarcomposedstate[id] = (obj)->state, 0))
734: #else
735: #define PetscObjectComposedDataSetScalarstar(obj,id,data) \
736:         PetscObjectComposedDataSetRealstar(obj,id,data)
737: #endif
738: /*MC
739:    PetscObjectComposedDataGetScalarstar - retrieve scalar array data
740:    attached to an object

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

746:    Not collective

748:    Input parameters:
749: +  obj - the object from which data is to be retrieved
750: -  id - the identifier for the data

752:    Output parameters
753: +  data - the data to be retrieved
754: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

758:    Level: developer
759: M*/
760: #if defined(PETSC_USE_COMPLEX)
761: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)                                 \
762:   ((((obj)->scalarstarcomposedstate && ((obj)->scalarstarcomposedstate[id] == (obj)->state)) ? \
763:        (data = (obj)->scalarstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
764: #else
765: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)         \
766:         PetscObjectComposedDataGetRealstar(obj,id,data,flag)
767: #endif

769: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);

771: PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscBool *);

773: PETSC_EXTERN PetscMPIInt Petsc_Counter_keyval;
774: PETSC_EXTERN PetscMPIInt Petsc_InnerComm_keyval;
775: PETSC_EXTERN PetscMPIInt Petsc_OuterComm_keyval;
776: PETSC_EXTERN PetscMPIInt Petsc_Seq_keyval;
777: PETSC_EXTERN PetscMPIInt Petsc_Shared_keyval;

779: /*
780:   PETSc communicators have this attribute, see
781:   PetscCommDuplicate(), PetscCommDestroy(), PetscCommGetNewTag(), PetscObjectGetName()
782: */
783: typedef struct {
784:   PetscMPIInt tag;              /* next free tag value */
785:   PetscInt    refcount;         /* number of references, communicator can be freed when this reaches 0 */
786:   PetscInt    namecount;        /* used to generate the next name, as in Vec_0, Mat_1, ... */
787: } PetscCommCounter;

789: /*E
790:     PetscOffloadFlag - indicates which memory (CPU, GPU, or none contains valid vector

792:    PETSC_OFFLOAD_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
793:    PETSC_OFFLOAD_GPU - GPU has valid vector/matrix entries
794:    PETSC_OFFLOAD_CPU - CPU has valid vector/matrix entries
795:    PETSC_OFFLOAD_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

797:    Level: developer
798: E*/
799: typedef enum {PETSC_OFFLOAD_UNALLOCATED,PETSC_OFFLOAD_GPU,PETSC_OFFLOAD_CPU,PETSC_OFFLOAD_BOTH} PetscOffloadFlag;

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

803: #define REDUCE_SUM  0
804: #define REDUCE_MAX  1
805: #define REDUCE_MIN  2

807: typedef struct {
808:   MPI_Comm    comm;
809:   MPI_Request request;
810:   PetscBool   async;
811:   PetscScalar *lvalues;     /* this are the reduced values before call to MPI_Allreduce() */
812:   PetscScalar *gvalues;     /* values after call to MPI_Allreduce() */
813:   void        **invecs;     /* for debugging only, vector/memory used with each op */
814:   PetscInt    *reducetype;  /* is particular value to be summed or maxed? */
815:   SRState     state;        /* are we calling xxxBegin() or xxxEnd()? */
816:   PetscInt    maxops;       /* total amount of space we have for requests */
817:   PetscInt    numopsbegin;  /* number of requests that have been queued in */
818:   PetscInt    numopsend;    /* number of requests that have been gotten by user */
819: } PetscSplitReduction;

821: PETSC_EXTERN PetscErrorCode PetscSplitReductionGet(MPI_Comm,PetscSplitReduction**);
822: PETSC_EXTERN PetscErrorCode PetscSplitReductionEnd(PetscSplitReduction*);
823: PETSC_EXTERN PetscErrorCode PetscSplitReductionExtend(PetscSplitReduction*);

825: #if !defined(PETSC_SKIP_SPINLOCK)
826: #if defined(PETSC_HAVE_THREADSAFETY)
827: #  if defined(PETSC_HAVE_CONCURRENCYKIT)
828: #if defined(__cplusplus)
829: /*  CK does not have extern "C" protection in their include files */
830: extern "C" {
831: #endif
832: #include <ck_spinlock.h>
833: #if defined(__cplusplus)
834: }
835: #endif
836: typedef ck_spinlock_t PetscSpinlock;
837: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockCreate(PetscSpinlock *ck_spinlock)
838: {
839:   ck_spinlock_init(ck_spinlock);
840:   return 0;
841: }
842: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockLock(PetscSpinlock *ck_spinlock)
843: {
844:   ck_spinlock_lock(ck_spinlock);
845:   return 0;
846: }
847: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockUnlock(PetscSpinlock *ck_spinlock)
848: {
849:   ck_spinlock_unlock(ck_spinlock);
850:   return 0;
851: }
852: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockDestroy(PetscSpinlock *ck_spinlock)
853: {
854:   return 0;
855: }
856: #  elif defined(PETSC_HAVE_OPENMP)

858: #include <omp.h>
859: typedef omp_lock_t PetscSpinlock;
860: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockCreate(PetscSpinlock *omp_lock)
861: {
862:   omp_init_lock(omp_lock);
863:   return 0;
864: }
865: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockLock(PetscSpinlock *omp_lock)
866: {
867:   omp_set_lock(omp_lock);
868:   return 0;
869: }
870: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockUnlock(PetscSpinlock *omp_lock)
871: {
872:   omp_unset_lock(omp_lock);
873:   return 0;
874: }
875: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockDestroy(PetscSpinlock *omp_lock)
876: {
877:   omp_destroy_lock(omp_lock);
878:   return 0;
879: }
880: #else
881: Thread safety requires either --with-openmp or --download-concurrencykit
882: #endif

884: #else
885: typedef int PetscSpinlock;
886: #define PetscSpinlockCreate(a)  0
887: #define PetscSpinlockLock(a)    0
888: #define PetscSpinlockUnlock(a)  0
889: #define PetscSpinlockDestroy(a) 0
890: #endif

892: #if defined(PETSC_HAVE_THREADSAFETY)
893: extern PetscSpinlock PetscViewerASCIISpinLockOpen;
894: extern PetscSpinlock PetscViewerASCIISpinLockStdout;
895: extern PetscSpinlock PetscViewerASCIISpinLockStderr;
896: extern PetscSpinlock PetscCommSpinLock;
897: #endif
898: #endif

900: #endif /* _PETSCHEAD_H */