Actual source code: dapreallocate.c
petsc-3.14.6 2021-03-30
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 DM object
10: - preallocCenterDim - The dimension of points which connect adjacent entries
12: Level: developer
14: Notes:
15: $ FEM: Two points p and q are adjacent if q \in closure(star(p)), preallocCenterDim = dim
16: $ FVM: Two points p and q are adjacent if q \in star(cone(p)), preallocCenterDim = dim-1
17: $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), preallocCenterDim = 0
19: .seealso: DMCreateMatrix(), DMDAPreallocateOperator()
20: @*/
21: PetscErrorCode DMDASetPreallocationCenterDimension(DM dm, PetscInt preallocCenterDim)
22: {
23: DM_DA *mesh = (DM_DA*) dm->data;
27: mesh->preallocCenterDim = preallocCenterDim;
28: return(0);
29: }
31: /*@
32: DMDAGetPreallocationCenterDimension - Return the topology used to determine adjacency
34: Input Parameter:
35: . dm - The DM object
37: Output Parameter:
38: . preallocCenterDim - The dimension of points which connect adjacent entries
40: Level: developer
42: Notes:
43: $ FEM: Two points p and q are adjacent if q \in closure(star(p)), preallocCenterDim = dim
44: $ FVM: Two points p and q are adjacent if q \in star(cone(p)), preallocCenterDim = dim-1
45: $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), preallocCenterDim = 0
47: .seealso: DMCreateMatrix(), DMDAPreallocateOperator(), DMDASetPreallocationCenterDimension()
48: @*/
49: PetscErrorCode DMDAGetPreallocationCenterDimension(DM dm, PetscInt *preallocCenterDim)
50: {
51: DM_DA *mesh = (DM_DA*) dm->data;
56: *preallocCenterDim = mesh->preallocCenterDim;
57: return(0);
58: }