Actual source code: petscimpl.h

petsc-3.4.5 2014-06-29
  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: /*
 11:    All major PETSc data structures have a common core; this is defined
 12:    below by PETSCHEADER.

 14:    PetscHeaderCreate() should be used whenever creating a PETSc structure.
 15: */

 17: /*
 18:    PetscOps: structure of core operations that all PETSc objects support.

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

 37: typedef struct {
 38:    PetscErrorCode (*getcomm)(PetscObject,MPI_Comm *);
 39:    PetscErrorCode (*view)(PetscObject,PetscViewer);
 40:    PetscErrorCode (*destroy)(PetscObject*);
 41:    PetscErrorCode (*compose)(PetscObject,const char[],PetscObject);
 42:    PetscErrorCode (*query)(PetscObject,const char[],PetscObject *);
 43:    PetscErrorCode (*composefunction)(PetscObject,const char[],void (*)(void));
 44:    PetscErrorCode (*queryfunction)(PetscObject,const char[],void (**)(void));
 45: } PetscOps;

 47: typedef enum {PETSC_FORTRAN_CALLBACK_CLASS,PETSC_FORTRAN_CALLBACK_SUBTYPE,PETSC_FORTRAN_CALLBACK_MAXTYPE} PetscFortranCallbackType;
 48: typedef int PetscFortranCallbackId;
 49: #define PETSC_SMALLEST_FORTRAN_CALLBACK ((PetscFortranCallbackId)1000)
 50: PETSC_EXTERN PetscErrorCode PetscFortranCallbackRegister(PetscClassId,const char*,PetscFortranCallbackId*);
 51: PETSC_EXTERN PetscErrorCode PetscFortranCallbackGetSizes(PetscClassId,PetscInt*,PetscInt*);

 53: typedef struct {
 54:   void (*func)(void);
 55:   void *ctx;
 56: } PetscFortranCallback;

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

103:   PetscInt             noptionhandler;
104:   PetscErrorCode       (*optionhandler[PETSC_MAX_OPTIONS_HANDLER])(PetscObject,void*);
105:   PetscErrorCode       (*optiondestroy[PETSC_MAX_OPTIONS_HANDLER])(PetscObject,void*);
106:   void                 *optionctx[PETSC_MAX_OPTIONS_HANDLER];
107:   PetscPrecision       precision;
108:   PetscBool            optionsprinted;
109: #if defined(PETSC_HAVE_AMS)
110:   PetscInt             amsmem;
111:   PetscBool            amspublishblock; /* if PETSC_TRUE and publishing objects then will block at PetscObjectAMSBlock() */
112:   PetscBool            amsblock;
113: #endif
114: } _p_PetscObject;

116: #define PETSCHEADER(ObjectOps) \
117:   _p_PetscObject hdr;          \
118:   ObjectOps      *ops

120: #define  PETSCFREEDHEADER -1

122: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscObjectFunction)(PetscObject*); /* force cast in next macro to NEVER use extern "C" style */
123: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscObjectViewerFunction)(PetscObject,PetscViewer);

125: /*@C
126:     PetscHeaderCreate - Creates a PETSc object of a particular class, indicated by tp

128:     Input Parameters:
129: +   tp - the data structure type of the object (for example _p_Vec)
130: .   pops - the data structure type of the objects operations (for example VecOps)
131: .   classid - the classid associated with this object (for example VEC_CLASSID)
132: .   class_name - string name of class; should be static (for example "Vec")
133: .   com - the MPI Communicator
134: .   des - the destroy routine for this object (for example VecDestroy())
135: -   vie - the view routine for this object (for example VecView())

137:     Output Parameter:
138: .   h - the newly created object

140:     Level: developer

142:    Developer Note: This currently is a CPP macro because it takes the types (for example _p_Vec and VecOps) as arguments

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

146: @*/
147: #define PetscHeaderCreate(h,tp,pops,classid,class_name,descr,mansec,com,des,vie) \
148:   (PetscNew(struct tp,&(h)) ||                                                  \
149:    PetscNew(PetscOps,&(((PetscObject)(h))->bops)) ||                            \
150:    PetscNew(pops,&((h)->ops)) ||                                                \
151:    PetscHeaderCreate_Private((PetscObject)h,classid,class_name,descr,mansec,com,(PetscObjectFunction)des,(PetscObjectViewerFunction)vie) || \
152:    PetscLogObjectCreate(h) ||                                                   \
153:    PetscLogObjectMemory(h, sizeof(struct tp) + sizeof(PetscOps) + sizeof(pops)))

