Actual source code: petscimpl.h

petsc-3.6.4 2016-04-12
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[1];
 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[1]

119: #define  PETSCFREEDHEADER -1

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

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

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

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

139:     Level: developer

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

143: @*/
144: #define PetscHeaderCreate(h,classid,class_name,descr,mansec,comm,destroy,view) \
145:   (PetscNew(&(h)) || \
146:    PetscHeaderCreate_Private((PetscObject)h,classid,class_name,descr,mansec,comm,(PetscObjectDestroyFunction)destroy,(PetscObjectViewFunction)view) || \
147:    PetscLogObjectCreate(h) || \
148:    PetscLogObjectMemory((PetscObject)h,sizeof(*(h))))

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

153: /*@C
154:     PetscHeaderDestroy - Final step in destroying a PetscObject

156:     Input Parameters:
157: .   h - the header created with PetscHeaderCreate()

159:     Level: developer

161: .seealso: PetscHeaderCreate()
162: @*/
163: #define PetscHeaderDestroy(h) (PetscHeaderDestroy_Private((PetscObject)(*h)) || PetscFree(*h))

165: PETSC_EXTERN PetscErrorCode PetscHeaderDestroy_Private(PetscObject);
166: PETSC_EXTERN PetscErrorCode PetscObjectCopyFortranFunctionPointers(PetscObject,PetscObject);
167: PETSC_EXTERN PetscErrorCode PetscObjectSetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId*,void(*)(void),void *ctx);
168: PETSC_EXTERN PetscErrorCode PetscObjectGetFortranCallback(PetscObject,PetscFortranCallbackType,PetscFortranCallbackId,void(**)(void),void **ctx);

170: PETSC_INTERN PetscErrorCode PetscCitationsInitialize(void);
171: PETSC_INTERN PetscErrorCode PetscOptionsFindPair_Private(const char[],const char[],char**,PetscBool*);


175: /*
176:     Macros to test if a PETSc object is valid and if pointers are valid
177: */
178: #if !defined(PETSC_USE_DEBUG)


189: #else

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

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

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

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

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

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

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

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

244: #endif

246: #if !defined(PETSC_USE_DEBUG)


259: #else

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

283:   do {                                                  \
286:   } while (0)

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

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

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

369:    Synopsis:
370:    #include "petsc/private/petscimpl.h"
371:    PetscErrorCode PetscObjectStateIncrease(PetscObject obj)

373:    Logically Collective

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

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

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

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

392:    Level: developer

394:    seealso: PetscObjectStateGet()

396:    Concepts: state

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

