Actual source code: petscimpl.h

petsc-3.3-p7 2013-05-11
  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.
 19:    
 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:       destroy()         - Is the routine for destroying the entire PETSc object; 
 24:                           for example,MatDestroy() is the general matrix 
 25:                           destruction routine.
 26:       compose()         - Associates a PETSc object with another PETSc object.
 27:       query()           - Returns a different PETSc object that has been associated
 28:                           with the first object.
 29:       composefunction() - Attaches an additional registered function.
 30:       queryfunction()   - Requests a registered function that has been registered.
 31:       publish()         - Not currently used
 32: */

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

 45: /*
 46:    All PETSc objects begin with the fields defined in PETSCHEADER.
 47:    The PetscObject is a way of examining these fields regardless of 
 48:    the specific object. In C++ this could be a base abstract class
 49:    from which all objects are derived.
 50: */
 51: #define PETSC_MAX_OPTIONS_HANDLER 5
 52: typedef struct _p_PetscObject {
 53:   PetscClassId   classid;
 54:   PetscOps       *bops;
 55:   MPI_Comm       comm;
 56:   PetscInt       type;
 57:   PetscLogDouble flops,time,mem;
 58:   PetscInt       id;
 59:   PetscInt       refct;
 60:   PetscMPIInt    tag;
 61:   PetscFList     qlist;
 62:   PetscOList     olist;
 63:   char           *class_name;
 64:   char           *description;
 65:   char           *mansec;
 66:   char           *type_name;
 67:   PetscObject    parent;
 68:   PetscInt       parentid;
 69:   char*          name;
 70:   char           *prefix;
 71:   PetscInt       tablevel;
 72:   void           *cpp;
 73:   PetscInt       amem;
 74:   PetscInt       state;
 75:   PetscInt       int_idmax,        intstar_idmax;
 76:   PetscInt       *intcomposedstate,*intstarcomposedstate;
 77:   PetscInt       *intcomposeddata, **intstarcomposeddata;
 78:   PetscInt       real_idmax,        realstar_idmax;
 79:   PetscInt       *realcomposedstate,*realstarcomposedstate;
 80:   PetscReal      *realcomposeddata, **realstarcomposeddata;
 81:   PetscInt       scalar_idmax,        scalarstar_idmax;
 82:   PetscInt       *scalarcomposedstate,*scalarstarcomposedstate;
 83:   PetscScalar    *scalarcomposeddata, **scalarstarcomposeddata;
 84:   void           (**fortran_func_pointers)(void);                  /* used by Fortran interface functions to stash user provided Fortran functions */
 85:   PetscInt       num_fortran_func_pointers;                        /* number of Fortran function pointers allocated */
 86:   void           *python_context;
 87:   PetscErrorCode (*python_destroy)(void*);

 89:   PetscInt       noptionhandler;
 90:   PetscErrorCode (*optionhandler[PETSC_MAX_OPTIONS_HANDLER])(PetscObject,void*);
 91:   PetscErrorCode (*optiondestroy[PETSC_MAX_OPTIONS_HANDLER])(PetscObject,void*);
 92:   void           *optionctx[PETSC_MAX_OPTIONS_HANDLER];
 93:   PetscPrecision precision;
 94:   PetscBool      optionsprinted;
 95: } _p_PetscObject;

 97: #define PETSCHEADER(ObjectOps) \
 98:   _p_PetscObject hdr;               \
 99:   ObjectOps      *ops

101: #define  PETSCFREEDHEADER -1

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

106: /*@C
107:     PetscHeaderCreate - Creates a PETSc object

109:     Input Parameters:
110: +   tp - the data structure type of the object
111: .   pops - the data structure type of the objects operations (for example VecOps)
112: .   cook - the classid associated with this object
113: .   t - type (no longer should be used)
114: .   class_name - string name of class; should be static
115: .   com - the MPI Communicator
116: .   des - the destroy routine for this object
117: -   vie - the view routine for this object

119:     Output Parameter:
120: .   h - the newly created object

122:     Level: developer

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

126: @*/
127: #define PetscHeaderCreate(h,tp,pops,cook,t,class_name,descr,mansec,com,des,vie) \
128:   (PetscNew(struct tp,&(h)) ||                                                \
129:    PetscNew(PetscOps,&(((PetscObject)(h))->bops)) ||                        \
130:    PetscNew(pops,&((h)->ops)) ||                                        \
131:    PetscHeaderCreate_Private((PetscObject)h,cook,t,class_name,descr,mansec,com,(PetscObjectFunction)des,(PetscObjectViewerFunction)vie) || \
132:    PetscLogObjectCreate(h) ||                                                \
133:    PetscLogObjectMemory(h, sizeof(struct tp) + sizeof(PetscOps) + sizeof(pops)))

