Actual source code: dmstagimpl.h
petsc-3.13.6 2020-09-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, set by DMStagInitialize().
15: Some may be adjusted by DMSetUp() */
16: PetscInt N[DMSTAG_MAX_DIM]; /* Global dimensions (elements) */
17: PetscInt n[DMSTAG_MAX_DIM]; /* Local dimensions (elements) */
18: PetscInt *l[DMSTAG_MAX_DIM]; /* Elements/rank in each direction */
19: PetscInt dof[DMSTAG_MAX_STRATA]; /* Dof per point for each stratum */
20: DMStagStencilType stencilType; /* Elementwise stencil type */
21: PetscInt stencilWidth; /* Elementwise ghost width */
22: DMBoundaryType boundaryType[DMSTAG_MAX_DIM]; /* Physical domain ghosting type */
23: PetscInt nRanks[DMSTAG_MAX_DIM]; /* Ranks in each direction */
25: /* Fields unrelated to setup */
26: DMType coordinateDMType; /* DM type to create for coordinates */
28: /* Data above is copied by DMStagDuplicateWithoutSetup(), while data below is not */
30: /* Fields populated by DMSetUp() */
31: PetscInt nGhost[DMSTAG_MAX_DIM]; /* Local dimensions (w/ ghosts) */
32: PetscInt start[DMSTAG_MAX_DIM]; /* First element number */
33: PetscInt startGhost[DMSTAG_MAX_DIM]; /* First element number (w/ ghosts) */
34: PetscMPIInt rank[DMSTAG_MAX_DIM]; /* Location in grid of ranks */
35: PetscMPIInt *neighbors; /* dim^3 local ranks */
36: VecScatter gtol; /* Global --> Local */
37: VecScatter ltog_injective; /* Local --> Global, injective */
38: PetscInt *locationOffsets; /* Offsets for points in loc. rep. */
40: /* Additional convenience fields populated by DMSetUp() (easily computed from the above) */
41: PetscInt entriesPerElement; /* Entries stored with each element */
42: PetscInt entries; /* Local number of entries */
43: PetscInt entriesGhost; /* Local numbers of entries w/ ghosts */
44: PetscBool firstRank[DMSTAG_MAX_DIM]; /* First rank in this dim? */
45: PetscBool lastRank[DMSTAG_MAX_DIM]; /* Last rank in this dim? */
47: } DM_Stag;
49: PETSC_INTERN PetscErrorCode DMStagDuplicateWithoutSetup(DM,MPI_Comm,DM*);
50: PETSC_INTERN PetscErrorCode DMStagInitialize(DMBoundaryType,DMBoundaryType,DMBoundaryType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,DMStagStencilType,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM);
51: PETSC_INTERN PetscErrorCode DMSetUp_Stag_1d(DM);
52: PETSC_INTERN PetscErrorCode DMSetUp_Stag_2d(DM);
53: PETSC_INTERN PetscErrorCode DMSetUp_Stag_3d(DM);
54: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective_1d(DM);
55: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective_2d(DM);
56: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective_3d(DM);
57: PETSC_INTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit_1d(DM,PetscReal,PetscReal);
58: PETSC_INTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit_2d(DM,PetscReal,PetscReal,PetscReal,PetscReal);
59: PETSC_INTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit_3d(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
61: #endif