Actual source code: dmstagimpl.h
petsc-3.12.5 2020-03-29
1: #if !defined(DMSTAGIMPL_H_)
2: #define DMSTAGIMPL_H_
4: #include <petscdmstag.h>
5: #include <petsc/private/dmimpl.h>
7: #define DMSTAG_MAX_DIM 3
8: #define DMSTAG_MAX_STRATA DMSTAG_MAX_DIM + 1
10: /* This value is 1 + 3^DMSTAG_MAX DIM */
11: #define DMSTAG_NUMBER_LOCATIONS 28
13: typedef struct {
14: /* Fields which may require being set before DMSetUp() is called */
15: PetscInt N[DMSTAG_MAX_DIM]; /* Global dimensions (elements) */
16: PetscInt n[DMSTAG_MAX_DIM]; /* Local dimensions (elements) */
17: PetscInt dof[DMSTAG_MAX_STRATA]; /* Dof per point for each stratum */
18: DMStagStencilType stencilType; /* Elementwise stencil type */
19: PetscInt stencilWidth; /* Elementwise ghost width */
20: DMBoundaryType boundaryType[DMSTAG_MAX_DIM]; /* Physical domain ghosting type */
21: PetscInt nRanks[DMSTAG_MAX_DIM]; /* Ranks in each direction */
23: /* Additional fields populated by DMSetUp() */
24: PetscInt nGhost[DMSTAG_MAX_DIM]; /* Local dimensions (w/ ghosts) */
25: PetscInt start[DMSTAG_MAX_DIM]; /* First element number */
26: PetscInt startGhost[DMSTAG_MAX_DIM]; /* First element number (w/ ghosts) */
27: PetscMPIInt rank[DMSTAG_MAX_DIM]; /* Location in grid of ranks */
28: PetscMPIInt *neighbors; /* dim^3 local ranks */
29: PetscInt *l[DMSTAG_MAX_DIM]; /* Elements/rank in each direction */
30: VecScatter gtol; /* Global --> Local */
31: VecScatter ltog_injective; /* Local --> Global, injective */
32: PetscInt *locationOffsets; /* Offsets for points in loc. rep. */
34: /* Coordinates */
35: DMType coordinateDMType; /* DM type to create for coordinates */
37: /* Convenience (easily computed from the above) */
38: PetscInt entriesPerElement; /* Entries stored with each element */
39: PetscInt entries; /* Local number of entries */
40: PetscInt entriesGhost; /* Local numbers of entries w/ ghosts */
41: PetscBool firstRank[DMSTAG_MAX_DIM]; /* First rank in this dim? */
42: PetscBool lastRank[DMSTAG_MAX_DIM]; /* Last rank in this dim? */
43: } DM_Stag;
45: PETSC_INTERN PetscErrorCode DMSetUp_Stag_1d(DM);
46: PETSC_INTERN PetscErrorCode DMSetUp_Stag_2d(DM);
47: PETSC_INTERN PetscErrorCode DMSetUp_Stag_3d(DM);
48: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective_1d(DM);
49: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective_2d(DM);
50: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective_3d(DM);
51: PETSC_INTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit_1d(DM,PetscReal,PetscReal);
52: PETSC_INTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit_2d(DM,PetscReal,PetscReal,PetscReal,PetscReal);
53: PETSC_INTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit_3d(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
55: #endif