155: PETSC_EXTERN PetscErrorCode PetscComposedQuantitiesDestroy(PetscObject obj);
156: PETSC_EXTERN PetscErrorCode PetscHeaderCreate_Private(PetscObject,PetscClassId,const char[],const char[],const char[],MPI_Comm,PetscErrorCode (*)(PetscObject*),PetscErrorCode (*)(PetscObject,PetscViewer));

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

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

164:     Level: developer

166:    Developer Note: This currently is a CPP macro because it accesses (*h)->ops which is a field in the derived class but not the PetscObject base class

168: .seealso: PetscHeaderCreate()
169: @*/
170: #define PetscHeaderDestroy(h)                         \
171:   (PetscHeaderDestroy_Private((PetscObject)(*h)) ||   \
172:    PetscFree((*h)->ops) ||                            \
173:    PetscFree(*h))

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

180: PETSC_INTERN PetscErrorCode PetscOptionsFindPair_Private(const char[],const char[],char**,PetscBool*);


184: /*
185:     Macros to test if a PETSc object is valid and if pointers are valid
186: */
187: #if !defined(PETSC_USE_DEBUG)


197: #else

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

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

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

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

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

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

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

247: #endif

249: #if !defined(PETSC_USE_DEBUG)


261: #else

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

285:   do {                                                  \
288:   } while (0)

