Actual source code: petscimpl.h

petsc-3.5.4 2015-05-23
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: /*
 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,memchildren;
 71:   PetscObjectId        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:   PetscObjectId        parentid;
 82:   char*                name;
 83:   char                 *prefix;
 84:   PetscInt             tablevel;
 85:   void                 *cpp;
 86:   PetscObjectState     state;
 87:   PetscInt             int_idmax,        intstar_idmax;
 88:   PetscObjectState     *intcomposedstate,*intstarcomposedstate;
 89:   PetscInt             *intcomposeddata, **intstarcomposeddata;
 90:   PetscInt             real_idmax,        realstar_idmax;
 91:   PetscObjectState     *realcomposedstate,*realstarcomposedstate;
 92:   PetscReal            *realcomposeddata, **realstarcomposeddata;
 93:   PetscInt             scalar_idmax,        scalarstar_idmax;
 94:   PetscObjectState     *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_SAWS)
110:   PetscBool            amsmem;          /* if PETSC_TRUE then this object is registered with SAWs and visible to clients */
111:   PetscBool            amspublishblock; /* if PETSC_TRUE and publishing objects then will block at PetscObjectSAWsBlock() */
112: #endif
113: } _p_PetscObject;

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

119: #define  PETSCFREEDHEADER -1

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

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

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

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

139:     Level: developer

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

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

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

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

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

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

163:     Level: developer

165:    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

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

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

179: PETSC_INTERN PetscErrorCode PetscCitationsInitialize(void);
180: PETSC_INTERN PetscErrorCode 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)


198: #else

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

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

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

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

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

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

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

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

253: #endif

255: #if !defined(PETSC_USE_DEBUG)


267: #else

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

291:   do {                                                  \
294:   } while (0)

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

306:   do {                                                                  \
307:     PetscErrorCode _7_ierr;                                             \
308:     PetscReal b1[2],b2[2];                                              \
309:     b1[0] = -b; b1[1] = b;                                              \
310:     _7_MPI_Allreduce(b1,b2,2,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));CHKERRQ(_7_ierr); \
311:     if (-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_MPI_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)

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

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

341: #endif

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

347:    Level: developer

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

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

361:    Level: developer

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

372: /*MC
373:    PetscObjectStateIncrease - Increases the state of any PetscObject

375:    Synopsis:
376:    #include "petsc-private/petscimpl.h"
377:    PetscErrorCode PetscObjectStateIncrease(PetscObject obj)

379:    Logically Collective

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

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

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

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

398:    Level: developer

400:    seealso: PetscObjectStateGet()

402:    Concepts: state

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

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

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

424:    Not collective

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

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

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

440: /*MC
441:    PetscObjectComposedDataGetInt - retrieve integer data attached to an object

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

447:    Not collective

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

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

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

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

465: /*MC
466:    PetscObjectComposedDataSetIntstar - attach integer array data to a PetscObject

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

472:    Not collective

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

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

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

489: /*MC
490:    PetscObjectComposedDataGetIntstar - retrieve integer array data
491:    attached to an object

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

497:    Not collective

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

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

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

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

515: /*MC
516:    PetscObjectComposedDataSetReal - attach real data to a PetscObject

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

522:    Not collective

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

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

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

539: /*MC
540:    PetscObjectComposedDataGetReal - retrieve real data attached to an object

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

546:    Not collective

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

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

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

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

564: /*MC
565:    PetscObjectComposedDataSetRealstar - attach real array data to a PetscObject

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

571:    Not collective

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

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

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

588: /*MC
589:    PetscObjectComposedDataGetRealstar - retrieve real array data
590:    attached to an object

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

596:    Not collective

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

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

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

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

614: /*MC
615:    PetscObjectComposedDataSetScalar - attach scalar data to a PetscObject

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

621:    Not collective

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

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

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

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

649:    Not collective

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

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

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

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

672: /*MC
673:    PetscObjectComposedDataSetScalarstar - attach scalar array data to a PetscObject

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

679:    Not collective

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

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

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

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

708:    Not collective

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

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

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

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

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

733: PETSC_EXTERN PetscMPIInt Petsc_Counter_keyval;
734: PETSC_EXTERN PetscMPIInt Petsc_InnerComm_keyval;
735: PETSC_EXTERN PetscMPIInt Petsc_OuterComm_keyval;

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

747: #if defined(PETSC_HAVE_CUSP)
748: /*E
749:     PetscCUSPFlag - indicates which memory (CPU, GPU, or none contains valid vector

751:    PETSC_CUSP_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
752:    PETSC_CUSP_GPU - GPU has valid vector/matrix entries
753:    PETSC_CUSP_CPU - CPU has valid vector/matrix entries
754:    PETSC_CUSP_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

756:    Level: developer
757: E*/
758: typedef enum {PETSC_CUSP_UNALLOCATED,PETSC_CUSP_GPU,PETSC_CUSP_CPU,PETSC_CUSP_BOTH} PetscCUSPFlag;
759: #endif

761: #if defined(PETSC_HAVE_VIENNACL)
762: /*E
763:     PetscViennaCLFlag - indicates which memory (CPU, GPU, or none contains valid vector

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

770:    Level: developer
771: E*/
772: typedef enum {PETSC_VIENNACL_UNALLOCATED,PETSC_VIENNACL_GPU,PETSC_VIENNACL_CPU,PETSC_VIENNACL_BOTH} PetscViennaCLFlag;
773: #endif

775: #endif /* _PETSCHEAD_H */