xref: /petsc/include/petscdm.h (revision bef01ed6f7ef6a42402475cc702347cd0deaa49f)
1e1589f56SBarry Smith /*
2e1589f56SBarry Smith       Objects to manage the interactions between the mesh data structures and the algebraic objects
3e1589f56SBarry Smith */
43c48a1e8SJed Brown #if !defined(__PETSCDM_H)
53c48a1e8SJed Brown #define __PETSCDM_H
62c8e378dSBarry Smith #include <petscmat.h>
71e25c274SJed Brown #include <petscdmtypes.h>
8e1589f56SBarry Smith 
9014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInitializePackage(const char[]);
10e1589f56SBarry Smith 
11014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID;
12e1589f56SBarry Smith 
1376bdecfbSBarry Smith /*J
14e1589f56SBarry Smith     DMType - String with the name of a PETSc DM or the creation function
15e1589f56SBarry Smith        with an optional dynamic library name, for example
166ae3a549SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mydmcreate()
17e1589f56SBarry Smith 
18e1589f56SBarry Smith    Level: beginner
19e1589f56SBarry Smith 
20e1589f56SBarry Smith .seealso: DMSetType(), DM
2176bdecfbSBarry Smith J*/
2219fd82e9SBarry Smith typedef const char* DMType;
23e1589f56SBarry Smith #define DMDA        "da"
24e1589f56SBarry Smith #define DMADDA      "adda"
25e1589f56SBarry Smith #define DMCOMPOSITE "composite"
26e1589f56SBarry Smith #define DMSLICED    "sliced"
27fe1899a2SJed Brown #define DMSHELL     "shell"
28b30b9b2eSMatthew G Knepley #define DMMESH      "mesh"
29ab7f58a0SBarry Smith #define DMPLEX      "plex"
30b30b9b2eSMatthew G Knepley #define DMCARTESIAN "cartesian"
318ac4e037SJed Brown #define DMREDUNDANT "redundant"
3201bc414fSDmitry Karpeev #define DMAKKT      "akkt"
333a19ef87SMatthew G Knepley #define DMPATCH     "patch"
34e1589f56SBarry Smith 
35140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList DMList;
36014dd563SJed Brown PETSC_EXTERN PetscBool         DMRegisterAllCalled;
37014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
3819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
3919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
40014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRegister(const char[],const char[],const char[],PetscErrorCode (*)(DM));
41014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRegisterAll(const char []);
42014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);
43e1589f56SBarry Smith 
44e1589f56SBarry Smith 
45e1589f56SBarry Smith /*MC
46e1589f56SBarry Smith   DMRegisterDynamic - Adds a new DM component implementation
47e1589f56SBarry Smith 
48e1589f56SBarry Smith   Synopsis:
49f2ba6396SBarry Smith   #include "petscdm.h"
50e1589f56SBarry Smith   PetscErrorCode DMRegisterDynamic(const char *name,const char *path,const char *func_name, PetscErrorCode (*create_func)(DM))
51e1589f56SBarry Smith 
52e1589f56SBarry Smith   Not Collective
53e1589f56SBarry Smith 
54e1589f56SBarry Smith   Input Parameters:
55e1589f56SBarry Smith + name        - The name of a new user-defined creation routine
56e1589f56SBarry Smith . path        - The path (either absolute or relative) of the library containing this routine
57e1589f56SBarry Smith . func_name   - The name of routine to create method context
58e1589f56SBarry Smith - create_func - The creation routine itself
59e1589f56SBarry Smith 
60e1589f56SBarry Smith   Notes:
61e1589f56SBarry Smith   DMRegisterDynamic() may be called multiple times to add several user-defined DMs
62e1589f56SBarry Smith 
63e1589f56SBarry Smith   If dynamic libraries are used, then the fourth input argument (routine_create) is ignored.
64e1589f56SBarry Smith 
65e1589f56SBarry Smith   Sample usage:
66e1589f56SBarry Smith .vb
67e1589f56SBarry Smith     DMRegisterDynamic("my_da","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyDMCreate", MyDMCreate);
68e1589f56SBarry Smith .ve
69e1589f56SBarry Smith 
70e1589f56SBarry Smith   Then, your DM type can be chosen with the procedural interface via
71e1589f56SBarry Smith .vb
72e1589f56SBarry Smith     DMCreate(MPI_Comm, DM *);
73e1589f56SBarry Smith     DMSetType(DM,"my_da_name");
74e1589f56SBarry Smith .ve
75e1589f56SBarry Smith    or at runtime via the option
76e1589f56SBarry Smith .vb
77e1589f56SBarry Smith     -da_type my_da_name
78e1589f56SBarry Smith .ve
79e1589f56SBarry Smith 
80e1589f56SBarry Smith   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
81e1589f56SBarry Smith          If your function is not being put into a shared library then use DMRegister() instead
82e1589f56SBarry Smith 
83e1589f56SBarry Smith   Level: advanced
84e1589f56SBarry Smith 
85e1589f56SBarry Smith .keywords: DM, register
86e1589f56SBarry Smith .seealso: DMRegisterAll(), DMRegisterDestroy(), DMRegister()
87e1589f56SBarry Smith M*/
88e1589f56SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
89e1589f56SBarry Smith #define DMRegisterDynamic(a,b,c,d) DMRegister(a,b,c,0)
90e1589f56SBarry Smith #else
91e1589f56SBarry Smith #define DMRegisterDynamic(a,b,c,d) DMRegister(a,b,c,d)
92e1589f56SBarry Smith #endif
93e1589f56SBarry Smith 
94014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
1062348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
1072348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMappingBlock(DM,ISLocalToGlobalMapping*);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
11219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,MatType,ISColoring*);
11319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,MatType,Mat*);
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,VecScatter*);
117aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*);
118aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
128ca266f36SBarry Smith PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,const char[]);
129ca266f36SBarry Smith 
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
133baf369e7SPeter Brune PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
13819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);
139e1589f56SBarry Smith 
1406636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
1416636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
1426636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
1436636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
1446636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
145e87bb0d3SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,IS*);
1466636e97aSMatthew G Knepley 
1475dbd56e3SPeter Brune /* block hook interface */
148be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
149be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);
1505dbd56e3SPeter Brune 
151014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
15219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
15319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
156014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
157014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
158014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
159b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
160014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);
16193d92d96SBarry Smith 
16281d26defSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *);
16316621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
1648d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
165e30e807fSPeter Brune PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
166e7c4fc90SDmitry Karpeev 
167014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
168014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
169014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
170e1589f56SBarry Smith 
1715f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
1725f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
173c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
174c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
1755f1ad066SMatthew G Knepley 
176e1589f56SBarry Smith typedef struct NLF_DAAD* NLF;
177e1589f56SBarry Smith 
178bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221
1797da65231SMatthew G Knepley 
1807da65231SMatthew G Knepley /* FEM support */
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
1837da65231SMatthew G Knepley 
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
187eaf8d80aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
191b21d0597SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
192057b4bcdSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
19388ed4aceSMatthew G Knepley 
194af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
195af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
196af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
197af122d2aSMatthew G Knepley 
1987da65231SMatthew G Knepley typedef struct {
1997da65231SMatthew G Knepley   PetscInt         numQuadPoints; /* The number of quadrature points on an element */
2007da65231SMatthew G Knepley   const PetscReal *quadPoints;    /* The quadrature point coordinates */
2017da65231SMatthew G Knepley   const PetscReal *quadWeights;   /* The quadrature weights */
2027da65231SMatthew G Knepley   PetscInt         numBasisFuncs; /* The number of finite element basis functions on an element */
2037da65231SMatthew G Knepley   PetscInt         numComponents; /* The number of components for each basis function */
2047da65231SMatthew G Knepley   const PetscReal *basis;         /* The basis functions tabulated at the quadrature points */
2057da65231SMatthew G Knepley   const PetscReal *basisDer;      /* The basis function derivatives tabulated at the quadrature points */
2067da65231SMatthew G Knepley } PetscQuadrature;
20783e5549fSMatthew G Knepley 
20892fc3d87SMatthew G Knepley typedef struct {
20992fc3d87SMatthew G Knepley   PetscQuadrature *quad;
21092fc3d87SMatthew G Knepley   void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
21192fc3d87SMatthew G Knepley   void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
21292fc3d87SMatthew G Knepley   void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */
21392fc3d87SMatthew G Knepley   void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */
21492fc3d87SMatthew G Knepley   void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */
21592fc3d87SMatthew G Knepley   void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */
216*bef01ed6SMatthew G. Knepley   void (**bcFuncs)(const PetscReal[], PetscScalar *); /* The boundary condition function for each field component */
217652b88e8SMatthew G. Knepley   PetscQuadrature *quadBd;
218652b88e8SMatthew G. Knepley   void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
219652b88e8SMatthew G. Knepley   void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
22092fc3d87SMatthew G Knepley } PetscFEM;
22192fc3d87SMatthew G Knepley 
22283e5549fSMatthew G Knepley typedef enum {PETSC_UNIT_LENGTH, PETSC_UNIT_MASS, PETSC_UNIT_TIME, PETSC_UNIT_CURRENT, PETSC_UNIT_TEMPERATURE, PETSC_UNIT_AMOUNT, PETSC_UNIT_LUMINOSITY, NUM_PETSC_UNITS} PetscUnit;
223e87bb0d3SMatthew G Knepley 
224e87bb0d3SMatthew G Knepley struct _DMInterpolationInfo {
225e87bb0d3SMatthew G Knepley   MPI_Comm   comm;
226e87bb0d3SMatthew G Knepley   PetscInt   dim;    /*1 The spatial dimension of points */
227e87bb0d3SMatthew G Knepley   PetscInt   nInput; /* The number of input points */
228e87bb0d3SMatthew G Knepley   PetscReal *points; /* The input point coordinates */
229e87bb0d3SMatthew G Knepley   PetscInt  *cells;  /* The cell containing each point */
230e87bb0d3SMatthew G Knepley   PetscInt   n;      /* The number of local points */
231e87bb0d3SMatthew G Knepley   Vec        coords; /* The point coordinates */
232e87bb0d3SMatthew G Knepley   PetscInt   dof;    /* The number of components to interpolate */
233e87bb0d3SMatthew G Knepley };
234e87bb0d3SMatthew G Knepley typedef struct _DMInterpolationInfo *DMInterpolationInfo;
235e87bb0d3SMatthew G Knepley 
236e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
237e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
238e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
239e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
240e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
241e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
242e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
243e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
244e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
245e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
246e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
247e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
248e1589f56SBarry Smith #endif
249