xref: /petsc/src/dm/impls/da/dapreallocate.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I      "petscdmda.h"   I*/
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>
33d183407SMatthew G. Knepley #include <petscsf.h>
43d183407SMatthew G. Knepley 
53d183407SMatthew G. Knepley /*@
63d183407SMatthew G. Knepley   DMDASetPreallocationCenterDimension - Determine the topology used to determine adjacency
73d183407SMatthew G. Knepley 
83d183407SMatthew G. Knepley   Input Parameters:
93d183407SMatthew G. Knepley + dm - The DM object
103d183407SMatthew G. Knepley - preallocCenterDim - The dimension of points which connect adjacent entries
113d183407SMatthew G. Knepley 
123d183407SMatthew G. Knepley   Level: developer
133d183407SMatthew G. Knepley 
143d183407SMatthew G. Knepley   Notes:
153d183407SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), preallocCenterDim = dim
163d183407SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    preallocCenterDim = dim-1
173d183407SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), preallocCenterDim = 0
183d183407SMatthew G. Knepley 
19db781477SPatrick Sanan .seealso: `DMCreateMatrix()`, `DMDAPreallocateOperator()`
203d183407SMatthew G. Knepley @*/
21*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetPreallocationCenterDimension(DM dm, PetscInt preallocCenterDim)
22*d71ae5a4SJacob Faibussowitsch {
233d183407SMatthew G. Knepley   DM_DA *mesh = (DM_DA *)dm->data;
243d183407SMatthew G. Knepley 
253d183407SMatthew G. Knepley   PetscFunctionBegin;
26a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
273d183407SMatthew G. Knepley   mesh->preallocCenterDim = preallocCenterDim;
283d183407SMatthew G. Knepley   PetscFunctionReturn(0);
293d183407SMatthew G. Knepley }
303d183407SMatthew G. Knepley 
313d183407SMatthew G. Knepley /*@
323d183407SMatthew G. Knepley   DMDAGetPreallocationCenterDimension - Return the topology used to determine adjacency
333d183407SMatthew G. Knepley 
343d183407SMatthew G. Knepley   Input Parameter:
353d183407SMatthew G. Knepley . dm - The DM object
363d183407SMatthew G. Knepley 
373d183407SMatthew G. Knepley   Output Parameter:
383d183407SMatthew G. Knepley . preallocCenterDim - The dimension of points which connect adjacent entries
393d183407SMatthew G. Knepley 
403d183407SMatthew G. Knepley   Level: developer
413d183407SMatthew G. Knepley 
423d183407SMatthew G. Knepley   Notes:
433d183407SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), preallocCenterDim = dim
443d183407SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    preallocCenterDim = dim-1
453d183407SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), preallocCenterDim = 0
463d183407SMatthew G. Knepley 
47db781477SPatrick Sanan .seealso: `DMCreateMatrix()`, `DMDAPreallocateOperator()`, `DMDASetPreallocationCenterDimension()`
483d183407SMatthew G. Knepley @*/
49*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetPreallocationCenterDimension(DM dm, PetscInt *preallocCenterDim)
50*d71ae5a4SJacob Faibussowitsch {
513d183407SMatthew G. Knepley   DM_DA *mesh = (DM_DA *)dm->data;
523d183407SMatthew G. Knepley 
533d183407SMatthew G. Knepley   PetscFunctionBegin;
54a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
553d183407SMatthew G. Knepley   PetscValidIntPointer(preallocCenterDim, 2);
563d183407SMatthew G. Knepley   *preallocCenterDim = mesh->preallocCenterDim;
573d183407SMatthew G. Knepley   PetscFunctionReturn(0);
583d183407SMatthew G. Knepley }
59