135: PETSC_EXTERN PetscErrorCode PetscComposedQuantitiesDestroy(PetscObject obj);
136: PETSC_EXTERN PetscErrorCode PetscHeaderCreate_Private(PetscObject,PetscClassId,PetscInt,const char[],const char[],const char[],MPI_Comm,PetscErrorCode (*)(PetscObject*),PetscErrorCode (*)(PetscObject,PetscViewer));

138: /*@C
139:     PetscHeaderDestroy - Final step in destroying a PetscObject

141:     Input Parameters:
142: .   h - the header created with PetscHeaderCreate()

144:     Level: developer

146: .seealso: PetscHeaderCreate()
147: @*/
148: #define PetscHeaderDestroy(h)                           \
149:   (PetscLogObjectDestroy((PetscObject)(*h)) ||           \
150:    PetscComposedQuantitiesDestroy((PetscObject)*h) || \
151:    PetscHeaderDestroy_Private((PetscObject)(*h)) || \
152:    PetscFree((*h)->ops) ||                           \
153:    PetscFree(*h))

155: PETSC_EXTERN PetscErrorCode PetscHeaderDestroy_Private(PetscObject);
156: PETSC_EXTERN PetscErrorCode PetscObjectCopyFortranFunctionPointers(PetscObject,PetscObject);

158: /* ---------------------------------------------------------------------------------------*/

160: #if !defined(PETSC_USE_DEBUG)


