Actual source code: petscimpl.h

petsc-3.7.5 2017-01-01
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: #if defined(PETSC_HAVE_CUDA)
 11: #include <cuda.h>
 12: #include <cublas_v2.h>
 13: #endif

 15: /*
 16:    All major PETSc data structures have a common core; this is defined
 17:    below by PETSCHEADER.

 19:    PetscHeaderCreate() should be used whenever creating a PETSc structure.
 20: */

 22: /*
 23:    PetscOps: structure of core operations that all PETSc objects support.

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

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

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

 58: typedef struct {
 59:   void (*func)(void);
 60:   void *ctx;
 61: } PetscFortranCallback;

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

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

121: #define PETSCHEADER(ObjectOps) \
122:   _p_PetscObject hdr;          \
123:   ObjectOps      ops[1]

125: #define  PETSCFREEDHEADER -1

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

130: /*@C
131:     PetscHeaderCreate - Creates a PETSc object of a particular class

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

142:     Output Parameter:
143: .   h - the newly created object

145:     Level: developer

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

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

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

159: /*@C
160:     PetscHeaderDestroy - Final step in destroying a PetscObject

162:     Input Parameters:
163: .   h - the header created with PetscHeaderCreate()

165:     Level: developer

167: .seealso: PetscHeaderCreate()
168: @*/
169: #define PetscHeaderDestroy(h) (PetscHeaderDestroy_Private((PetscObject)(*h)) || PetscFree(*h))

171: PETSC_EXTERN PetscErrorCode PetscHeaderDestroy_Private(PetscObject);
172: PETSC_EXTERN PetscErrorCode PetscObjectCopyFortranFunctionPointers(PetscObject,PetscObject);
173: PETSC_EXTERN PetscErrorCode PetscObjectSetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId*,void(*)(void),void *ctx);
174: PETSC_EXTERN PetscErrorCode PetscObjectGetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId,void(**)(void),void **ctx);

176: PETSC_INTERN PetscErrorCode PetscCitationsInitialize(void);
177: PETSC_INTERN PetscErrorCode PetscOptionsFindPair_Private(PetscOptions,const char[],const char[],char**,PetscBool*);


181: /*
182:     Macros to test if a PETSc object is valid and if pointers are valid
183: */
184: #if !defined(PETSC_USE_DEBUG)


195: #else

198:   do {                                                                  \
199:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object: Parameter # %d",arg); \
201:     if (((PetscObject)(h))->classid != ck) {                            \
202:       if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free: Parameter # %d",arg); \
203:       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong type of object: Parameter # %d",arg); \
204:     }                                                                   \
205:   } while (0)

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

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

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

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

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

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

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

250: #endif

252: #if !defined(PETSC_USE_DEBUG)


265: #else

267: /*
268:     For example, in the dot product between two vectors,
269:   both vectors must be either Seq or MPI, not one of each
270: */
272:   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);
273: /*
274:    Use this macro to check if the type is set
275: */
277:   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);
278: /*
279:    Sometimes object must live on same communicator to inter-operate
280: */
282:   do {                                                                  \
283:     PetscErrorCode _6_ierr,__flag;                                      \
284:     _6_MPI_Comm_compare(PetscObjectComm((PetscObject)a),PetscObjectComm((PetscObject)b),&__flag);CHKERRQ(_6_ierr);                                                   \
285:     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); \
286:   } while (0)

289:   do {                                                  \
292:   } while (0)

