xref: /petsc/src/dm/impls/da/dapreallocate.c (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
1 #include <petsc/private/dmdaimpl.h>   /*I      "petscdmda.h"   I*/
2 #include <petsc/private/isimpl.h>
3 #include <petscsf.h>
4 
5 /*@
6   DMDASetPreallocationCenterDimension - Determine the topology used to determine adjacency
7 
8   Input Parameters:
9 + dm - The DM object
10 - preallocCenterDim - The dimension of points which connect adjacent entries
11 
12   Level: developer
13 
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
18 
19 .seealso: DMCreateMatrix(), DMDAPreallocateOperator()
20 @*/
21 PetscErrorCode DMDASetPreallocationCenterDimension(DM dm, PetscInt preallocCenterDim)
22 {
23   DM_DA *mesh = (DM_DA*) dm->data;
24 
25   PetscFunctionBegin;
26   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA);
27   mesh->preallocCenterDim = preallocCenterDim;
28   PetscFunctionReturn(0);
29 }
30 
31 /*@
32   DMDAGetPreallocationCenterDimension - Return the topology used to determine adjacency
33 
34   Input Parameter:
35 . dm - The DM object
36 
37   Output Parameter:
38 . preallocCenterDim - The dimension of points which connect adjacent entries
39 
40   Level: developer
41 
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
46 
47 .seealso: DMCreateMatrix(), DMDAPreallocateOperator(), DMDASetPreallocationCenterDimension()
48 @*/
49 PetscErrorCode DMDAGetPreallocationCenterDimension(DM dm, PetscInt *preallocCenterDim)
50 {
51   DM_DA *mesh = (DM_DA*) dm->data;
52 
53   PetscFunctionBegin;
54   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA);
55   PetscValidIntPointer(preallocCenterDim, 2);
56   *preallocCenterDim = mesh->preallocCenterDim;
57   PetscFunctionReturn(0);
58 }
59