291:   do {                                                                  \
292:     PetscErrorCode _7_ierr;                                             \
293:     PetscReal b1[2],b2[2];                                              \
294:     b1[0] = -PetscRealPart(b); b1[1] = PetscRealPart(b);                \
295:     _7_MPI_Allreduce(b1,b2,2,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
296:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Scalar value must be same on all processes, argument # %d",c); \
297:   } while (0)

300:   do {                                                                  \
301:     PetscErrorCode _7_ierr;                                             \
302:     PetscReal b1[2],b2[2];                                              \
303:     b1[0] = -b; b1[1] = b;                                              \
304:     _7_MPI_Allreduce(b1,b2,2,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
305:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Real value must be same on all processes, argument # %d",c); \
306:   } while (0)

309:   do {                                                                  \
310:     PetscErrorCode _7_ierr;                                             \
311:     PetscInt b1[2],b2[2];                                               \
312:     b1[0] = -b; b1[1] = b;                                              \
313:     _7_MPI_Allreduce(b1,b2,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
314:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Int value must be same on all processes, argument # %d",c); \
315:   } while (0)

318:   do {                                                                  \
319:     PetscErrorCode _7_ierr;                                             \
320:     PetscMPIInt b1[2],b2[2];                                            \
321:     b1[0] = -(PetscMPIInt)b; b1[1] = (PetscMPIInt)b;                    \
322:     _7_MPI_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
323:     if (-b2[0] != b2[1]) SETERRQ1(PetscObjectComm((PetscObject)a),PETSC_ERR_ARG_WRONG,"Bool value must be same on all processes, argument # %d",c); \
324:   } while (0)

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

335: #endif

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

341:    Level: developer

343: .seealso: PetscUseMethod()
344: */
345: #define  PetscTryMethod(obj,A,B,C) \
346:   0;{ PetscErrorCode (*f)B, __ierr; \
347:     __PetscObjectQueryFunction((PetscObject)obj,A,&f);CHKERRQ(__ierr); \
348:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
349:   }

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

355:    Level: developer

357: .seealso: PetscTryMethod()
358: */
359: #define  PetscUseMethod(obj,A,B,C) \
360:   0;{ PetscErrorCode (*f)B, __ierr; \
361:     __PetscObjectQueryFunction((PetscObject)obj,A,&f);CHKERRQ(__ierr); \
362:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
363:     else SETERRQ1(PetscObjectComm((PetscObject)obj),PETSC_ERR_SUP,"Cannot locate function %s in object",A); \
364:   }

366: /*MC
367:    PetscObjectStateIncrease - Increases the state of any PetscObject,
368:    regardless of the type.

370:    Synopsis:
371:    #include "petscsys.h"
372:    PetscErrorCode PetscObjectStateIncrease(PetscObject obj)

374:    Not Collective

376:    Input Parameter:
377: .  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
378:          cast with a (PetscObject), for example,
379:          PetscObjectStateIncrease((PetscObject)mat);

381:    Notes: object state is an integer which gets increased every time
382:    the object is changed. By saving and later querying the object state
383:    one can determine whether information about the object is still current.
384:    Currently, state is maintained for Vec and Mat objects.

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

391:    Level: developer

393:    seealso: PetscObjectStateQuery(), PetscObjectStateDecrease()

395:    Concepts: state

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

400: /*MC
401:    PetscObjectStateDecrease - Decreases the state of any PetscObject,
402:    regardless of the type.

404:    Synopsis:
405:    #include "petscsys.h"
406:    PetscErrorCode PetscObjectStateDecrease(PetscObject obj)

408:    Not Collective

410:    Input Parameter:
411: .  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
412:          cast with a (PetscObject), for example,
413:          PetscObjectStateIncrease((PetscObject)mat);

415:    Notes: object state is an integer which gets increased every time
416:    the object is changed. By saving and later querying the object state
417:    one can determine whether information about the object is still current.
418:    Currently, state is maintained for Vec and Mat objects.

420:    Level: developer

422:    seealso: PetscObjectStateQuery(), PetscObjectStateIncrease()

424:    Concepts: state

426: M*/
427: #define PetscObjectStateDecrease(obj) ((obj)->state--,0)

429: PETSC_EXTERN PetscErrorCode PetscObjectStateQuery(PetscObject,PetscInt*);
430: PETSC_EXTERN PetscErrorCode PetscObjectSetState(PetscObject,PetscInt);
431: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataRegister(PetscInt*);
432: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseInt(PetscObject);
433: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseIntstar(PetscObject);
434: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseReal(PetscObject);
435: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseRealstar(PetscObject);
436: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalar(PetscObject);
437: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalarstar(PetscObject);
438: PETSC_EXTERN PetscInt         PetscObjectComposedDataMax;
439: /*MC
440:    PetscObjectComposedDataSetInt - attach integer data to a PetscObject

442:    Synopsis:
443:    #include "petscsys.h"
444:    PetscErrorCode PetscObjectComposedDataSetInt(PetscObject obj,int id,int data)

446:    Not collective

448:    Input parameters:
449: +  obj - the object to which data is to be attached
450: .  id - the identifier for the data
451: -  data - the data to  be attached

453:    Notes
454:    The data identifier can best be determined through a call to
455:    PetscObjectComposedDataRegister()

457:    Level: developer
458: M*/
459: #define PetscObjectComposedDataSetInt(obj,id,data)                                      \
460:   ((((obj)->int_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseInt(obj)) ||  \
461:    ((obj)->intcomposeddata[id] = data,(obj)->intcomposedstate[id] = (obj)->state, 0))

463: /*MC
464:    PetscObjectComposedDataGetInt - retrieve integer data attached to an object

466:    Synopsis:
467:    #include "petscsys.h"
468:    PetscErrorCode PetscObjectComposedDataGetInt(PetscObject obj,int id,int data,PetscBool  flag)

470:    Not collective

472:    Input parameters:
473: +  obj - the object from which data is to be retrieved
474: -  id - the identifier for the data

476:    Output parameters
477: +  data - the data to be retrieved
478: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

482:    Level: developer
483: M*/
484: #define PetscObjectComposedDataGetInt(obj,id,data,flag)                            \
485:   ((((obj)->intcomposedstate && ((obj)->intcomposedstate[id] == (obj)->state)) ?   \
486:    (data = (obj)->intcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

488: /*MC
489:    PetscObjectComposedDataSetIntstar - attach integer array data to a PetscObject

491:    Synopsis:
492:    #include "petscsys.h"
493:    PetscErrorCode PetscObjectComposedDataSetIntstar(PetscObject obj,int id,int *data)

495:    Not collective

497:    Input parameters:
498: +  obj - the object to which data is to be attached
499: .  id - the identifier for the data
500: -  data - the data to  be attached

502:    Notes
503:    The data identifier can best be determined through a call to
504:    PetscObjectComposedDataRegister()

506:    Level: developer
507: M*/
508: #define PetscObjectComposedDataSetIntstar(obj,id,data)                                          \
509:   ((((obj)->intstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseIntstar(obj)) ||  \
510:    ((obj)->intstarcomposeddata[id] = data,(obj)->intstarcomposedstate[id] = (obj)->state, 0))

512: /*MC
513:    PetscObjectComposedDataGetIntstar - retrieve integer array data
514:    attached to an object

516:    Synopsis:
517:    #include "petscsys.h"
518:    PetscErrorCode PetscObjectComposedDataGetIntstar(PetscObject obj,int id,int *data,PetscBool  flag)

520:    Not collective

522:    Input parameters:
523: +  obj - the object from which data is to be retrieved
524: -  id - the identifier for the data

526:    Output parameters
527: +  data - the data to be retrieved
528: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

532:    Level: developer
533: M*/
534: #define PetscObjectComposedDataGetIntstar(obj,id,data,flag)                               \
535:   ((((obj)->intstarcomposedstate && ((obj)->intstarcomposedstate[id] == (obj)->state)) ?  \
536:    (data = (obj)->intstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

538: /*MC
539:    PetscObjectComposedDataSetReal - attach real data to a PetscObject

541:    Synopsis:
542:    #include "petscsys.h"
543:    PetscErrorCode PetscObjectComposedDataSetReal(PetscObject obj,int id,PetscReal data)

545:    Not collective

547:    Input parameters:
548: +  obj - the object to which data is to be attached
549: .  id - the identifier for the data
550: -  data - the data to  be attached

552:    Notes
553:    The data identifier can best be determined through a call to
554:    PetscObjectComposedDataRegister()

556:    Level: developer
557: M*/
558: #define PetscObjectComposedDataSetReal(obj,id,data)                                       \
559:   ((((obj)->real_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseReal(obj)) ||  \
560:    ((obj)->realcomposeddata[id] = data,(obj)->realcomposedstate[id] = (obj)->state, 0))

562: /*MC
563:    PetscObjectComposedDataGetReal - retrieve real data attached to an object

565:    Synopsis:
566:    #include "petscsys.h"
567:    PetscErrorCode PetscObjectComposedDataGetReal(PetscObject obj,int id,PetscReal data,PetscBool  flag)

569:    Not collective

571:    Input parameters:
572: +  obj - the object from which data is to be retrieved
573: -  id - the identifier for the data

575:    Output parameters
576: +  data - the data to be retrieved
577: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

581:    Level: developer
582: M*/
583: #define PetscObjectComposedDataGetReal(obj,id,data,flag)                            \
584:   ((((obj)->realcomposedstate && ((obj)->realcomposedstate[id] == (obj)->state)) ?  \
585:    (data = (obj)->realcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

587: /*MC
588:    PetscObjectComposedDataSetRealstar - attach real array data to a PetscObject

590:    Synopsis:
591:    #include "petscsys.h"
592:    PetscErrorCode PetscObjectComposedDataSetRealstar(PetscObject obj,int id,PetscReal *data)

594:    Not collective

596:    Input parameters:
597: +  obj - the object to which data is to be attached
598: .  id - the identifier for the data
599: -  data - the data to  be attached

601:    Notes
602:    The data identifier can best be determined through a call to
603:    PetscObjectComposedDataRegister()

605:    Level: developer
606: M*/
607: #define PetscObjectComposedDataSetRealstar(obj,id,data)                                           \
608:   ((((obj)->realstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseRealstar(obj)) ||  \
609:    ((obj)->realstarcomposeddata[id] = data, (obj)->realstarcomposedstate[id] = (obj)->state, 0))

611: /*MC
612:    PetscObjectComposedDataGetRealstar - retrieve real array data
613:    attached to an object

615:    Synopsis:
616:    #include "petscsys.h"
617:    PetscErrorCode PetscObjectComposedDataGetRealstar(PetscObject obj,int id,PetscReal *data,PetscBool  flag)

619:    Not collective

621:    Input parameters:
622: +  obj - the object from which data is to be retrieved
623: -  id - the identifier for the data

625:    Output parameters
626: +  data - the data to be retrieved
627: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

631:    Level: developer
632: M*/
633: #define PetscObjectComposedDataGetRealstar(obj,id,data,flag)                                \
634:   ((((obj)->realstarcomposedstate && ((obj)->realstarcomposedstate[id] == (obj)->state)) ?  \
635:    (data = (obj)->realstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

637: /*MC
638:    PetscObjectComposedDataSetScalar - attach scalar data to a PetscObject

640:    Synopsis:
641:    #include "petscsys.h"
642:    PetscErrorCode PetscObjectComposedDataSetScalar(PetscObject obj,int id,PetscScalar data)

644:    Not collective

646:    Input parameters:
647: +  obj - the object to which data is to be attached
648: .  id - the identifier for the data
649: -  data - the data to  be attached

651:    Notes
652:    The data identifier can best be determined through a call to
653:    PetscObjectComposedDataRegister()

655:    Level: developer
656: M*/
657: #if defined(PETSC_USE_COMPLEX)
658: #define PetscObjectComposedDataSetScalar(obj,id,data)                                        \
659:   ((((obj)->scalar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalar(obj)) || \
660:    ((obj)->scalarcomposeddata[id] = data,(obj)->scalarcomposedstate[id] = (obj)->state, 0))
661: #else
662: #define PetscObjectComposedDataSetScalar(obj,id,data) \
663:         PetscObjectComposedDataSetReal(obj,id,data)
664: #endif
665: /*MC
666:    PetscObjectComposedDataGetScalar - retrieve scalar data attached to an object

668:    Synopsis:
669:    #include "petscsys.h"
670:    PetscErrorCode PetscObjectComposedDataGetScalar(PetscObject obj,int id,PetscScalar 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: #if defined(PETSC_USE_COMPLEX)
687: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                              \
688:   ((((obj)->scalarcomposedstate && ((obj)->scalarcomposedstate[id] == (obj)->state) ) ? \
689:    (data = (obj)->scalarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
690: #else
691: #define PetscObjectComposedDataGetScalar(obj,id,data,flag)                             \
692:         PetscObjectComposedDataGetReal(obj,id,data,flag)
693: #endif

695: /*MC
696:    PetscObjectComposedDataSetScalarstar - attach scalar array data to a PetscObject

698:    Synopsis:
699:    #include "petscsys.h"
700:    PetscErrorCode PetscObjectComposedDataSetScalarstar(PetscObject obj,int id,PetscScalar *data)

702:    Not collective

704:    Input parameters:
705: +  obj - the object to which data is to be attached
706: .  id - the identifier for the data
707: -  data - the data to  be attached

709:    Notes
710:    The data identifier can best be determined through a call to
711:    PetscObjectComposedDataRegister()

713:    Level: developer
714: M*/
715: #if defined(PETSC_USE_COMPLEX)
716: #define PetscObjectComposedDataSetScalarstar(obj,id,data)                                             \
717:   ((((obj)->scalarstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseScalarstar(obj)) ||  \
718:    ((obj)->scalarstarcomposeddata[id] = data,(obj)->scalarstarcomposedstate[id] = (obj)->state, 0))
719: #else
720: #define PetscObjectComposedDataSetScalarstar(obj,id,data) \
721:         PetscObjectComposedDataSetRealstar(obj,id,data)
722: #endif
723: /*MC
724:    PetscObjectComposedDataGetScalarstar - retrieve scalar array data
725:    attached to an object

727:    Synopsis:
728:    #include "petscsys.h"
729:    PetscErrorCode PetscObjectComposedDataGetScalarstar(PetscObject obj,int id,PetscScalar *data,PetscBool  flag)

731:    Not collective

733:    Input parameters:
734: +  obj - the object from which data is to be retrieved
735: -  id - the identifier for the data

737:    Output parameters
738: +  data - the data to be retrieved
739: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

743:    Level: developer
744: M*/
745: #if defined(PETSC_USE_COMPLEX)
746: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)                                 \
747:   ((((obj)->scalarstarcomposedstate && ((obj)->scalarstarcomposedstate[id] == (obj)->state)) ? \
748:        (data = (obj)->scalarstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
749: #else
750: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)         \
751:         PetscObjectComposedDataGetRealstar(obj,id,data,flag)
752: #endif

754: PETSC_EXTERN PetscMPIInt Petsc_Counter_keyval;
755: PETSC_EXTERN PetscMPIInt Petsc_InnerComm_keyval;
756: PETSC_EXTERN PetscMPIInt Petsc_OuterComm_keyval;

758: /*
759:   PETSc communicators have this attribute, see
760:   PetscCommDuplicate(), PetscCommDestroy(), PetscCommGetNewTag(), PetscObjectGetName()
761: */
762: typedef struct {
763:   PetscMPIInt tag;              /* next free tag value */
764:   PetscInt    refcount;         /* number of references, communicator can be freed when this reaches 0 */
765:   PetscInt    namecount;        /* used to generate the next name, as in Vec_0, Mat_1, ... */
766: } PetscCommCounter;

768: #if defined(PETSC_HAVE_CUSP)
769: /*E
770:     PetscCUSPFlag - indicates which memory (CPU, GPU, or none contains valid vector

772:    PETSC_CUSP_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
773:    PETSC_CUSP_GPU - GPU has valid vector/matrix entries
774:    PETSC_CUSP_CPU - CPU has valid vector/matrix entries
775:    PETSC_CUSP_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

777:    Level: developer
778: E*/
779: typedef enum {PETSC_CUSP_UNALLOCATED,PETSC_CUSP_GPU,PETSC_CUSP_CPU,PETSC_CUSP_BOTH} PetscCUSPFlag;
780: #endif

782: #endif /* _PETSCHEAD_H */