295:   do {                                                                  \
296:     PetscErrorCode _7_ierr;                                             \
297:     PetscReal b1[6],b2[6];                                              \
298:     if (PetscIsNanScalar(b)) {b1[4] = -1; b1[5] = 1;} else b1[4] = b1[5] = 0; \
299:     b1[0] = -PetscRealPart(b); b1[1] = PetscRealPart(b);b1[2] = -PetscImaginaryPart(b); b1[3] = PetscImaginaryPart(b);         \
300:     _7_MPI_Allreduce(b1,b2,6,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
301:     if (!(b2[4] == -1 && b2[5] == 1) && (-b2[0] != b2[1] || -b2[2] != b2[3])) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Scalar value must be same on all processes, argument # %d",c); \
302:   } while (0)

305:   do {                                                                  \
306:     PetscErrorCode _7_ierr;                                             \
307:     PetscReal b1[4],b2[4];                                              \
308:     if (PetscIsNanReal(b)) {b1[2] = -1; b1[3] = 1;} else b1[2] = b1[3] = 0; \
309:     b1[0] = -b; b1[1] = b;                                                  \
310:     _7_MPI_Allreduce(b1,b2,4,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
311:     if (!(b2[2] == -1 && b2[3] == 1) && -b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Real value must be same on all processes, argument # %d",c); \
312:   } while (0)

315:   do {                                                                  \
316:     PetscErrorCode _7_ierr;                                             \
317:     PetscInt b1[2],b2[2];                                               \
318:     b1[0] = -b; b1[1] = b;                                              \
319:     _7_MPIU_Allreduce(b1,b2,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
320:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Int value must be same on all processes, argument # %d",c); \
321:   } while (0)


326:   do {                                                                  \
327:     PetscErrorCode _7_ierr;                                             \
328:     PetscMPIInt b1[2],b2[2];                                            \
329:     b1[0] = -(PetscMPIInt)b; b1[1] = (PetscMPIInt)b;                    \
330:     _7_MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
331:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Bool value must be same on all processes, argument # %d",c); \
332:   } while (0)

335:   do {                                                                  \
336:     PetscErrorCode _7_ierr;                                             \
337:     PetscMPIInt b1[2],b2[2];                                            \
338:     b1[0] = -(PetscMPIInt)b; b1[1] = (PetscMPIInt)b;                    \
339:     _7_MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
340:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes, argument # %d",c); \
341:   } while (0)

343: #endif

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

349:    Level: developer

351: .seealso: PetscUseMethod()
352: */
353: #define  PetscTryMethod(obj,A,B,C) \
354:   0;{ PetscErrorCode (*f)B, __ierr; \
355:     __PetscObjectQueryFunction((PetscObject)obj,A,&f);CHKERRQ(__ierr); \
356:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
357:   }

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

363:    Level: developer

365: .seealso: PetscTryMethod()
366: */
367: #define  PetscUseMethod(obj,A,B,C) \
368:   0;{ PetscErrorCode (*f)B, __ierr; \
369:     __PetscObjectQueryFunction((PetscObject)obj,A,&f);CHKERRQ(__ierr); \
370:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
371:     else SETERRQ1(PetscObjectComm((PetscObject)obj),PETSC_ERR_SUP,"Cannot locate function %s in object",A); \
372:   }

374: /*MC
375:    PetscObjectStateIncrease - Increases the state of any PetscObject

377:    Synopsis:
378:    #include "petsc/private/petscimpl.h"
379:    PetscErrorCode PetscObjectStateIncrease(PetscObject obj)

381:    Logically Collective

383:    Input Parameter:
384: .  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
385:          cast with a (PetscObject), for example,
386:          PetscObjectStateIncrease((PetscObject)mat);

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

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

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

400:    Level: developer

402:    seealso: PetscObjectStateGet()

404:    Concepts: state

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

409: PETSC_EXTERN PetscErrorCode PetscObjectStateGet(PetscObject,PetscObjectState*);
410: PETSC_EXTERN PetscErrorCode PetscObjectStateSet(PetscObject,PetscObjectState);
411: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataRegister(PetscInt*);
412: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseInt(PetscObject);
413: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseIntstar(PetscObject);
414: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseReal(PetscObject);
415: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseRealstar(PetscObject);
416: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalar(PetscObject);
417: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalarstar(PetscObject);
418: PETSC_EXTERN PetscInt         PetscObjectComposedDataMax;
419: /*MC
420:    PetscObjectComposedDataSetInt - attach integer data to a PetscObject

422:    Synopsis:
423:    #include "petsc/private/petscimpl.h"
424:    PetscErrorCode PetscObjectComposedDataSetInt(PetscObject obj,int id,int data)

426:    Not collective

428:    Input parameters:
429: +  obj - the object to which data is to be attached
430: .  id - the identifier for the data
431: -  data - the data to  be attached

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

436:    Level: developer
437: M*/
438: #define PetscObjectComposedDataSetInt(obj,id,data)                                      \
439:   ((((obj)->int_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseInt(obj)) ||  \
440:    ((obj)->intcomposeddata[id] = data,(obj)->intcomposedstate[id] = (obj)->state, 0))

442: /*MC
443:    PetscObjectComposedDataGetInt - retrieve integer data attached to an object

445:    Synopsis:
446:    #include "petsc/private/petscimpl.h"
447:    PetscErrorCode PetscObjectComposedDataGetInt(PetscObject obj,int id,int data,PetscBool  flag)

449:    Not collective

451:    Input parameters:
452: +  obj - the object from which data is to be retrieved
453: -  id - the identifier for the data

455:    Output parameters
456: +  data - the data to be retrieved
457: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

461:    Level: developer
462: M*/
463: #define PetscObjectComposedDataGetInt(obj,id,data,flag)                            \
464:   ((((obj)->intcomposedstate && ((obj)->intcomposedstate[id] == (obj)->state)) ?   \
465:    (data = (obj)->intcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

467: /*MC
468:    PetscObjectComposedDataSetIntstar - attach integer array data to a PetscObject

470:    Synopsis:
471:    #include "petsc/private/petscimpl.h"
472:    PetscErrorCode PetscObjectComposedDataSetIntstar(PetscObject obj,int id,int *data)

474:    Not collective

476:    Input parameters:
477: +  obj - the object to which data is to be attached
478: .  id - the identifier for the data
479: -  data - the data to  be attached

481:    Notes
482:    The data identifier can best be determined through a call to
483:    PetscObjectComposedDataRegister()

485:    Level: developer
486: M*/
487: #define PetscObjectComposedDataSetIntstar(obj,id,data)                                          \
488:   ((((obj)->intstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseIntstar(obj)) ||  \
489:    ((obj)->intstarcomposeddata[id] = data,(obj)->intstarcomposedstate[id] = (obj)->state, 0))

491: /*MC
492:    PetscObjectComposedDataGetIntstar - retrieve integer array data
493:    attached to an object

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

499:    Not collective

501:    Input parameters:
502: +  obj - the object from which data is to be retrieved
503: -  id - the identifier for the data

505:    Output parameters
506: +  data - the data to be retrieved
507: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

511:    Level: developer
512: M*/
513: #define PetscObjectComposedDataGetIntstar(obj,id,data,flag)                               \
514:   ((((obj)->intstarcomposedstate && ((obj)->intstarcomposedstate[id] == (obj)->state)) ?  \
515:    (data = (obj)->intstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

517: /*MC
518:    PetscObjectComposedDataSetReal - attach real data to a PetscObject

520:    Synopsis:
521:    #include "petsc/private/petscimpl.h"
522:    PetscErrorCode PetscObjectComposedDataSetReal(PetscObject obj,int id,PetscReal data)

524:    Not collective

526:    Input parameters:
527: +  obj - the object to which data is to be attached
528: .  id - the identifier for the data
529: -  data - the data to  be attached

531:    Notes
532:    The data identifier can best be determined through a call to
533:    PetscObjectComposedDataRegister()

535:    Level: developer
536: M*/
537: #define PetscObjectComposedDataSetReal(obj,id,data)                                       \
538:   ((((obj)->real_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseReal(obj)) ||  \
539:    ((obj)->realcomposeddata[id] = data,(obj)->realcomposedstate[id] = (obj)->state, 0))

541: /*MC
542:    PetscObjectComposedDataGetReal - retrieve real data attached to an object

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

548:    Not collective

550:    Input parameters:
551: +  obj - the object from which data is to be retrieved
552: -  id - the identifier for the data

554:    Output parameters
555: +  data - the data to be retrieved
556: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

560:    Level: developer
561: M*/
562: #define PetscObjectComposedDataGetReal(obj,id,data,flag)                            \
563:   ((((obj)->realcomposedstate && ((obj)->realcomposedstate[id] == (obj)->state)) ?  \
564:    (data = (obj)->realcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

566: /*MC
567:    PetscObjectComposedDataSetRealstar - attach real array data to a PetscObject

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

573:    Not collective

575:    Input parameters:
576: +  obj - the object to which data is to be attached
577: .  id - the identifier for the data
578: -  data - the data to  be attached

580:    Notes
581:    The data identifier can best be determined through a call to
582:    PetscObjectComposedDataRegister()

584:    Level: developer
585: M*/
586: #define PetscObjectComposedDataSetRealstar(obj,id,data)                                           \
587:   ((((obj)->realstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseRealstar(obj)) ||  \
588:    ((obj)->realstarcomposeddata[id] = data, (obj)->realstarcomposedstate[id] = (obj)->state, 0))

590: /*MC
591:    PetscObjectComposedDataGetRealstar - retrieve real array data
592:    attached to an object

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

598:    Not collective

600:    Input parameters:
601: +  obj - the object from which data is to be retrieved
602: -  id - the identifier for the data

604:    Output parameters
605: +  data - the data to be retrieved
606: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

610:    Level: developer
611: M*/
612: #define PetscObjectComposedDataGetRealstar(obj,id,data,flag)                                \
613:   ((((obj)->realstarcomposedstate && ((obj)->realstarcomposedstate[id] == (obj)->state)) ?  \
614:    (data = (obj)->realstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

616: /*MC
617:    PetscObjectComposedDataSetScalar - attach scalar data to a PetscObject

619:    Synopsis:
620:    #include "petsc/private/petscimpl.h"
621:    PetscErrorCode PetscObjectComposedDataSetScalar(PetscObject obj,int id,PetscScalar data)

623:    Not collective

625:    Input parameters:
626: +  obj - the object to which data is to be attached
627: .  id - the identifier for the data
628: -  data - the data to  be attached

630:    Notes
631:    The data identifier can best be determined through a call to
632:    PetscObjectComposedDataRegister()

634:    Level: developer
635: M*/
636: #if defined(PETSC_USE_COMPLEX)
637: #define PetscObjectComposedDataSetScalar(obj,id,data)                                        \
638:   ((((obj)->scalar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalar(obj)) || \
639:    ((obj)->scalarcomposeddata[id] = data,(obj)->scalarcomposedstate[id] = (obj)->state, 0))
640: #else
641: #define PetscObjectComposedDataSetScalar(obj,id,data) \
642:         PetscObjectComposedDataSetReal(obj,id,data)
643: #endif
644: /*MC
645:    PetscObjectComposedDataGetScalar - retrieve scalar data attached to an object

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

651:    Not collective

653:    Input parameters:
654: +  obj - the object from which data is to be retrieved
655: -  id - the identifier for the data

657:    Output parameters
658: +  data - the data to be retrieved
659: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

663:    Level: developer
664: M*/
665: #if defined(PETSC_USE_COMPLEX)
666: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                              \
667:   ((((obj)->scalarcomposedstate && ((obj)->scalarcomposedstate[id] == (obj)->state) ) ? \
668:    (data = (obj)->scalarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
669: #else
670: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                             \
671:         PetscObjectComposedDataGetReal(obj,id,data,flag)
672: #endif

674: /*MC
675:    PetscObjectComposedDataSetScalarstar - attach scalar array data to a PetscObject

677:    Synopsis:
678:    #include "petsc/private/petscimpl.h"
679:    PetscErrorCode PetscObjectComposedDataSetScalarstar(PetscObject obj,int id,PetscScalar *data)

681:    Not collective

683:    Input parameters:
684: +  obj - the object to which data is to be attached
685: .  id - the identifier for the data
686: -  data - the data to  be attached

688:    Notes
689:    The data identifier can best be determined through a call to
690:    PetscObjectComposedDataRegister()

692:    Level: developer
693: M*/
694: #if defined(PETSC_USE_COMPLEX)
695: #define PetscObjectComposedDataSetScalarstar(obj,id,data)                                             \
696:   ((((obj)->scalarstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalarstar(obj)) ||  \
697:    ((obj)->scalarstarcomposeddata[id] = data,(obj)->scalarstarcomposedstate[id] = (obj)->state, 0))
698: #else
699: #define PetscObjectComposedDataSetScalarstar(obj,id,data) \
700:         PetscObjectComposedDataSetRealstar(obj,id,data)
701: #endif
702: /*MC
703:    PetscObjectComposedDataGetScalarstar - retrieve scalar array data
704:    attached to an object

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

710:    Not collective

712:    Input parameters:
713: +  obj - the object from which data is to be retrieved
714: -  id - the identifier for the data

716:    Output parameters
717: +  data - the data to be retrieved
718: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

722:    Level: developer
723: M*/
724: #if defined(PETSC_USE_COMPLEX)
725: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)                                 \
726:   ((((obj)->scalarstarcomposedstate && ((obj)->scalarstarcomposedstate[id] == (obj)->state)) ? \
727:        (data = (obj)->scalarstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
728: #else
729: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)         \
730:         PetscObjectComposedDataGetRealstar(obj,id,data,flag)
731: #endif

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

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

737: PETSC_EXTERN PetscMPIInt Petsc_Counter_keyval;
738: PETSC_EXTERN PetscMPIInt Petsc_InnerComm_keyval;
739: PETSC_EXTERN PetscMPIInt Petsc_OuterComm_keyval;

741: /*
742:   PETSc communicators have this attribute, see
743:   PetscCommDuplicate(), PetscCommDestroy(), PetscCommGetNewTag(), PetscObjectGetName()
744: */
745: typedef struct {
746:   PetscMPIInt tag;              /* next free tag value */
747:   PetscInt    refcount;         /* number of references, communicator can be freed when this reaches 0 */
748:   PetscInt    namecount;        /* used to generate the next name, as in Vec_0, Mat_1, ... */
749: } PetscCommCounter;

751: #if defined(PETSC_HAVE_CUSP)
752: /*E
753:     PetscCUSPFlag - indicates which memory (CPU, GPU, or none contains valid vector

755:    PETSC_CUSP_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
756:    PETSC_CUSP_GPU - GPU has valid vector/matrix entries
757:    PETSC_CUSP_CPU - CPU has valid vector/matrix entries
758:    PETSC_CUSP_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

760:    Level: developer
761: E*/
762: typedef enum {PETSC_CUSP_UNALLOCATED,PETSC_CUSP_GPU,PETSC_CUSP_CPU,PETSC_CUSP_BOTH} PetscCUSPFlag;
763: #elif defined(PETSC_HAVE_VIENNACL)
764: /*E
765:     PetscViennaCLFlag - indicates which memory (CPU, GPU, or none contains valid vector

767:    PETSC_VIENNACL_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
768:    PETSC_VIENNACL_GPU - GPU has valid vector/matrix entries
769:    PETSC_VIENNACL_CPU - CPU has valid vector/matrix entries
770:    PETSC_VIENNACL_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

772:    Level: developer
773: E*/
774: typedef enum {PETSC_VIENNACL_UNALLOCATED,PETSC_VIENNACL_GPU,PETSC_VIENNACL_CPU,PETSC_VIENNACL_BOTH} PetscViennaCLFlag;
775: #elif defined(PETSC_HAVE_VECCUDA)
776: /*E
777:     PetscCUDAFlag - indicates which memory (CPU, GPU, or none contains valid vector

779:    PETSC_CUDA_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
780:    PETSC_CUDA_GPU - GPU has valid vector/matrix entries
781:    PETSC_CUDA_CPU - CPU has valid vector/matrix entries
782:    PETSC_CUDA_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

784:    Level: developer
785: E*/
786: typedef enum {PETSC_CUDA_UNALLOCATED,PETSC_CUDA_GPU,PETSC_CUDA_CPU,PETSC_CUDA_BOTH} PetscCUDAFlag;
787: #endif

789: #if defined(PETSC_HAVE_CUSP) || defined(PETSC_HAVE_VECCUDA)
790: PETSC_EXTERN cublasHandle_t cublasv2handle;
791: #endif

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

795: #define REDUCE_SUM  0
796: #define REDUCE_MAX  1
797: #define REDUCE_MIN  2

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

813: PETSC_EXTERN PetscErrorCode PetscSplitReductionGet(MPI_Comm,PetscSplitReduction**);
814: PETSC_EXTERN PetscErrorCode PetscSplitReductionEnd(PetscSplitReduction*);
815: PETSC_EXTERN PetscErrorCode PetscSplitReductionExtend(PetscSplitReduction*);

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

850: #include <omp.h>
851: typedef omp_lock_t PetscSpinlock;
852: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockCreate(PetscSpinlock *omp_lock)
853: {
854:   omp_init_lock(omp_lock);
855:   return 0;
856: }
857: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockLock(PetscSpinlock *omp_lock)
858: {
859:   omp_set_lock(omp_lock);
860:   return 0;
861: }
862: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockUnlock(PetscSpinlock *omp_lock)
863: {
864:   omp_unset_lock(omp_lock);
865:   return 0;
866: }
867: PETSC_STATIC_INLINE PetscErrorCode PetscSpinlockDestroy(PetscSpinlock *omp_lock)
868: {
869:   omp_destroy_lock(omp_lock);
870:   return 0;
871: }
872: #else
873: Thread safety requires either --with-openmp or --download-concurrencykit
874: #endif

876: #else
877: typedef int PetscSpinlock;
878: #define PetscSpinlockCreate(a)  0
879: #define PetscSpinlockLock(a)    0
880: #define PetscSpinlockUnlock(a)  0
881: #define PetscSpinlockDestroy(a) 0
882: #endif

884: #if defined(PETSC_HAVE_THREADSAFETY)
885: extern PetscSpinlock PetscViewerASCIISpinLockOpen;
886: extern PetscSpinlock PetscViewerASCIISpinLockStdout;
887: extern PetscSpinlock PetscViewerASCIISpinLockStderr;
888: extern PetscSpinlock PetscCommSpinLock;
889: #endif
890: #endif

892: #endif /* _PETSCHEAD_H */