401: PETSC_EXTERN PetscErrorCode PetscObjectStateGet(PetscObject,PetscObjectState*);
402: PETSC_EXTERN PetscErrorCode PetscObjectStateSet(PetscObject,PetscObjectState);
403: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataRegister(PetscInt*);
404: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseInt(PetscObject);
405: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseIntstar(PetscObject);
406: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseReal(PetscObject);
407: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseRealstar(PetscObject);
408: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalar(PetscObject);
409: PETSC_EXTERN PetscErrorCode PetscObjectComposedDataIncreaseScalarstar(PetscObject);
410: PETSC_EXTERN PetscInt         PetscObjectComposedDataMax;
411: /*MC
412:    PetscObjectComposedDataSetInt - attach integer data to a PetscObject

414:    Synopsis:
415:    #include "petsc/private/petscimpl.h"
416:    PetscErrorCode PetscObjectComposedDataSetInt(PetscObject obj,int id,int data)

418:    Not collective

420:    Input parameters:
421: +  obj - the object to which data is to be attached
422: .  id - the identifier for the data
423: -  data - the data to  be attached

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

428:    Level: developer
429: M*/
430: #define PetscObjectComposedDataSetInt(obj,id,data)                                      \
431:   ((((obj)->int_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseInt(obj)) ||  \
432:    ((obj)->intcomposeddata[id] = data,(obj)->intcomposedstate[id] = (obj)->state, 0))

434: /*MC
435:    PetscObjectComposedDataGetInt - retrieve integer data attached to an object

437:    Synopsis:
438:    #include "petsc/private/petscimpl.h"
439:    PetscErrorCode PetscObjectComposedDataGetInt(PetscObject obj,int id,int data,PetscBool  flag)

441:    Not collective

443:    Input parameters:
444: +  obj - the object from which data is to be retrieved
445: -  id - the identifier for the data

447:    Output parameters
448: +  data - the data to be retrieved
449: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

453:    Level: developer
454: M*/
455: #define PetscObjectComposedDataGetInt(obj,id,data,flag)                            \
456:   ((((obj)->intcomposedstate && ((obj)->intcomposedstate[id] == (obj)->state)) ?   \
457:    (data = (obj)->intcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

459: /*MC
460:    PetscObjectComposedDataSetIntstar - attach integer array data to a PetscObject

462:    Synopsis:
463:    #include "petsc/private/petscimpl.h"
464:    PetscErrorCode PetscObjectComposedDataSetIntstar(PetscObject obj,int id,int *data)

466:    Not collective

468:    Input parameters:
469: +  obj - the object to which data is to be attached
470: .  id - the identifier for the data
471: -  data - the data to  be attached

473:    Notes
474:    The data identifier can best be determined through a call to
475:    PetscObjectComposedDataRegister()

477:    Level: developer
478: M*/
479: #define PetscObjectComposedDataSetIntstar(obj,id,data)                                          \
480:   ((((obj)->intstar_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseIntstar(obj)) ||  \
481:    ((obj)->intstarcomposeddata[id] = data,(obj)->intstarcomposedstate[id] = (obj)->state, 0))

483: /*MC
484:    PetscObjectComposedDataGetIntstar - retrieve integer array data
485:    attached to an object

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

491:    Not collective

493:    Input parameters:
494: +  obj - the object from which data is to be retrieved
495: -  id - the identifier for the data

497:    Output parameters
498: +  data - the data to be retrieved
499: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

503:    Level: developer
504: M*/
505: #define PetscObjectComposedDataGetIntstar(obj,id,data,flag)                               \
506:   ((((obj)->intstarcomposedstate && ((obj)->intstarcomposedstate[id] == (obj)->state)) ?  \
507:    (data = (obj)->intstarcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

509: /*MC
510:    PetscObjectComposedDataSetReal - attach real data to a PetscObject

512:    Synopsis:
513:    #include "petsc/private/petscimpl.h"
514:    PetscErrorCode PetscObjectComposedDataSetReal(PetscObject obj,int id,PetscReal data)

516:    Not collective

518:    Input parameters:
519: +  obj - the object to which data is to be attached
520: .  id - the identifier for the data
521: -  data - the data to  be attached

523:    Notes
524:    The data identifier can best be determined through a call to
525:    PetscObjectComposedDataRegister()

527:    Level: developer
528: M*/
529: #define PetscObjectComposedDataSetReal(obj,id,data)                                       \
530:   ((((obj)->real_idmax < PetscObjectComposedDataMax) && PetscObjectComposedDataIncreaseReal(obj)) ||  \
531:    ((obj)->realcomposeddata[id] = data,(obj)->realcomposedstate[id] = (obj)->state, 0))

533: /*MC
534:    PetscObjectComposedDataGetReal - retrieve real data attached to an object

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

540:    Not collective

542:    Input parameters:
543: +  obj - the object from which data is to be retrieved
544: -  id - the identifier for the data

546:    Output parameters
547: +  data - the data to be retrieved
548: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

552:    Level: developer
553: M*/
554: #define PetscObjectComposedDataGetReal(obj,id,data,flag)                            \
555:   ((((obj)->realcomposedstate && ((obj)->realcomposedstate[id] == (obj)->state)) ?  \
556:    (data = (obj)->realcomposeddata[id],flag = PETSC_TRUE) : (flag = PETSC_FALSE)),0)

558: /*MC
559:    PetscObjectComposedDataSetRealstar - attach real array data to a PetscObject

561:    Synopsis:
562:    #include "petsc/private/petscimpl.h"
563:    PetscErrorCode PetscObjectComposedDataSetRealstar(PetscObject obj,int id,PetscReal *data)

565:    Not collective

567:    Input parameters:
568: +  obj - the object to which data is to be attached
569: .  id - the identifier for the data
570: -  data - the data to  be attached

572:    Notes
573:    The data identifier can best be determined through a call to
574:    PetscObjectComposedDataRegister()

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

582: /*MC
583:    PetscObjectComposedDataGetRealstar - retrieve real array data
584:    attached to an object

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

590:    Not collective

592:    Input parameters:
593: +  obj - the object from which data is to be retrieved
594: -  id - the identifier for the data

596:    Output parameters
597: +  data - the data to be retrieved
598: -  flag - PETSC_TRUE if the data item exists and is valid, PETSC_FALSE otherwise

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

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

608: /*MC
609:    PetscObjectComposedDataSetScalar - attach scalar data to a PetscObject

611:    Synopsis:
612:    #include "petsc/private/petscimpl.h"
613:    PetscErrorCode PetscObjectComposedDataSetScalar(PetscObject obj,int id,PetscScalar data)

615:    Not collective

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

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

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

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

643:    Not collective

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

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

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

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

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

669:    Synopsis:
670:    #include "petsc/private/petscimpl.h"
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:    #include "petsc/private/petscimpl.h"
700:    PetscErrorCode PetscObjectComposedDataGetScalarstar(PetscObject obj,int id,PetscScalar *data,PetscBool  flag)

702:    Not collective

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

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

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

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

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

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

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

741: #if defined(PETSC_HAVE_CUSP)
742: /*E
743:     PetscCUSPFlag - indicates which memory (CPU, GPU, or none contains valid vector

745:    PETSC_CUSP_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
746:    PETSC_CUSP_GPU - GPU has valid vector/matrix entries
747:    PETSC_CUSP_CPU - CPU has valid vector/matrix entries
748:    PETSC_CUSP_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

750:    Level: developer
751: E*/
752: typedef enum {PETSC_CUSP_UNALLOCATED,PETSC_CUSP_GPU,PETSC_CUSP_CPU,PETSC_CUSP_BOTH} PetscCUSPFlag;
753: #endif

755: #if defined(PETSC_HAVE_VIENNACL)
756: /*E
757:     PetscViennaCLFlag - indicates which memory (CPU, GPU, or none contains valid vector

759:    PETSC_VIENNACL_UNALLOCATED  - no memory contains valid matrix entries; NEVER used for vectors
760:    PETSC_VIENNACL_GPU - GPU has valid vector/matrix entries
761:    PETSC_VIENNACL_CPU - CPU has valid vector/matrix entries
762:    PETSC_VIENNACL_BOTH - Both GPU and CPU have valid vector/matrix entries and they match

764:    Level: developer
765: E*/
766: typedef enum {PETSC_VIENNACL_UNALLOCATED,PETSC_VIENNACL_GPU,PETSC_VIENNACL_CPU,PETSC_VIENNACL_BOTH} PetscViennaCLFlag;
767: #endif

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

771: #define REDUCE_SUM  0
772: #define REDUCE_MAX  1
773: #define REDUCE_MIN  2

775: typedef struct {
776:   MPI_Comm    comm;
777:   MPI_Request request;
778:   PetscBool   async;
779:   PetscScalar *lvalues;     /* this are the reduced values before call to MPI_Allreduce() */
780:   PetscScalar *gvalues;     /* values after call to MPI_Allreduce() */
781:   void        **invecs;     /* for debugging only, vector/memory used with each op */
782:   PetscInt    *reducetype;  /* is particular value to be summed or maxed? */
783:   SRState     state;        /* are we calling xxxBegin() or xxxEnd()? */
784:   PetscInt    maxops;       /* total amount of space we have for requests */
785:   PetscInt    numopsbegin;  /* number of requests that have been queued in */
786:   PetscInt    numopsend;    /* number of requests that have been gotten by user */
787: } PetscSplitReduction;

789: PETSC_EXTERN PetscErrorCode PetscSplitReductionGet(MPI_Comm,PetscSplitReduction**);
790: PETSC_EXTERN PetscErrorCode PetscSplitReductionEnd(PetscSplitReduction*);
791: PETSC_EXTERN PetscErrorCode PetscSplitReductionExtend(PetscSplitReduction*);

793: #endif /* _PETSCHEAD_H */