Actual source code: dapreallocate.c

  1: #include <petsc/private/dmdaimpl.h>
  2: #include <petsc/private/isimpl.h>
  3: #include <petscsf.h>

  5: /*@
  6:   DMDASetPreallocationCenterDimension - Determine the topology used to determine adjacency

  8:   Input Parameters:
  9: + dm                - The `DMDA` object
 10: - preallocCenterDim - The dimension of points which connect adjacent entries

 12:   Level: developer

 14:   Notes:
 15: .vb
 16:      FEM:   Two points p and q are adjacent if q \in closure(star(p)), preallocCenterDim = dim
 17:      FVM:   Two points p and q are adjacent if q \in star(cone(p)),    preallocCenterDim = dim-1
 18:      FVM++: Two points p and q are adjacent if q \in star(closure(p)), preallocCenterDim = 0
 19: .ve

 21: .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateMatrix()`, `DMDAPreallocateOperator()`
 22: @*/
 23: PetscErrorCode DMDASetPreallocationCenterDimension(DM dm, PetscInt preallocCenterDim)
 24: {
 25:   DM_DA *mesh = (DM_DA *)dm->data;

 27:   PetscFunctionBegin;
 29:   mesh->preallocCenterDim = preallocCenterDim;
 30:   PetscFunctionReturn(PETSC_SUCCESS);
 31: }

 33: /*@
 34:   DMDAGetPreallocationCenterDimension - Return the topology used to determine adjacency

 36:   Input Parameter:
 37: . dm - The `DMDA` object

 39:   Output Parameter:
 40: . preallocCenterDim - The dimension of points which connect adjacent entries

 42:   Level: developer

 44:   Notes:
 45: .vb
 46:      FEM:   Two points p and q are adjacent if q \in closure(star(p)), preallocCenterDim = dim
 47:      FVM:   Two points p and q are adjacent if q \in star(cone(p)),    preallocCenterDim = dim-1
 48:      FVM++: Two points p and q are adjacent if q \in star(closure(p)), preallocCenterDim = 0
 49: .ve

 51: .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateMatrix()`, `DMDAPreallocateOperator()`, `DMDASetPreallocationCenterDimension()`
 52: @*/
 53: PetscErrorCode DMDAGetPreallocationCenterDimension(DM dm, PetscInt *preallocCenterDim)
 54: {
 55:   DM_DA *mesh = (DM_DA *)dm->data;

 57:   PetscFunctionBegin;
 59:   PetscAssertPointer(preallocCenterDim, 2);
 60:   *preallocCenterDim = mesh->preallocCenterDim;
 61:   PetscFunctionReturn(PETSC_SUCCESS);
 62: }