Actual source code: petscdmstag.h

  1: #ifndef PETSC_DMSTAG_H
  2: #define PETSC_DMSTAG_H

  4: #include <petscdm.h>
  5: #include <petscdmproduct.h>

  7: /* SUBMANSEC = DMStag */

  9: /*E
 10:   DMStagStencilLocation - enumerated type denoting a location relative to an element in a `DMSTAG` grid

 12:   The interpretation of these values is dimension-dependent.

 14:   Level: beginner

 16:   Developer Note:
 17:    The order of the enum entries is significant, as it corresponds to the canonical numbering
 18:    of DOFs, and the fact that the numbering starts at 0 may also be used by the implementation.

 20: .seealso: [](chapter_stag), `DMSTAG`, `DMStagStencil`, `DMStagGetLocationSlot()`, `DMStagStencilType`
 21: E*/
 22: typedef enum {
 23:   DMSTAG_NULL_LOCATION = 0,
 24:   DMSTAG_BACK_DOWN_LEFT,
 25:   DMSTAG_BACK_DOWN,
 26:   DMSTAG_BACK_DOWN_RIGHT,
 27:   DMSTAG_BACK_LEFT,
 28:   DMSTAG_BACK,
 29:   DMSTAG_BACK_RIGHT,
 30:   DMSTAG_BACK_UP_LEFT,
 31:   DMSTAG_BACK_UP,
 32:   DMSTAG_BACK_UP_RIGHT,
 33:   DMSTAG_DOWN_LEFT,
 34:   DMSTAG_DOWN,
 35:   DMSTAG_DOWN_RIGHT,
 36:   DMSTAG_LEFT,
 37:   DMSTAG_ELEMENT,
 38:   DMSTAG_RIGHT,
 39:   DMSTAG_UP_LEFT,
 40:   DMSTAG_UP,
 41:   DMSTAG_UP_RIGHT,
 42:   DMSTAG_FRONT_DOWN_LEFT,
 43:   DMSTAG_FRONT_DOWN,
 44:   DMSTAG_FRONT_DOWN_RIGHT,
 45:   DMSTAG_FRONT_LEFT,
 46:   DMSTAG_FRONT,
 47:   DMSTAG_FRONT_RIGHT,
 48:   DMSTAG_FRONT_UP_LEFT,
 49:   DMSTAG_FRONT_UP,
 50:   DMSTAG_FRONT_UP_RIGHT
 51: } DMStagStencilLocation;
 52: PETSC_EXTERN const char *const DMStagStencilLocations[]; /* Corresponding strings (see stagstencil.c) */

 54: /*S
 55:   DMStagStencil - data structure representing a degree of freedom on a `DMSTAG` grid

 57:   Data structure (C struct), analogous to describing a degree of freedom associated with a `DMSTAG` object,
 58:   in terms of a global element index in each of up to three directions, a "location" as defined by `DMStagStencilLocation`,
 59:   and a component number. Primarily for use with `DMStagMatSetValuesStencil()` (compare with use of `MatStencil` with `MatSetValuesStencil()`).

 61:   Note:
 62:   The component (c) field must always be set, even if there is a single component at a given location (in which case c should be set to 0).

 64: Level: beginner

 66: .seealso: [](chapter_stag), `DMSTAG`, `DMStagMatSetValuesStencil()`, `DMStagVecSetValuesStencil()`, `DMStagStencilLocation`, `DMStagSetStencilWidth()`,
 67:           `DMStagSetStencilType()`, `DMStagVecGetValuesStencil()`, `DMStagStencilLocation`
 68: S*/
 69: typedef struct {
 70:   DMStagStencilLocation loc;
 71:   PetscInt              i, j, k, c;
 72: } DMStagStencil;

 74: /*E
 75:   DMStagStencilType - Elementwise stencil type, determining which neighbors participate in communication

 77:   Level: beginner

 79: .seealso: [](chapter_stag), `DMSTAG`, `DMStagCreate1d()`, `DMStagCreate2d()`, `DMStagCreate3d()`, `DMStagStencil`, `DMDAStencilType`, `DMStagStencilLocation`
 80: E*/

 82: typedef enum {
 83:   DMSTAG_STENCIL_NONE = 0,
 84:   DMSTAG_STENCIL_STAR,
 85:   DMSTAG_STENCIL_BOX
 86: } DMStagStencilType;
 87: PETSC_EXTERN const char *const DMStagStencilTypes[]; /* Corresponding strings (see stagstencil.c) */

 89: PETSC_EXTERN PetscErrorCode DMCreate_Stag(DM);
 90: PETSC_EXTERN PetscErrorCode DMStagCreate1d(MPI_Comm, DMBoundaryType, PetscInt, PetscInt, PetscInt, DMStagStencilType, PetscInt, const PetscInt[], DM *);
 91: PETSC_EXTERN PetscErrorCode DMStagCreate2d(MPI_Comm, DMBoundaryType, DMBoundaryType, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, DMStagStencilType, PetscInt, const PetscInt[], const PetscInt[], DM *);
 92: PETSC_EXTERN PetscErrorCode DMStagCreate3d(MPI_Comm, DMBoundaryType, DMBoundaryType, DMBoundaryType, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, DMStagStencilType, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], DM *);
 93: PETSC_EXTERN PetscErrorCode DMStagCreateCompatibleDMStag(DM, PetscInt, PetscInt, PetscInt, PetscInt, DM *);
 94: PETSC_EXTERN PetscErrorCode DMStagCreateISFromStencils(DM, PetscInt, DMStagStencil *, IS *);
 95: PETSC_EXTERN PetscErrorCode DMStagGetBoundaryTypes(DM, DMBoundaryType *, DMBoundaryType *, DMBoundaryType *);
 96: PETSC_EXTERN PetscErrorCode DMStagGetCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
 97: PETSC_EXTERN PetscErrorCode DMStagGetDOF(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
 98: PETSC_EXTERN PetscErrorCode DMStagGetEntries(DM, PetscInt *);
 99: PETSC_EXTERN PetscErrorCode DMStagGetEntriesLocal(DM, PetscInt *);
100: PETSC_EXTERN PetscErrorCode DMStagGetEntriesPerElement(DM, PetscInt *);
101: PETSC_EXTERN PetscErrorCode DMStagGetGhostCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
102: PETSC_EXTERN PetscErrorCode DMStagGetGlobalSizes(DM, PetscInt *, PetscInt *, PetscInt *);
103: PETSC_EXTERN PetscErrorCode DMStagGetIsFirstRank(DM, PetscBool *, PetscBool *, PetscBool *);
104: PETSC_EXTERN PetscErrorCode DMStagGetIsLastRank(DM, PetscBool *, PetscBool *, PetscBool *);
105: PETSC_EXTERN PetscErrorCode DMStagGetLocalSizes(DM, PetscInt *, PetscInt *, PetscInt *);
106: PETSC_EXTERN PetscErrorCode DMStagGetLocationDOF(DM, DMStagStencilLocation, PetscInt *);
107: PETSC_EXTERN PetscErrorCode DMStagGetLocationSlot(DM, DMStagStencilLocation, PetscInt, PetscInt *);
108: PETSC_EXTERN PetscErrorCode DMStagGetNumRanks(DM, PetscInt *, PetscInt *, PetscInt *);
109: PETSC_EXTERN PetscErrorCode DMStagGetOwnershipRanges(DM, const PetscInt **, const PetscInt **, const PetscInt **);
110: PETSC_EXTERN PetscErrorCode DMStagGetProductCoordinateArrays(DM, void *, void *, void *);
111: PETSC_EXTERN PetscErrorCode DMStagGetProductCoordinateArraysRead(DM, void *, void *, void *);
112: PETSC_EXTERN PetscErrorCode DMStagGetProductCoordinateLocationSlot(DM, DMStagStencilLocation, PetscInt *);
113: PETSC_EXTERN PetscErrorCode DMStagGetStencilType(DM, DMStagStencilType *);
114: PETSC_EXTERN PetscErrorCode DMStagGetStencilWidth(DM, PetscInt *);
115: PETSC_EXTERN PetscErrorCode DMStagMatGetValuesStencil(DM, Mat, PetscInt, const DMStagStencil *, PetscInt, const DMStagStencil *, PetscScalar *);
116: PETSC_EXTERN PetscErrorCode DMStagMatSetValuesStencil(DM, Mat, PetscInt, const DMStagStencil *, PetscInt, const DMStagStencil *, const PetscScalar *, InsertMode);
117: PETSC_EXTERN PetscErrorCode DMStagMigrateVec(DM, Vec, DM, Vec);
118: PETSC_EXTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective(DM);
119: PETSC_EXTERN PetscErrorCode DMStagRestoreProductCoordinateArrays(DM, void *, void *, void *);
120: PETSC_EXTERN PetscErrorCode DMStagRestoreProductCoordinateArraysRead(DM, void *, void *, void *);
121: PETSC_EXTERN PetscErrorCode DMStagRestrictSimple(DM, Vec, DM, Vec);
122: PETSC_EXTERN PetscErrorCode DMStagSetBoundaryTypes(DM, DMBoundaryType, DMBoundaryType, DMBoundaryType);
123: PETSC_EXTERN PetscErrorCode DMStagSetCoordinateDMType(DM, DMType);
124: PETSC_EXTERN PetscErrorCode DMStagSetDOF(DM, PetscInt, PetscInt, PetscInt, PetscInt);
125: PETSC_EXTERN PetscErrorCode DMStagSetGlobalSizes(DM, PetscInt, PetscInt, PetscInt);
126: PETSC_EXTERN PetscErrorCode DMStagSetNumRanks(DM, PetscInt, PetscInt, PetscInt);
127: PETSC_EXTERN PetscErrorCode DMStagSetOwnershipRanges(DM, PetscInt const *, PetscInt const *, PetscInt const *);
128: PETSC_EXTERN PetscErrorCode DMStagSetStencilType(DM, DMStagStencilType);
129: PETSC_EXTERN PetscErrorCode DMStagSetStencilWidth(DM, PetscInt);
130: PETSC_EXTERN PetscErrorCode DMStagSetUniformCoordinates(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
131: PETSC_EXTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
132: PETSC_EXTERN PetscErrorCode DMStagSetUniformCoordinatesProduct(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
133: PETSC_EXTERN PetscErrorCode DMStagStencilToIndexLocal(DM, PetscInt, PetscInt, const DMStagStencil *, PetscInt *);
134: PETSC_EXTERN PetscErrorCode DMStagVecGetArray(DM, Vec, void *);
135: PETSC_EXTERN PetscErrorCode DMStagVecGetArrayRead(DM, Vec, void *);
136: PETSC_EXTERN PetscErrorCode DMStagVecGetValuesStencil(DM, Vec, PetscInt, const DMStagStencil *, PetscScalar *);
137: PETSC_EXTERN PetscErrorCode DMStagVecRestoreArray(DM, Vec, void *);
138: PETSC_EXTERN PetscErrorCode DMStagVecRestoreArrayRead(DM, Vec, void *);
139: PETSC_EXTERN PetscErrorCode DMStagVecSetValuesStencil(DM, Vec, PetscInt, const DMStagStencil *, const PetscScalar *, InsertMode);
140: PETSC_EXTERN PetscErrorCode DMStagVecSplitToDMDA(DM, Vec, DMStagStencilLocation, PetscInt, DM *, Vec *);

142: PETSC_DEPRECATED_FUNCTION("Use DMStagGetProductCoordinateArraysRead() (since version 3.13") static inline PetscErrorCode DMStagGet1dCoordinateArraysDOFRead(DM dm, void *ax, void *ay, void *az)
143: {
144:   return DMStagGetProductCoordinateArraysRead(dm, ax, ay, az);
145: }
146: PETSC_DEPRECATED_FUNCTION("Use DMStagGetProductCoordinateLocationSlot() (since version 3.13") static inline PetscErrorCode DMStagGet1dCoordinateLocationSlot(DM dm, DMStagStencilLocation loc, PetscInt *s)
147: {
148:   return DMStagGetProductCoordinateLocationSlot(dm, loc, s);
149: }
150: PETSC_DEPRECATED_FUNCTION("Use DMStagGetStencilType() (since version 3.11)") static inline PetscErrorCode DMStagGetGhostType(DM dm, DMStagStencilType *s)
151: {
152:   return DMStagGetStencilType(dm, s);
153: }
154: PETSC_DEPRECATED_FUNCTION("Use DMStagRestoreProductCoordinateArraysRead() (since version 3.13") static inline PetscErrorCode DMStagRestore1dCoordinateArraysDOFRead(DM dm, void *ax, void *ay, void *az)
155: {
156:   return DMStagRestoreProductCoordinateArraysRead(dm, ax, ay, az);
157: }
158: PETSC_DEPRECATED_FUNCTION("Use DMStagSetStencilType() (since version 3.11)") static inline PetscErrorCode DMStagSetGhostType(DM dm, DMStagStencilType *s)
159: {
160:   return DMStagGetStencilType(dm, s);
161: }
162: PETSC_DEPRECATED_FUNCTION("Use DMStagVecGetArray() (since version 3.13") static inline PetscErrorCode DMStagVecGetArrayDOF(DM dm, Vec v, void *a)
163: {
164:   return DMStagVecGetArray(dm, v, a);
165: }
166: PETSC_DEPRECATED_FUNCTION("Use DMStagVecGetArrayRead() (since version 3.13") static inline PetscErrorCode DMStagVecGetArrayDOFRead(DM dm, Vec v, void *a)
167: {
168:   return DMStagVecGetArrayRead(dm, v, a);
169: }
170: PETSC_DEPRECATED_FUNCTION("Use DMStagVecRestoreArray() (since version 3.13") static inline PetscErrorCode DMStagVecRestoreArrayDOF(DM dm, Vec v, void *a)
171: {
172:   return DMStagVecRestoreArray(dm, v, a);
173: }
174: PETSC_DEPRECATED_FUNCTION("Use DMStagVecRestoreArrayRead() (since version 3.13") static inline PetscErrorCode DMStagVecRestoreArrayDOFRead(DM dm, Vec v, void *a)
175: {
176:   return DMStagVecRestoreArrayRead(dm, v, a);
177: }

179: #endif // PETSC_DMSTAG_H