169: #elif !defined(PETSC_HAVE_UNALIGNED_POINTERS)
170: /* 
171:     Macros to test if a PETSc object is valid and if pointers are
172: valid

174: */
176:   do {                                                                  \
177:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object: Parameter # %d",arg); \
178:     if ((unsigned long)(h) & (unsigned long)3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Invalid Pointer to Object: Parameter # %d",arg); \
179:     if (((PetscObject)(h))->classid != ck) {                            \
180:       if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free: Parameter # %d",arg); \
181:       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong type of object: Parameter # %d",arg); \
182:     }                                                                   \
183:   } while (0)

186:   do {                                                                  \
187:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object: Parameter # %d",arg); \
188:     if ((unsigned long)(h) & (unsigned long)3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Invalid Pointer to Object: Parameter # %d",arg); \
189:     else if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free: Parameter # %d",arg); \
190:     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); \
191:   } while (0)

194:   do {                                                                  \
195:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer: Parameter # %d",arg); \
196:     if ((unsigned long)(h) & (unsigned long)3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Invalid Pointer: Parameter # %d",arg); \
197:   } while (0)

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

203:   do {                                                                  \
204:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Null Pointer: Parameter # %d",arg); \
205:     if ((unsigned long)(h) & (unsigned long)3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Invalid Pointer to Int: Parameter # %d",arg); \
206:   } while (0)

208: #if !defined(PETSC_HAVE_DOUBLES_ALIGNED)
210:   do {                                                                  \
211:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer: Parameter # %d",arg); \
212:     if ((unsigned long)(h) & (unsigned long)3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Invalid Pointer to PetscScalar: Parameter # %d",arg); \
213:   } while (0)
214: #else
216:   do {                                                                  \
217:     if (!h) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer: Parameter # %d",arg); \
218:     if ((unsigned long)(h) & (unsigned long)7) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Invalid Pointer to PetscScalar: Parameter # %d",arg); \
219:   } while (0)
220: #endif

222: #else
223: /*
224:      Version where pointers don't have any particular alignment
225: */
227:   do {                                                                  \
228:     if (!h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object");  \
229:     if (((PetscObject)(h))->classid != ck) {                            \
230:       if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free"); \
231:       else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong Object");    \
232:     }                                                                   \
233:   } while (0)

236:   do {                                                                  \
237:     if (!h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Object");  \
238:     if (((PetscObject)(h))->classid == PETSCFREEDHEADER) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Object already free"); \
239:     else if (((PetscObject)(h))->classid < PETSC_SMALLEST_CLASSID || ((PetscObject)(h))->classid > PETSC_LARGEST_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Invalid type of object"); \
240:   } while (0)

243:   do {if (!h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer");} while (0)

246:   do {if (!h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer");} while (0)

249:   do {if (!h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer");} while (0)

251: #if !defined(PETSC_HAVE_DOUBLES_ALIGNED)
253:   do {if (!h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer");} while (0)
254: #else
256:   do {if (!h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Null Pointer");} while (0)
257: #endif

259: #endif

262: #if !defined(PETSC_USE_DEBUG)


274: #else

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

298:   do {                                                  \
301:   } while (0)

304:   do {                                                                  \
305:     PetscErrorCode _7_ierr;                                             \
306:     PetscReal b1[2],b2[2];                                              \
307:     b1[0] = -PetscRealPart(b); b1[1] = PetscRealPart(b);                \
308:     _7_MPI_Allreduce(b1,b2,2,MPIU_REAL,MPIU_MAX,((PetscObject)a)->comm);CHKERRQ(_7_ierr); \
309:     if (-b2[0] != b2[1]) SETERRQ1(((PetscObject)a)->comm,PETSC_ERR_ARG_WRONG,"Scalar value must be same on all processes, argument # %d",c); \
310:   } while (0)

313:   do {                                                                  \
314:     PetscErrorCode _7_ierr;                                             \
315:     PetscReal b1[2],b2[2];                                              \
316:     b1[0] = -b; b1[1] = b;                                              \
317:     _7_MPI_Allreduce(b1,b2,2,MPIU_REAL,MPIU_MAX,((PetscObject)a)->comm);CHKERRQ(_7_ierr); \
318:     if (-b2[0] != b2[1]) SETERRQ1(((PetscObject)a)->comm,PETSC_ERR_ARG_WRONG,"Real value must be same on all processes, argument # %d",c); \
319:   } while (0)

322:   do {                                                                  \
323:     PetscErrorCode _7_ierr;                                             \
324:     PetscInt b1[2],b2[2];                                               \
325:     b1[0] = -b; b1[1] = b;                                              \
326:     _7_MPI_Allreduce(b1,b2,2,MPIU_INT,MPI_MAX,((PetscObject)a)->comm);CHKERRQ(_7_ierr); \
327:     if (-b2[0] != b2[1]) SETERRQ1(((PetscObject)a)->comm,PETSC_ERR_ARG_WRONG,"Int value must be same on all processes, argument # %d",c); \
328:   } while (0)

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

340:   do {                                                                  \
341:     PetscErrorCode _7_ierr;                                             \
342:     PetscMPIInt b1[2],b2[2];                                            \
343:     b1[0] = -(PetscMPIInt)b; b1[1] = (PetscMPIInt)b;                    \
344:     _7_MPI_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,((PetscObject)a)->comm);CHKERRQ(_7_ierr); \
345:     if (-b2[0] != b2[1]) SETERRQ1(((PetscObject)a)->comm,PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes, argument # %d",c); \
346:   } while (0)

348: #endif

350: /*MC
351:    PetscObjectStateIncrease - Increases the state of any PetscObject, 
352:    regardless of the type.

354:    Synopsis:
355:    PetscErrorCode PetscObjectStateIncrease(PetscObject obj)

357:    Not Collective

359:    Input Parameter:
360: .  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
361:          cast with a (PetscObject), for example, 
362:          PetscObjectStateIncrease((PetscObject)mat);

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

369:    This routine is mostly for internal use by PETSc; a developer need only
370:    call it after explicit access to an object's internals. Routines such
371:    as VecSet or MatScale already call this routine. It is also called, as a 
372:    precaution, in VecRestoreArray, MatRestoreRow, MatRestoreArray.

374:    Level: developer

376:    seealso: PetscObjectStateQuery(), PetscObjectStateDecrease()

378:    Concepts: state

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

383: /*MC
384:    PetscObjectStateDecrease - Decreases the state of any PetscObject, 
385:    regardless of the type.

387:    Synopsis:
388:    PetscErrorCode PetscObjectStateDecrease(PetscObject obj)

390:    Not Collective

392:    Input Parameter:
393: .  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
394:          cast with a (PetscObject), for example, 
395:          PetscObjectStateIncrease((PetscObject)mat);

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

402:    Level: developer

404:    seealso: PetscObjectStateQuery(), PetscObjectStateIncrease()

406:    Concepts: state

408: M*/
409: #define PetscObjectStateDecrease(obj) ((obj)->state--,0)

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

424:    Synopsis:
425:    PetscErrorCode PetscObjectComposedDataSetInt(PetscObject obj,int id,int data)

427:    Not collective

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

434:    Notes
435:    The data identifier can best be determined through a call to
436:    PetscObjectComposedDataRegister()

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

444: /*MC
445:    PetscObjectComposedDataGetInt - retrieve integer data attached to an object

447:    Synopsis:
448:    PetscErrorCode PetscObjectComposedDataGetInt(PetscObject obj,int id,int data,PetscBool  flag)

450:    Not collective

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

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

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

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

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

471:    Synopsis:
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:    PetscErrorCode PetscObjectComposedDataGetIntstar(PetscObject obj,int id,int *data,PetscBool  flag)

498:    Not collective

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

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

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

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

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

519:    Synopsis:
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:    PetscErrorCode PetscObjectComposedDataGetReal(PetscObject obj,int id,PetscReal data,PetscBool  flag)

545:    Not collective

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

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

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

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

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

566:    Synopsis:
567:    PetscErrorCode PetscObjectComposedDataSetRealstar(PetscObject obj,int id,PetscReal *data)

569:    Not collective

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

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

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

586: /*MC
587:    PetscObjectComposedDataGetRealstar - retrieve real array data
588:    attached to an object

590:    Synopsis:
591:    PetscErrorCode PetscObjectComposedDataGetRealstar(PetscObject obj,int id,PetscReal *data,PetscBool  flag)

593:    Not collective

595:    Input parameters:
596: +  obj - the object from which data is to be retrieved
597: -  id - the identifier for the data

599:    Output parameters
600: +  data - the data to be retrieved
601: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

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

611: /*MC
612:    PetscObjectComposedDataSetScalar - attach scalar data to a PetscObject

614:    Synopsis:
615:    PetscErrorCode PetscObjectComposedDataSetScalar(PetscObject obj,int id,PetscScalar data)

617:    Not collective

619:    Input parameters:
620: +  obj - the object to which data is to be attached
621: .  id - the identifier for the data
622: -  data - the data to  be attached

624:    Notes
625:    The data identifier can best be determined through a call to
626:    PetscObjectComposedDataRegister()

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

641:    Synopsis:
642:    PetscErrorCode PetscObjectComposedDataGetScalar(PetscObject obj,int id,PetscScalar data,PetscBool  flag)

644:    Not collective

646:    Input parameters:
647: +  obj - the object from which data is to be retrieved
648: -  id - the identifier for the data

650:    Output parameters
651: +  data - the data to be retrieved
652: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

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

667: /*MC
668:    PetscObjectComposedDataSetScalarstar - attach scalar array data to a PetscObject 

670:    Synopsis:
671:    PetscErrorCode PetscObjectComposedDataSetScalarstar(PetscObject obj,int id,PetscScalar *data)

673:    Not collective

675:    Input parameters:
676: +  obj - the object to which data is to be attached
677: .  id - the identifier for the data
678: -  data - the data to  be attached

680:    Notes
681:    The data identifier can best be determined through a call to
682:    PetscObjectComposedDataRegister()

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

698:    Synopsis:
699:    PetscErrorCode PetscObjectComposedDataGetScalarstar(PetscObject obj,int id,PetscScalar *data,PetscBool  flag)

701:    Not collective

703:    Input parameters:
704: +  obj - the object from which data is to be retrieved
705: -  id - the identifier for the data

707:    Output parameters
708: +  data - the data to be retrieved
709: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

713:    Level: developer
714: M*/
715: #if defined(PETSC_USE_COMPLEX)
716: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)                                 \
717:   ((((obj)->scalarstarcomposedstate && ((obj)->scalarstarcomposedstate[id] == (obj)->state)) ? \
718:        (data = (obj)->scalarstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)
719: #else
720: #define PetscObjectComposedDataGetScalarstar(obj,id,data,flag)                 \
721:         PetscObjectComposedDataGetRealstar(obj,id,data,flag)
722: #endif

724: /* some vars for logging */
725: PETSC_EXTERN PetscBool PetscPreLoadingUsed;       /* true if we are or have done preloading */
726: PETSC_EXTERN PetscBool PetscPreLoadingOn;         /* true if we are currently in a preloading calculation */

728: PETSC_EXTERN PetscMPIInt Petsc_Counter_keyval;
729: PETSC_EXTERN PetscMPIInt Petsc_InnerComm_keyval;
730: PETSC_EXTERN PetscMPIInt Petsc_OuterComm_keyval;

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


743: #endif /* _PETSCHEAD_H */