xref: /petsc/src/dm/interface/dm.c (revision 9435951e5d58a9a8b01c7ece939244d29b452105)
1b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/
20c312b8eSJed Brown #include <petscsf.h>
32764a2aaSMatthew G. Knepley #include <petscds.h>
447c6ae99SBarry Smith 
5732e2eb9SMatthew G Knepley PetscClassId  DM_CLASSID;
6f089877aSRichard Tran Mills PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocalToLocal;
767a56275SMatthew G Knepley 
8bff4a2f0SMatthew G. Knepley const char *const DMBoundaryTypes[] = {"NONE","GHOSTED","MIRROR","PERIODIC","TWIST","DM_BOUNDARY_",0};
9bff4a2f0SMatthew G. Knepley 
1047c6ae99SBarry Smith #undef __FUNCT__
11a4121054SBarry Smith #define __FUNCT__ "DMCreate"
12a4121054SBarry Smith /*@
13de043629SMatthew G Knepley   DMCreate - Creates an empty DM object. The type can then be set with DMSetType().
14a4121054SBarry Smith 
15a4121054SBarry Smith    If you never  call DMSetType()  it will generate an
16a4121054SBarry Smith    error when you try to use the vector.
17a4121054SBarry Smith 
18a4121054SBarry Smith   Collective on MPI_Comm
19a4121054SBarry Smith 
20a4121054SBarry Smith   Input Parameter:
21a4121054SBarry Smith . comm - The communicator for the DM object
22a4121054SBarry Smith 
23a4121054SBarry Smith   Output Parameter:
24a4121054SBarry Smith . dm - The DM object
25a4121054SBarry Smith 
26a4121054SBarry Smith   Level: beginner
27a4121054SBarry Smith 
28a4121054SBarry Smith .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE
29a4121054SBarry Smith @*/
307087cfbeSBarry Smith PetscErrorCode  DMCreate(MPI_Comm comm,DM *dm)
31a4121054SBarry Smith {
32a4121054SBarry Smith   DM             v;
33a4121054SBarry Smith   PetscErrorCode ierr;
34a4121054SBarry Smith 
35a4121054SBarry Smith   PetscFunctionBegin;
361411c6eeSJed Brown   PetscValidPointer(dm,2);
370298fd71SBarry Smith   *dm = NULL;
388b984314SMatthew G. Knepley   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
39607a6623SBarry Smith   ierr = VecInitializePackage();CHKERRQ(ierr);
40607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
41607a6623SBarry Smith   ierr = DMInitializePackage();CHKERRQ(ierr);
42a4121054SBarry Smith 
4367c2884eSBarry Smith   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
44a4121054SBarry Smith   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
451411c6eeSJed Brown 
46e7c4fc90SDmitry Karpeev 
470298fd71SBarry Smith   v->ltogmap              = NULL;
481411c6eeSJed Brown   v->bs                   = 1;
49171400e9SBarry Smith   v->coloringtype         = IS_COLORING_GLOBAL;
5088ed4aceSMatthew G Knepley   ierr                    = PetscSFCreate(comm, &v->sf);CHKERRQ(ierr);
5188ed4aceSMatthew G Knepley   ierr                    = PetscSFCreate(comm, &v->defaultSF);CHKERRQ(ierr);
520298fd71SBarry Smith   v->defaultSection       = NULL;
530298fd71SBarry Smith   v->defaultGlobalSection = NULL;
54c6b900c6SMatthew G. Knepley   v->L                    = NULL;
55c6b900c6SMatthew G. Knepley   v->maxCell              = NULL;
569a9a41abSToby Isaac   v->dimEmbed             = PETSC_DEFAULT;
57435a35e8SMatthew G Knepley   {
58435a35e8SMatthew G Knepley     PetscInt i;
59435a35e8SMatthew G Knepley     for (i = 0; i < 10; ++i) {
600298fd71SBarry Smith       v->nullspaceConstructors[i] = NULL;
61435a35e8SMatthew G Knepley     }
62435a35e8SMatthew G Knepley   }
632764a2aaSMatthew G. Knepley   ierr = PetscDSCreate(comm, &v->prob);CHKERRQ(ierr);
6414f150ffSMatthew G. Knepley   v->dmBC = NULL;
65f4d763aaSMatthew G. Knepley   v->outputSequenceNum = -1;
66cdb7a50dSMatthew G. Knepley   v->outputSequenceVal = 0.0;
67c0dedaeaSBarry Smith   ierr = DMSetVecType(v,VECSTANDARD);CHKERRQ(ierr);
68b412c318SBarry Smith   ierr = DMSetMatType(v,MATAIJ);CHKERRQ(ierr);
691411c6eeSJed Brown   *dm = v;
70a4121054SBarry Smith   PetscFunctionReturn(0);
71a4121054SBarry Smith }
72a4121054SBarry Smith 
73a4121054SBarry Smith #undef __FUNCT__
7438221697SMatthew G. Knepley #define __FUNCT__ "DMClone"
7538221697SMatthew G. Knepley /*@
7638221697SMatthew G. Knepley   DMClone - Creates a DM object with the same topology as the original.
7738221697SMatthew G. Knepley 
7838221697SMatthew G. Knepley   Collective on MPI_Comm
7938221697SMatthew G. Knepley 
8038221697SMatthew G. Knepley   Input Parameter:
8138221697SMatthew G. Knepley . dm - The original DM object
8238221697SMatthew G. Knepley 
8338221697SMatthew G. Knepley   Output Parameter:
8438221697SMatthew G. Knepley . newdm  - The new DM object
8538221697SMatthew G. Knepley 
8638221697SMatthew G. Knepley   Level: beginner
8738221697SMatthew G. Knepley 
8838221697SMatthew G. Knepley .keywords: DM, topology, create
8938221697SMatthew G. Knepley @*/
9038221697SMatthew G. Knepley PetscErrorCode DMClone(DM dm, DM *newdm)
9138221697SMatthew G. Knepley {
9238221697SMatthew G. Knepley   PetscSF        sf;
9338221697SMatthew G. Knepley   Vec            coords;
9438221697SMatthew G. Knepley   void          *ctx;
951de53e9aSMatthew G. Knepley   PetscInt       dim;
9638221697SMatthew G. Knepley   PetscErrorCode ierr;
9738221697SMatthew G. Knepley 
9838221697SMatthew G. Knepley   PetscFunctionBegin;
9938221697SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10038221697SMatthew G. Knepley   PetscValidPointer(newdm,2);
10138221697SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), newdm);CHKERRQ(ierr);
1021de53e9aSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1031de53e9aSMatthew G. Knepley   ierr = DMSetDimension(*newdm, dim);CHKERRQ(ierr);
10438221697SMatthew G. Knepley   if (dm->ops->clone) {
10538221697SMatthew G. Knepley     ierr = (*dm->ops->clone)(dm, newdm);CHKERRQ(ierr);
10638221697SMatthew G. Knepley   }
107cd27c7c9SMatthew G. Knepley   (*newdm)->setupcalled = PETSC_TRUE;
10838221697SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
10938221697SMatthew G. Knepley   ierr = DMSetPointSF(*newdm, sf);CHKERRQ(ierr);
11038221697SMatthew G. Knepley   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
11138221697SMatthew G. Knepley   ierr = DMSetApplicationContext(*newdm, ctx);CHKERRQ(ierr);
11238221697SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr);
11338221697SMatthew G. Knepley   if (coords) {
11438221697SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*newdm, coords);CHKERRQ(ierr);
11538221697SMatthew G. Knepley   } else {
11638221697SMatthew G. Knepley     ierr = DMGetCoordinates(dm, &coords);CHKERRQ(ierr);
11738221697SMatthew G. Knepley     if (coords) {ierr = DMSetCoordinates(*newdm, coords);CHKERRQ(ierr);}
11838221697SMatthew G. Knepley   }
119c6b900c6SMatthew G. Knepley   if (dm->maxCell) {
120c6b900c6SMatthew G. Knepley     const PetscReal *maxCell, *L;
121c6b900c6SMatthew G. Knepley     ierr = DMGetPeriodicity(dm,     &maxCell, &L);CHKERRQ(ierr);
122c6b900c6SMatthew G. Knepley     ierr = DMSetPeriodicity(*newdm,  maxCell,  L);CHKERRQ(ierr);
123c6b900c6SMatthew G. Knepley   }
12438221697SMatthew G. Knepley   PetscFunctionReturn(0);
12538221697SMatthew G. Knepley }
12638221697SMatthew G. Knepley 
12738221697SMatthew G. Knepley #undef __FUNCT__
1289a42bb27SBarry Smith #define __FUNCT__ "DMSetVecType"
1299a42bb27SBarry Smith /*@C
130564755cdSBarry Smith        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
1319a42bb27SBarry Smith 
132aa219208SBarry Smith    Logically Collective on DMDA
1339a42bb27SBarry Smith 
1349a42bb27SBarry Smith    Input Parameter:
1359a42bb27SBarry Smith +  da - initial distributed array
1368154be41SBarry Smith .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
1379a42bb27SBarry Smith 
1389a42bb27SBarry Smith    Options Database:
139dd85299cSBarry Smith .   -dm_vec_type ctype
1409a42bb27SBarry Smith 
1419a42bb27SBarry Smith    Level: intermediate
1429a42bb27SBarry Smith 
143c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType, DMGetVecType()
1449a42bb27SBarry Smith @*/
14519fd82e9SBarry Smith PetscErrorCode  DMSetVecType(DM da,VecType ctype)
1469a42bb27SBarry Smith {
1479a42bb27SBarry Smith   PetscErrorCode ierr;
1489a42bb27SBarry Smith 
1499a42bb27SBarry Smith   PetscFunctionBegin;
1509a42bb27SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1519a42bb27SBarry Smith   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
15219fd82e9SBarry Smith   ierr = PetscStrallocpy(ctype,(char**)&da->vectype);CHKERRQ(ierr);
1539a42bb27SBarry Smith   PetscFunctionReturn(0);
1549a42bb27SBarry Smith }
1559a42bb27SBarry Smith 
1569a42bb27SBarry Smith #undef __FUNCT__
157c0dedaeaSBarry Smith #define __FUNCT__ "DMGetVecType"
158c0dedaeaSBarry Smith /*@C
159c0dedaeaSBarry Smith        DMGetVecType - Gets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
160c0dedaeaSBarry Smith 
161c0dedaeaSBarry Smith    Logically Collective on DMDA
162c0dedaeaSBarry Smith 
163c0dedaeaSBarry Smith    Input Parameter:
164c0dedaeaSBarry Smith .  da - initial distributed array
165c0dedaeaSBarry Smith 
166c0dedaeaSBarry Smith    Output Parameter:
167c0dedaeaSBarry Smith .  ctype - the vector type
168c0dedaeaSBarry Smith 
169c0dedaeaSBarry Smith    Level: intermediate
170c0dedaeaSBarry Smith 
171c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
172c0dedaeaSBarry Smith @*/
173c0dedaeaSBarry Smith PetscErrorCode  DMGetVecType(DM da,VecType *ctype)
174c0dedaeaSBarry Smith {
175c0dedaeaSBarry Smith   PetscFunctionBegin;
176c0dedaeaSBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
177c0dedaeaSBarry Smith   *ctype = da->vectype;
178c0dedaeaSBarry Smith   PetscFunctionReturn(0);
179c0dedaeaSBarry Smith }
180c0dedaeaSBarry Smith 
181c0dedaeaSBarry Smith #undef __FUNCT__
1825f1ad066SMatthew G Knepley #define __FUNCT__ "VecGetDM"
1835f1ad066SMatthew G Knepley /*@
18434f98d34SBarry Smith   VecGetDM - Gets the DM defining the data layout of the vector
1855f1ad066SMatthew G Knepley 
1865f1ad066SMatthew G Knepley   Not collective
1875f1ad066SMatthew G Knepley 
1885f1ad066SMatthew G Knepley   Input Parameter:
1895f1ad066SMatthew G Knepley . v - The Vec
1905f1ad066SMatthew G Knepley 
1915f1ad066SMatthew G Knepley   Output Parameter:
1925f1ad066SMatthew G Knepley . dm - The DM
1935f1ad066SMatthew G Knepley 
1945f1ad066SMatthew G Knepley   Level: intermediate
1955f1ad066SMatthew G Knepley 
1965f1ad066SMatthew G Knepley .seealso: VecSetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType()
1975f1ad066SMatthew G Knepley @*/
1985f1ad066SMatthew G Knepley PetscErrorCode VecGetDM(Vec v, DM *dm)
1995f1ad066SMatthew G Knepley {
2005f1ad066SMatthew G Knepley   PetscErrorCode ierr;
2015f1ad066SMatthew G Knepley 
2025f1ad066SMatthew G Knepley   PetscFunctionBegin;
2035f1ad066SMatthew G Knepley   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
2045f1ad066SMatthew G Knepley   PetscValidPointer(dm,2);
2055f1ad066SMatthew G Knepley   ierr = PetscObjectQuery((PetscObject) v, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr);
2065f1ad066SMatthew G Knepley   PetscFunctionReturn(0);
2075f1ad066SMatthew G Knepley }
2085f1ad066SMatthew G Knepley 
2095f1ad066SMatthew G Knepley #undef __FUNCT__
2105f1ad066SMatthew G Knepley #define __FUNCT__ "VecSetDM"
2115f1ad066SMatthew G Knepley /*@
2125f1ad066SMatthew G Knepley   VecSetDM - Sets the DM defining the data layout of the vector
2135f1ad066SMatthew G Knepley 
2145f1ad066SMatthew G Knepley   Not collective
2155f1ad066SMatthew G Knepley 
2165f1ad066SMatthew G Knepley   Input Parameters:
2175f1ad066SMatthew G Knepley + v - The Vec
2185f1ad066SMatthew G Knepley - dm - The DM
2195f1ad066SMatthew G Knepley 
2205f1ad066SMatthew G Knepley   Level: intermediate
2215f1ad066SMatthew G Knepley 
2225f1ad066SMatthew G Knepley .seealso: VecGetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType()
2235f1ad066SMatthew G Knepley @*/
2245f1ad066SMatthew G Knepley PetscErrorCode VecSetDM(Vec v, DM dm)
2255f1ad066SMatthew G Knepley {
2265f1ad066SMatthew G Knepley   PetscErrorCode ierr;
2275f1ad066SMatthew G Knepley 
2285f1ad066SMatthew G Knepley   PetscFunctionBegin;
2295f1ad066SMatthew G Knepley   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
230d7f50e27SLisandro Dalcin   if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2);
2315f1ad066SMatthew G Knepley   ierr = PetscObjectCompose((PetscObject) v, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr);
2325f1ad066SMatthew G Knepley   PetscFunctionReturn(0);
2335f1ad066SMatthew G Knepley }
2345f1ad066SMatthew G Knepley 
2355f1ad066SMatthew G Knepley #undef __FUNCT__
236521d9a4cSLisandro Dalcin #define __FUNCT__ "DMSetMatType"
237521d9a4cSLisandro Dalcin /*@C
238521d9a4cSLisandro Dalcin        DMSetMatType - Sets the type of matrix created with DMCreateMatrix()
239521d9a4cSLisandro Dalcin 
240521d9a4cSLisandro Dalcin    Logically Collective on DM
241521d9a4cSLisandro Dalcin 
242521d9a4cSLisandro Dalcin    Input Parameter:
243521d9a4cSLisandro Dalcin +  dm - the DM context
244521d9a4cSLisandro Dalcin .  ctype - the matrix type
245521d9a4cSLisandro Dalcin 
246521d9a4cSLisandro Dalcin    Options Database:
247521d9a4cSLisandro Dalcin .   -dm_mat_type ctype
248521d9a4cSLisandro Dalcin 
249521d9a4cSLisandro Dalcin    Level: intermediate
250521d9a4cSLisandro Dalcin 
251c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMGetMatType()
252521d9a4cSLisandro Dalcin @*/
25319fd82e9SBarry Smith PetscErrorCode  DMSetMatType(DM dm,MatType ctype)
254521d9a4cSLisandro Dalcin {
255521d9a4cSLisandro Dalcin   PetscErrorCode ierr;
25688f0584fSBarry Smith 
257521d9a4cSLisandro Dalcin   PetscFunctionBegin;
258521d9a4cSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
259521d9a4cSLisandro Dalcin   ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
26019fd82e9SBarry Smith   ierr = PetscStrallocpy(ctype,(char**)&dm->mattype);CHKERRQ(ierr);
261521d9a4cSLisandro Dalcin   PetscFunctionReturn(0);
262521d9a4cSLisandro Dalcin }
263521d9a4cSLisandro Dalcin 
264521d9a4cSLisandro Dalcin #undef __FUNCT__
265c0dedaeaSBarry Smith #define __FUNCT__ "DMGetMatType"
266c0dedaeaSBarry Smith /*@C
267c0dedaeaSBarry Smith        DMGetMatType - Gets the type of matrix created with DMCreateMatrix()
268c0dedaeaSBarry Smith 
269c0dedaeaSBarry Smith    Logically Collective on DM
270c0dedaeaSBarry Smith 
271c0dedaeaSBarry Smith    Input Parameter:
272c0dedaeaSBarry Smith .  dm - the DM context
273c0dedaeaSBarry Smith 
274c0dedaeaSBarry Smith    Output Parameter:
275c0dedaeaSBarry Smith .  ctype - the matrix type
276c0dedaeaSBarry Smith 
277c0dedaeaSBarry Smith    Options Database:
278c0dedaeaSBarry Smith .   -dm_mat_type ctype
279c0dedaeaSBarry Smith 
280c0dedaeaSBarry Smith    Level: intermediate
281c0dedaeaSBarry Smith 
282c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMSetMatType()
283c0dedaeaSBarry Smith @*/
284c0dedaeaSBarry Smith PetscErrorCode  DMGetMatType(DM dm,MatType *ctype)
285c0dedaeaSBarry Smith {
286c0dedaeaSBarry Smith   PetscFunctionBegin;
287c0dedaeaSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
288c0dedaeaSBarry Smith   *ctype = dm->mattype;
289c0dedaeaSBarry Smith   PetscFunctionReturn(0);
290c0dedaeaSBarry Smith }
291c0dedaeaSBarry Smith 
292c0dedaeaSBarry Smith #undef __FUNCT__
293c688c046SMatthew G Knepley #define __FUNCT__ "MatGetDM"
294c688c046SMatthew G Knepley /*@
29534f98d34SBarry Smith   MatGetDM - Gets the DM defining the data layout of the matrix
296c688c046SMatthew G Knepley 
297c688c046SMatthew G Knepley   Not collective
298c688c046SMatthew G Knepley 
299c688c046SMatthew G Knepley   Input Parameter:
300c688c046SMatthew G Knepley . A - The Mat
301c688c046SMatthew G Knepley 
302c688c046SMatthew G Knepley   Output Parameter:
303c688c046SMatthew G Knepley . dm - The DM
304c688c046SMatthew G Knepley 
305c688c046SMatthew G Knepley   Level: intermediate
306c688c046SMatthew G Knepley 
307c688c046SMatthew G Knepley .seealso: MatSetDM(), DMCreateMatrix(), DMSetMatType()
308c688c046SMatthew G Knepley @*/
309c688c046SMatthew G Knepley PetscErrorCode MatGetDM(Mat A, DM *dm)
310c688c046SMatthew G Knepley {
311c688c046SMatthew G Knepley   PetscErrorCode ierr;
312c688c046SMatthew G Knepley 
313c688c046SMatthew G Knepley   PetscFunctionBegin;
314c688c046SMatthew G Knepley   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
315c688c046SMatthew G Knepley   PetscValidPointer(dm,2);
316c688c046SMatthew G Knepley   ierr = PetscObjectQuery((PetscObject) A, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr);
317c688c046SMatthew G Knepley   PetscFunctionReturn(0);
318c688c046SMatthew G Knepley }
319c688c046SMatthew G Knepley 
320c688c046SMatthew G Knepley #undef __FUNCT__
321c688c046SMatthew G Knepley #define __FUNCT__ "MatSetDM"
322c688c046SMatthew G Knepley /*@
323c688c046SMatthew G Knepley   MatSetDM - Sets the DM defining the data layout of the matrix
324c688c046SMatthew G Knepley 
325c688c046SMatthew G Knepley   Not collective
326c688c046SMatthew G Knepley 
327c688c046SMatthew G Knepley   Input Parameters:
328c688c046SMatthew G Knepley + A - The Mat
329c688c046SMatthew G Knepley - dm - The DM
330c688c046SMatthew G Knepley 
331c688c046SMatthew G Knepley   Level: intermediate
332c688c046SMatthew G Knepley 
333c688c046SMatthew G Knepley .seealso: MatGetDM(), DMCreateMatrix(), DMSetMatType()
334c688c046SMatthew G Knepley @*/
335c688c046SMatthew G Knepley PetscErrorCode MatSetDM(Mat A, DM dm)
336c688c046SMatthew G Knepley {
337c688c046SMatthew G Knepley   PetscErrorCode ierr;
338c688c046SMatthew G Knepley 
339c688c046SMatthew G Knepley   PetscFunctionBegin;
340c688c046SMatthew G Knepley   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
3418865f1eaSKarl Rupp   if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2);
342c688c046SMatthew G Knepley   ierr = PetscObjectCompose((PetscObject) A, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr);
343c688c046SMatthew G Knepley   PetscFunctionReturn(0);
344c688c046SMatthew G Knepley }
345c688c046SMatthew G Knepley 
346c688c046SMatthew G Knepley #undef __FUNCT__
3479a42bb27SBarry Smith #define __FUNCT__ "DMSetOptionsPrefix"
3489a42bb27SBarry Smith /*@C
3499a42bb27SBarry Smith    DMSetOptionsPrefix - Sets the prefix used for searching for all
350aa219208SBarry Smith    DMDA options in the database.
3519a42bb27SBarry Smith 
352aa219208SBarry Smith    Logically Collective on DMDA
3539a42bb27SBarry Smith 
3549a42bb27SBarry Smith    Input Parameter:
355aa219208SBarry Smith +  da - the DMDA context
3569a42bb27SBarry Smith -  prefix - the prefix to prepend to all option names
3579a42bb27SBarry Smith 
3589a42bb27SBarry Smith    Notes:
3599a42bb27SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
3609a42bb27SBarry Smith    The first character of all runtime options is AUTOMATICALLY the hyphen.
3619a42bb27SBarry Smith 
3629a42bb27SBarry Smith    Level: advanced
3639a42bb27SBarry Smith 
364aa219208SBarry Smith .keywords: DMDA, set, options, prefix, database
3659a42bb27SBarry Smith 
3669a42bb27SBarry Smith .seealso: DMSetFromOptions()
3679a42bb27SBarry Smith @*/
3687087cfbeSBarry Smith PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
3699a42bb27SBarry Smith {
3709a42bb27SBarry Smith   PetscErrorCode ierr;
3719a42bb27SBarry Smith 
3729a42bb27SBarry Smith   PetscFunctionBegin;
3739a42bb27SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3749a42bb27SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
3759a42bb27SBarry Smith   PetscFunctionReturn(0);
3769a42bb27SBarry Smith }
3779a42bb27SBarry Smith 
3789a42bb27SBarry Smith #undef __FUNCT__
37947c6ae99SBarry Smith #define __FUNCT__ "DMDestroy"
38047c6ae99SBarry Smith /*@
381aa219208SBarry Smith     DMDestroy - Destroys a vector packer or DMDA.
38247c6ae99SBarry Smith 
38347c6ae99SBarry Smith     Collective on DM
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith     Input Parameter:
38647c6ae99SBarry Smith .   dm - the DM object to destroy
38747c6ae99SBarry Smith 
38847c6ae99SBarry Smith     Level: developer
38947c6ae99SBarry Smith 
390e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
39147c6ae99SBarry Smith 
39247c6ae99SBarry Smith @*/
393fcfd50ebSBarry Smith PetscErrorCode  DMDestroy(DM *dm)
39447c6ae99SBarry Smith {
395c83e3748SMatthew G. Knepley   PetscInt       i, cnt = 0;
396dfe15315SJed Brown   DMNamedVecLink nlink,nnext;
39747c6ae99SBarry Smith   PetscErrorCode ierr;
39847c6ae99SBarry Smith 
39947c6ae99SBarry Smith   PetscFunctionBegin;
4006bf464f9SBarry Smith   if (!*dm) PetscFunctionReturn(0);
4016bf464f9SBarry Smith   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
40287e657c6SBarry Smith 
40387e657c6SBarry Smith   /* count all the circular references of DM and its contained Vecs */
404732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
4058865f1eaSKarl Rupp     if ((*dm)->localin[i])  cnt++;
4068865f1eaSKarl Rupp     if ((*dm)->globalin[i]) cnt++;
407732e2eb9SMatthew G Knepley   }
408dfe15315SJed Brown   for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++;
4092348bcf4SPeter Brune   for (nlink=(*dm)->namedlocal; nlink; nlink=nlink->next) cnt++;
4102348bcf4SPeter Brune   if ((*dm)->x) {
4112348bcf4SPeter Brune     DM obj;
4122348bcf4SPeter Brune     ierr = VecGetDM((*dm)->x, &obj);CHKERRQ(ierr);
4132348bcf4SPeter Brune     if (obj == *dm) cnt++;
4142348bcf4SPeter Brune   }
415732e2eb9SMatthew G Knepley 
4166bf464f9SBarry Smith   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
417732e2eb9SMatthew G Knepley   /*
418732e2eb9SMatthew G Knepley      Need this test because the dm references the vectors that
419732e2eb9SMatthew G Knepley      reference the dm, so destroying the dm calls destroy on the
420732e2eb9SMatthew G Knepley      vectors that cause another destroy on the dm
421732e2eb9SMatthew G Knepley   */
4226bf464f9SBarry Smith   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
4236bf464f9SBarry Smith   ((PetscObject) (*dm))->refct = 0;
424732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
4256bf464f9SBarry Smith     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
4266bf464f9SBarry Smith     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
427732e2eb9SMatthew G Knepley   }
428f490541aSPeter Brune   nnext=(*dm)->namedglobal;
4290298fd71SBarry Smith   (*dm)->namedglobal = NULL;
430f490541aSPeter Brune   for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named vectors */
4312348bcf4SPeter Brune     nnext = nlink->next;
4322348bcf4SPeter Brune     if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name);
4332348bcf4SPeter Brune     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
4342348bcf4SPeter Brune     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
4352348bcf4SPeter Brune     ierr = PetscFree(nlink);CHKERRQ(ierr);
4362348bcf4SPeter Brune   }
437f490541aSPeter Brune   nnext=(*dm)->namedlocal;
4380298fd71SBarry Smith   (*dm)->namedlocal = NULL;
439f490541aSPeter Brune   for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named local vectors */
440f490541aSPeter Brune     nnext = nlink->next;
441f490541aSPeter Brune     if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name);
442f490541aSPeter Brune     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
443f490541aSPeter Brune     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
444f490541aSPeter Brune     ierr = PetscFree(nlink);CHKERRQ(ierr);
445f490541aSPeter Brune   }
4462348bcf4SPeter Brune 
447b17ce1afSJed Brown   /* Destroy the list of hooks */
448c833c3b5SJed Brown   {
449c833c3b5SJed Brown     DMCoarsenHookLink link,next;
450b17ce1afSJed Brown     for (link=(*dm)->coarsenhook; link; link=next) {
451b17ce1afSJed Brown       next = link->next;
452b17ce1afSJed Brown       ierr = PetscFree(link);CHKERRQ(ierr);
453b17ce1afSJed Brown     }
4540298fd71SBarry Smith     (*dm)->coarsenhook = NULL;
455c833c3b5SJed Brown   }
456c833c3b5SJed Brown   {
457c833c3b5SJed Brown     DMRefineHookLink link,next;
458c833c3b5SJed Brown     for (link=(*dm)->refinehook; link; link=next) {
459c833c3b5SJed Brown       next = link->next;
460c833c3b5SJed Brown       ierr = PetscFree(link);CHKERRQ(ierr);
461c833c3b5SJed Brown     }
4620298fd71SBarry Smith     (*dm)->refinehook = NULL;
463c833c3b5SJed Brown   }
464be081cd6SPeter Brune   {
465be081cd6SPeter Brune     DMSubDomainHookLink link,next;
466be081cd6SPeter Brune     for (link=(*dm)->subdomainhook; link; link=next) {
467be081cd6SPeter Brune       next = link->next;
468be081cd6SPeter Brune       ierr = PetscFree(link);CHKERRQ(ierr);
469be081cd6SPeter Brune     }
4700298fd71SBarry Smith     (*dm)->subdomainhook = NULL;
471be081cd6SPeter Brune   }
472baf369e7SPeter Brune   {
473baf369e7SPeter Brune     DMGlobalToLocalHookLink link,next;
474baf369e7SPeter Brune     for (link=(*dm)->gtolhook; link; link=next) {
475baf369e7SPeter Brune       next = link->next;
476baf369e7SPeter Brune       ierr = PetscFree(link);CHKERRQ(ierr);
477baf369e7SPeter Brune     }
4780298fd71SBarry Smith     (*dm)->gtolhook = NULL;
479baf369e7SPeter Brune   }
480d4d07f1eSToby Isaac   {
481d4d07f1eSToby Isaac     DMLocalToGlobalHookLink link,next;
482d4d07f1eSToby Isaac     for (link=(*dm)->ltoghook; link; link=next) {
483d4d07f1eSToby Isaac       next = link->next;
484d4d07f1eSToby Isaac       ierr = PetscFree(link);CHKERRQ(ierr);
485d4d07f1eSToby Isaac     }
486d4d07f1eSToby Isaac     (*dm)->ltoghook = NULL;
487d4d07f1eSToby Isaac   }
488aa1993deSMatthew G Knepley   /* Destroy the work arrays */
489aa1993deSMatthew G Knepley   {
490aa1993deSMatthew G Knepley     DMWorkLink link,next;
491aa1993deSMatthew G Knepley     if ((*dm)->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out");
492aa1993deSMatthew G Knepley     for (link=(*dm)->workin; link; link=next) {
493aa1993deSMatthew G Knepley       next = link->next;
494aa1993deSMatthew G Knepley       ierr = PetscFree(link->mem);CHKERRQ(ierr);
495aa1993deSMatthew G Knepley       ierr = PetscFree(link);CHKERRQ(ierr);
496aa1993deSMatthew G Knepley     }
4970298fd71SBarry Smith     (*dm)->workin = NULL;
498aa1993deSMatthew G Knepley   }
499b17ce1afSJed Brown 
50052536dc3SBarry Smith   ierr = PetscObjectDestroy(&(*dm)->dmksp);CHKERRQ(ierr);
50152536dc3SBarry Smith   ierr = PetscObjectDestroy(&(*dm)->dmsnes);CHKERRQ(ierr);
50252536dc3SBarry Smith   ierr = PetscObjectDestroy(&(*dm)->dmts);CHKERRQ(ierr);
50352536dc3SBarry Smith 
5041a266240SBarry Smith   if ((*dm)->ctx && (*dm)->ctxdestroy) {
5051a266240SBarry Smith     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
5061a266240SBarry Smith   }
50787e657c6SBarry Smith   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
50871cd77b2SBarry Smith   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
5094dcab191SBarry Smith   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
5106bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
5116bf464f9SBarry Smith   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
512073dac72SJed Brown   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
51388ed4aceSMatthew G Knepley 
51488ed4aceSMatthew G Knepley   ierr = PetscSectionDestroy(&(*dm)->defaultSection);CHKERRQ(ierr);
51588ed4aceSMatthew G Knepley   ierr = PetscSectionDestroy(&(*dm)->defaultGlobalSection);CHKERRQ(ierr);
5168b1ab98fSJed Brown   ierr = PetscLayoutDestroy(&(*dm)->map);CHKERRQ(ierr);
51788ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&(*dm)->sf);CHKERRQ(ierr);
51888ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&(*dm)->defaultSF);CHKERRQ(ierr);
519af122d2aSMatthew G Knepley 
5206636e97aSMatthew G Knepley   ierr = DMDestroy(&(*dm)->coordinateDM);CHKERRQ(ierr);
5216636e97aSMatthew G Knepley   ierr = VecDestroy(&(*dm)->coordinates);CHKERRQ(ierr);
5226636e97aSMatthew G Knepley   ierr = VecDestroy(&(*dm)->coordinatesLocal);CHKERRQ(ierr);
523c6b900c6SMatthew G. Knepley   ierr = PetscFree((*dm)->maxCell);CHKERRQ(ierr);
524c6b900c6SMatthew G. Knepley   ierr = PetscFree((*dm)->L);CHKERRQ(ierr);
5256636e97aSMatthew G Knepley 
5262764a2aaSMatthew G. Knepley   ierr = PetscDSDestroy(&(*dm)->prob);CHKERRQ(ierr);
52714f150ffSMatthew G. Knepley   ierr = DMDestroy(&(*dm)->dmBC);CHKERRQ(ierr);
528e04113cfSBarry Smith   /* if memory was published with SAWs then destroy it */
529e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*dm);CHKERRQ(ierr);
530732e2eb9SMatthew G Knepley 
5316bf464f9SBarry Smith   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
532435a35e8SMatthew G Knepley   /* We do not destroy (*dm)->data here so that we can reference count backend objects */
533732e2eb9SMatthew G Knepley   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
53447c6ae99SBarry Smith   PetscFunctionReturn(0);
53547c6ae99SBarry Smith }
53647c6ae99SBarry Smith 
53747c6ae99SBarry Smith #undef __FUNCT__
538d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp"
539d7bf68aeSBarry Smith /*@
540d7bf68aeSBarry Smith     DMSetUp - sets up the data structures inside a DM object
541d7bf68aeSBarry Smith 
542d7bf68aeSBarry Smith     Collective on DM
543d7bf68aeSBarry Smith 
544d7bf68aeSBarry Smith     Input Parameter:
545d7bf68aeSBarry Smith .   dm - the DM object to setup
546d7bf68aeSBarry Smith 
547d7bf68aeSBarry Smith     Level: developer
548d7bf68aeSBarry Smith 
549e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
550d7bf68aeSBarry Smith 
551d7bf68aeSBarry Smith @*/
5527087cfbeSBarry Smith PetscErrorCode  DMSetUp(DM dm)
553d7bf68aeSBarry Smith {
554d7bf68aeSBarry Smith   PetscErrorCode ierr;
555d7bf68aeSBarry Smith 
556d7bf68aeSBarry Smith   PetscFunctionBegin;
557171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5588387afaaSJed Brown   if (dm->setupcalled) PetscFunctionReturn(0);
559d7bf68aeSBarry Smith   if (dm->ops->setup) {
560d7bf68aeSBarry Smith     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
561d7bf68aeSBarry Smith   }
5628387afaaSJed Brown   dm->setupcalled = PETSC_TRUE;
563d7bf68aeSBarry Smith   PetscFunctionReturn(0);
564d7bf68aeSBarry Smith }
565d7bf68aeSBarry Smith 
566d7bf68aeSBarry Smith #undef __FUNCT__
567d7bf68aeSBarry Smith #define __FUNCT__ "DMSetFromOptions"
568d7bf68aeSBarry Smith /*@
569d7bf68aeSBarry Smith     DMSetFromOptions - sets parameters in a DM from the options database
570d7bf68aeSBarry Smith 
571d7bf68aeSBarry Smith     Collective on DM
572d7bf68aeSBarry Smith 
573d7bf68aeSBarry Smith     Input Parameter:
574d7bf68aeSBarry Smith .   dm - the DM object to set options for
575d7bf68aeSBarry Smith 
576732e2eb9SMatthew G Knepley     Options Database:
577dd85299cSBarry Smith +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
578dd85299cSBarry Smith .   -dm_vec_type <type>  type of vector to create inside DM
579171400e9SBarry Smith .   -dm_mat_type <type>  type of matrix to create inside DM
580171400e9SBarry Smith -   -dm_coloring_type <global or ghosted>
581732e2eb9SMatthew G Knepley 
582d7bf68aeSBarry Smith     Level: developer
583d7bf68aeSBarry Smith 
584e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
585d7bf68aeSBarry Smith 
586d7bf68aeSBarry Smith @*/
5877087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions(DM dm)
588d7bf68aeSBarry Smith {
5897781c08eSBarry Smith   char           typeName[256];
590ca266f36SBarry Smith   PetscBool      flg;
591d7bf68aeSBarry Smith   PetscErrorCode ierr;
592d7bf68aeSBarry Smith 
593d7bf68aeSBarry Smith   PetscFunctionBegin;
594171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5953194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
5960298fd71SBarry Smith   ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,NULL);CHKERRQ(ierr);
597a264d7a6SBarry Smith   ierr = PetscOptionsFList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
598f9ba7244SBarry Smith   if (flg) {
599f9ba7244SBarry Smith     ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
600f9ba7244SBarry Smith   }
601a264d7a6SBarry Smith   ierr = PetscOptionsFList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype ? dm->mattype : typeName,typeName,sizeof(typeName),&flg);CHKERRQ(ierr);
602073dac72SJed Brown   if (flg) {
603521d9a4cSLisandro Dalcin     ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr);
604073dac72SJed Brown   }
6050298fd71SBarry Smith   ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,NULL);CHKERRQ(ierr);
606f9ba7244SBarry Smith   if (dm->ops->setfromoptions) {
607f9ba7244SBarry Smith     ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
608f9ba7244SBarry Smith   }
609f9ba7244SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
610f9ba7244SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
61182fcb398SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
612d7bf68aeSBarry Smith   PetscFunctionReturn(0);
613d7bf68aeSBarry Smith }
614d7bf68aeSBarry Smith 
615d7bf68aeSBarry Smith #undef __FUNCT__
61647c6ae99SBarry Smith #define __FUNCT__ "DMView"
617fc9bc008SSatish Balay /*@C
618aa219208SBarry Smith     DMView - Views a vector packer or DMDA.
61947c6ae99SBarry Smith 
62047c6ae99SBarry Smith     Collective on DM
62147c6ae99SBarry Smith 
62247c6ae99SBarry Smith     Input Parameter:
62347c6ae99SBarry Smith +   dm - the DM object to view
62447c6ae99SBarry Smith -   v - the viewer
62547c6ae99SBarry Smith 
62647c6ae99SBarry Smith     Level: developer
62747c6ae99SBarry Smith 
628e727c939SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
62947c6ae99SBarry Smith 
63047c6ae99SBarry Smith @*/
6317087cfbeSBarry Smith PetscErrorCode  DMView(DM dm,PetscViewer v)
63247c6ae99SBarry Smith {
63347c6ae99SBarry Smith   PetscErrorCode ierr;
63432c0f0efSBarry Smith   PetscBool      isbinary;
63547c6ae99SBarry Smith 
63647c6ae99SBarry Smith   PetscFunctionBegin;
637171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6383014e516SBarry Smith   if (!v) {
639ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm),&v);CHKERRQ(ierr);
6403014e516SBarry Smith   }
64198c3331eSBarry Smith   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)dm,v);CHKERRQ(ierr);
64232c0f0efSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
64332c0f0efSBarry Smith   if (isbinary) {
64455849f57SBarry Smith     PetscInt classid = DM_FILE_CLASSID;
64532c0f0efSBarry Smith     char     type[256];
64632c0f0efSBarry Smith 
64732c0f0efSBarry Smith     ierr = PetscViewerBinaryWrite(v,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
64832c0f0efSBarry Smith     ierr = PetscStrncpy(type,((PetscObject)dm)->type_name,256);CHKERRQ(ierr);
64932c0f0efSBarry Smith     ierr = PetscViewerBinaryWrite(v,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
65032c0f0efSBarry Smith   }
6510c010503SBarry Smith   if (dm->ops->view) {
6520c010503SBarry Smith     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
65347c6ae99SBarry Smith   }
65447c6ae99SBarry Smith   PetscFunctionReturn(0);
65547c6ae99SBarry Smith }
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith #undef __FUNCT__
65847c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector"
65947c6ae99SBarry Smith /*@
660aa219208SBarry Smith     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
66147c6ae99SBarry Smith 
66247c6ae99SBarry Smith     Collective on DM
66347c6ae99SBarry Smith 
66447c6ae99SBarry Smith     Input Parameter:
66547c6ae99SBarry Smith .   dm - the DM object
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith     Output Parameter:
66847c6ae99SBarry Smith .   vec - the global vector
66947c6ae99SBarry Smith 
670073dac72SJed Brown     Level: beginner
67147c6ae99SBarry Smith 
672e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
67347c6ae99SBarry Smith 
67447c6ae99SBarry Smith @*/
6757087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
67647c6ae99SBarry Smith {
67747c6ae99SBarry Smith   PetscErrorCode ierr;
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith   PetscFunctionBegin;
680171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
68147c6ae99SBarry Smith   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
68247c6ae99SBarry Smith   PetscFunctionReturn(0);
68347c6ae99SBarry Smith }
68447c6ae99SBarry Smith 
68547c6ae99SBarry Smith #undef __FUNCT__
68647c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector"
68747c6ae99SBarry Smith /*@
688aa219208SBarry Smith     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
68947c6ae99SBarry Smith 
69047c6ae99SBarry Smith     Not Collective
69147c6ae99SBarry Smith 
69247c6ae99SBarry Smith     Input Parameter:
69347c6ae99SBarry Smith .   dm - the DM object
69447c6ae99SBarry Smith 
69547c6ae99SBarry Smith     Output Parameter:
69647c6ae99SBarry Smith .   vec - the local vector
69747c6ae99SBarry Smith 
698073dac72SJed Brown     Level: beginner
69947c6ae99SBarry Smith 
700e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
70147c6ae99SBarry Smith 
70247c6ae99SBarry Smith @*/
7037087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
70447c6ae99SBarry Smith {
70547c6ae99SBarry Smith   PetscErrorCode ierr;
70647c6ae99SBarry Smith 
70747c6ae99SBarry Smith   PetscFunctionBegin;
708171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
70947c6ae99SBarry Smith   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
71047c6ae99SBarry Smith   PetscFunctionReturn(0);
71147c6ae99SBarry Smith }
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith #undef __FUNCT__
7141411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping"
7151411c6eeSJed Brown /*@
7161411c6eeSJed Brown    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
7171411c6eeSJed Brown 
7181411c6eeSJed Brown    Collective on DM
7191411c6eeSJed Brown 
7201411c6eeSJed Brown    Input Parameter:
7211411c6eeSJed Brown .  dm - the DM that provides the mapping
7221411c6eeSJed Brown 
7231411c6eeSJed Brown    Output Parameter:
7241411c6eeSJed Brown .  ltog - the mapping
7251411c6eeSJed Brown 
7261411c6eeSJed Brown    Level: intermediate
7271411c6eeSJed Brown 
7281411c6eeSJed Brown    Notes:
7291411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMapping() or
7301411c6eeSJed Brown    MatSetLocalToGlobalMapping().
7311411c6eeSJed Brown 
732fc31e74dSBarry Smith .seealso: DMCreateLocalVector()
7331411c6eeSJed Brown @*/
7347087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
7351411c6eeSJed Brown {
7361411c6eeSJed Brown   PetscErrorCode ierr;
7371411c6eeSJed Brown 
7381411c6eeSJed Brown   PetscFunctionBegin;
7391411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7401411c6eeSJed Brown   PetscValidPointer(ltog,2);
7411411c6eeSJed Brown   if (!dm->ltogmap) {
74237d0c07bSMatthew G Knepley     PetscSection section, sectionGlobal;
74337d0c07bSMatthew G Knepley 
74437d0c07bSMatthew G Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
74537d0c07bSMatthew G Knepley     if (section) {
74637d0c07bSMatthew G Knepley       PetscInt *ltog;
74737d0c07bSMatthew G Knepley       PetscInt pStart, pEnd, size, p, l;
74837d0c07bSMatthew G Knepley 
74937d0c07bSMatthew G Knepley       ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
75037d0c07bSMatthew G Knepley       ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
75137d0c07bSMatthew G Knepley       ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
752785e854fSJed Brown       ierr = PetscMalloc1(size, &ltog);CHKERRQ(ierr); /* We want the local+overlap size */
75337d0c07bSMatthew G Knepley       for (p = pStart, l = 0; p < pEnd; ++p) {
75437d0c07bSMatthew G Knepley         PetscInt dof, off, c;
75537d0c07bSMatthew G Knepley 
75637d0c07bSMatthew G Knepley         /* Should probably use constrained dofs */
75737d0c07bSMatthew G Knepley         ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
75837d0c07bSMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &off);CHKERRQ(ierr);
75937d0c07bSMatthew G Knepley         for (c = 0; c < dof; ++c, ++l) {
76037d0c07bSMatthew G Knepley           ltog[l] = off+c;
76137d0c07bSMatthew G Knepley         }
76237d0c07bSMatthew G Knepley       }
763f0413b6fSBarry Smith       ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, 1,size, ltog, PETSC_OWN_POINTER, &dm->ltogmap);CHKERRQ(ierr);
7643bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)dm, (PetscObject)dm->ltogmap);CHKERRQ(ierr);
76537d0c07bSMatthew G Knepley     } else {
766184d77edSJed Brown       if (!dm->ops->getlocaltoglobalmapping) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
767184d77edSJed Brown       ierr = (*dm->ops->getlocaltoglobalmapping)(dm);CHKERRQ(ierr);
7681411c6eeSJed Brown     }
76937d0c07bSMatthew G Knepley   }
7701411c6eeSJed Brown   *ltog = dm->ltogmap;
7711411c6eeSJed Brown   PetscFunctionReturn(0);
7721411c6eeSJed Brown }
7731411c6eeSJed Brown 
7741411c6eeSJed Brown #undef __FUNCT__
7751411c6eeSJed Brown #define __FUNCT__ "DMGetBlockSize"
7761411c6eeSJed Brown /*@
7771411c6eeSJed Brown    DMGetBlockSize - Gets the inherent block size associated with a DM
7781411c6eeSJed Brown 
7791411c6eeSJed Brown    Not Collective
7801411c6eeSJed Brown 
7811411c6eeSJed Brown    Input Parameter:
7821411c6eeSJed Brown .  dm - the DM with block structure
7831411c6eeSJed Brown 
7841411c6eeSJed Brown    Output Parameter:
7851411c6eeSJed Brown .  bs - the block size, 1 implies no exploitable block structure
7861411c6eeSJed Brown 
7871411c6eeSJed Brown    Level: intermediate
7881411c6eeSJed Brown 
789fc31e74dSBarry Smith .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMapping()
7901411c6eeSJed Brown @*/
7917087cfbeSBarry Smith PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
7921411c6eeSJed Brown {
7931411c6eeSJed Brown   PetscFunctionBegin;
7941411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7951411c6eeSJed Brown   PetscValidPointer(bs,2);
7961411c6eeSJed Brown   if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet");
7971411c6eeSJed Brown   *bs = dm->bs;
7981411c6eeSJed Brown   PetscFunctionReturn(0);
7991411c6eeSJed Brown }
8001411c6eeSJed Brown 
8011411c6eeSJed Brown #undef __FUNCT__
802e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation"
80347c6ae99SBarry Smith /*@
804e727c939SJed Brown     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
80547c6ae99SBarry Smith 
80647c6ae99SBarry Smith     Collective on DM
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith     Input Parameter:
80947c6ae99SBarry Smith +   dm1 - the DM object
81047c6ae99SBarry Smith -   dm2 - the second, finer DM object
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith     Output Parameter:
81347c6ae99SBarry Smith +  mat - the interpolation
81447c6ae99SBarry Smith -  vec - the scaling (optional)
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith     Level: developer
81747c6ae99SBarry Smith 
81885afcc9aSBarry Smith     Notes:  For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by
81985afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
820d52bd9f3SBarry Smith 
8211f588964SMatthew G Knepley         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMGetCoordinateDM()) to interpolate the mesh coordinate vectors
822d52bd9f3SBarry Smith         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
82385afcc9aSBarry Smith 
82485afcc9aSBarry Smith 
825e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith @*/
828e727c939SJed Brown PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
82947c6ae99SBarry Smith {
83047c6ae99SBarry Smith   PetscErrorCode ierr;
83147c6ae99SBarry Smith 
83247c6ae99SBarry Smith   PetscFunctionBegin;
833171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
834171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
83525296bd5SBarry Smith   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
83647c6ae99SBarry Smith   PetscFunctionReturn(0);
83747c6ae99SBarry Smith }
83847c6ae99SBarry Smith 
83947c6ae99SBarry Smith #undef __FUNCT__
840e727c939SJed Brown #define __FUNCT__ "DMCreateInjection"
84147c6ae99SBarry Smith /*@
842e727c939SJed Brown     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith     Collective on DM
84547c6ae99SBarry Smith 
84647c6ae99SBarry Smith     Input Parameter:
84747c6ae99SBarry Smith +   dm1 - the DM object
84847c6ae99SBarry Smith -   dm2 - the second, finer DM object
84947c6ae99SBarry Smith 
85047c6ae99SBarry Smith     Output Parameter:
85147c6ae99SBarry Smith .   ctx - the injection
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith     Level: developer
85447c6ae99SBarry Smith 
85585afcc9aSBarry Smith    Notes:  For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by
85685afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
85785afcc9aSBarry Smith 
858e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith @*/
861e727c939SJed Brown PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
86247c6ae99SBarry Smith {
86347c6ae99SBarry Smith   PetscErrorCode ierr;
86447c6ae99SBarry Smith 
86547c6ae99SBarry Smith   PetscFunctionBegin;
866171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
867171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
86847c6ae99SBarry Smith   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
86947c6ae99SBarry Smith   PetscFunctionReturn(0);
87047c6ae99SBarry Smith }
87147c6ae99SBarry Smith 
87247c6ae99SBarry Smith #undef __FUNCT__
873e727c939SJed Brown #define __FUNCT__ "DMCreateColoring"
874b412c318SBarry Smith /*@
875b412c318SBarry Smith     DMCreateColoring - Gets coloring for a DM
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith     Collective on DM
87847c6ae99SBarry Smith 
87947c6ae99SBarry Smith     Input Parameter:
88047c6ae99SBarry Smith +   dm - the DM object
881b412c318SBarry Smith -   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith     Output Parameter:
88447c6ae99SBarry Smith .   coloring - the coloring
88547c6ae99SBarry Smith 
88647c6ae99SBarry Smith     Level: developer
88747c6ae99SBarry Smith 
888b412c318SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix(), DMSetMatType()
88947c6ae99SBarry Smith 
890aab9d709SJed Brown @*/
891b412c318SBarry Smith PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,ISColoring *coloring)
89247c6ae99SBarry Smith {
89347c6ae99SBarry Smith   PetscErrorCode ierr;
89447c6ae99SBarry Smith 
89547c6ae99SBarry Smith   PetscFunctionBegin;
896171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
897ce94432eSBarry Smith   if (!dm->ops->getcoloring) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No coloring for this type of DM yet");
898b412c318SBarry Smith   ierr = (*dm->ops->getcoloring)(dm,ctype,coloring);CHKERRQ(ierr);
89947c6ae99SBarry Smith   PetscFunctionReturn(0);
90047c6ae99SBarry Smith }
90147c6ae99SBarry Smith 
90247c6ae99SBarry Smith #undef __FUNCT__
903950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix"
904b412c318SBarry Smith /*@
905950540a4SJed Brown     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
90647c6ae99SBarry Smith 
90747c6ae99SBarry Smith     Collective on DM
90847c6ae99SBarry Smith 
90947c6ae99SBarry Smith     Input Parameter:
910b412c318SBarry Smith .   dm - the DM object
91147c6ae99SBarry Smith 
91247c6ae99SBarry Smith     Output Parameter:
91347c6ae99SBarry Smith .   mat - the empty Jacobian
91447c6ae99SBarry Smith 
915073dac72SJed Brown     Level: beginner
91647c6ae99SBarry Smith 
91794013140SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
91894013140SBarry Smith        do not need to do it yourself.
91994013140SBarry Smith 
92094013140SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
921aa219208SBarry Smith        the nonzero pattern call DMDASetMatPreallocateOnly()
92294013140SBarry Smith 
92394013140SBarry Smith        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
92494013140SBarry Smith        internally by PETSc.
92594013140SBarry Smith 
92694013140SBarry Smith        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
927aa219208SBarry Smith        the indices for the global numbering for DMDAs which is complicated.
92894013140SBarry Smith 
929b412c318SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMSetMatType()
93047c6ae99SBarry Smith 
931aab9d709SJed Brown @*/
932b412c318SBarry Smith PetscErrorCode  DMCreateMatrix(DM dm,Mat *mat)
93347c6ae99SBarry Smith {
93447c6ae99SBarry Smith   PetscErrorCode ierr;
93547c6ae99SBarry Smith 
93647c6ae99SBarry Smith   PetscFunctionBegin;
937171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
938607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
939c7b7c8a4SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
940c7b7c8a4SJed Brown   PetscValidPointer(mat,3);
941b412c318SBarry Smith   ierr = (*dm->ops->creatematrix)(dm,mat);CHKERRQ(ierr);
94247c6ae99SBarry Smith   PetscFunctionReturn(0);
94347c6ae99SBarry Smith }
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith #undef __FUNCT__
946732e2eb9SMatthew G Knepley #define __FUNCT__ "DMSetMatrixPreallocateOnly"
947732e2eb9SMatthew G Knepley /*@
948950540a4SJed Brown   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
949732e2eb9SMatthew G Knepley     preallocated but the nonzero structure and zero values will not be set.
950732e2eb9SMatthew G Knepley 
951732e2eb9SMatthew G Knepley   Logically Collective on DMDA
952732e2eb9SMatthew G Knepley 
953732e2eb9SMatthew G Knepley   Input Parameter:
954732e2eb9SMatthew G Knepley + dm - the DM
955732e2eb9SMatthew G Knepley - only - PETSC_TRUE if only want preallocation
956732e2eb9SMatthew G Knepley 
957732e2eb9SMatthew G Knepley   Level: developer
958950540a4SJed Brown .seealso DMCreateMatrix()
959732e2eb9SMatthew G Knepley @*/
960732e2eb9SMatthew G Knepley PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
961732e2eb9SMatthew G Knepley {
962732e2eb9SMatthew G Knepley   PetscFunctionBegin;
963732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
964732e2eb9SMatthew G Knepley   dm->prealloc_only = only;
965732e2eb9SMatthew G Knepley   PetscFunctionReturn(0);
966732e2eb9SMatthew G Knepley }
967732e2eb9SMatthew G Knepley 
968732e2eb9SMatthew G Knepley #undef __FUNCT__
969a89ea682SMatthew G Knepley #define __FUNCT__ "DMGetWorkArray"
970a89ea682SMatthew G Knepley /*@C
971aa1993deSMatthew G Knepley   DMGetWorkArray - Gets a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray()
972a89ea682SMatthew G Knepley 
973a89ea682SMatthew G Knepley   Not Collective
974a89ea682SMatthew G Knepley 
975a89ea682SMatthew G Knepley   Input Parameters:
976a89ea682SMatthew G Knepley + dm - the DM object
977aa1993deSMatthew G Knepley . count - The minium size
978aa1993deSMatthew G Knepley - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT)
979a89ea682SMatthew G Knepley 
980a89ea682SMatthew G Knepley   Output Parameter:
981a89ea682SMatthew G Knepley . array - the work array
982a89ea682SMatthew G Knepley 
983a89ea682SMatthew G Knepley   Level: developer
984a89ea682SMatthew G Knepley 
985a89ea682SMatthew G Knepley .seealso DMDestroy(), DMCreate()
986a89ea682SMatthew G Knepley @*/
987aa1993deSMatthew G Knepley PetscErrorCode DMGetWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem)
988a89ea682SMatthew G Knepley {
989a89ea682SMatthew G Knepley   PetscErrorCode ierr;
990aa1993deSMatthew G Knepley   DMWorkLink     link;
991aa1993deSMatthew G Knepley   size_t         size;
992a89ea682SMatthew G Knepley 
993a89ea682SMatthew G Knepley   PetscFunctionBegin;
994a89ea682SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
995aa1993deSMatthew G Knepley   PetscValidPointer(mem,4);
996aa1993deSMatthew G Knepley   if (dm->workin) {
997aa1993deSMatthew G Knepley     link       = dm->workin;
998aa1993deSMatthew G Knepley     dm->workin = dm->workin->next;
999aa1993deSMatthew G Knepley   } else {
1000b00a9115SJed Brown     ierr = PetscNewLog(dm,&link);CHKERRQ(ierr);
1001a89ea682SMatthew G Knepley   }
1002aa1993deSMatthew G Knepley   ierr = PetscDataTypeGetSize(dtype,&size);CHKERRQ(ierr);
1003aa1993deSMatthew G Knepley   if (size*count > link->bytes) {
1004aa1993deSMatthew G Knepley     ierr        = PetscFree(link->mem);CHKERRQ(ierr);
1005aa1993deSMatthew G Knepley     ierr        = PetscMalloc(size*count,&link->mem);CHKERRQ(ierr);
1006aa1993deSMatthew G Knepley     link->bytes = size*count;
1007aa1993deSMatthew G Knepley   }
1008aa1993deSMatthew G Knepley   link->next   = dm->workout;
1009aa1993deSMatthew G Knepley   dm->workout  = link;
1010aa1993deSMatthew G Knepley   *(void**)mem = link->mem;
1011a89ea682SMatthew G Knepley   PetscFunctionReturn(0);
1012a89ea682SMatthew G Knepley }
1013a89ea682SMatthew G Knepley 
1014aa1993deSMatthew G Knepley #undef __FUNCT__
1015aa1993deSMatthew G Knepley #define __FUNCT__ "DMRestoreWorkArray"
1016aa1993deSMatthew G Knepley /*@C
1017aa1993deSMatthew G Knepley   DMRestoreWorkArray - Restores a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray()
1018aa1993deSMatthew G Knepley 
1019aa1993deSMatthew G Knepley   Not Collective
1020aa1993deSMatthew G Knepley 
1021aa1993deSMatthew G Knepley   Input Parameters:
1022aa1993deSMatthew G Knepley + dm - the DM object
1023aa1993deSMatthew G Knepley . count - The minium size
1024aa1993deSMatthew G Knepley - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT)
1025aa1993deSMatthew G Knepley 
1026aa1993deSMatthew G Knepley   Output Parameter:
1027aa1993deSMatthew G Knepley . array - the work array
1028aa1993deSMatthew G Knepley 
1029aa1993deSMatthew G Knepley   Level: developer
1030aa1993deSMatthew G Knepley 
1031aa1993deSMatthew G Knepley .seealso DMDestroy(), DMCreate()
1032aa1993deSMatthew G Knepley @*/
1033aa1993deSMatthew G Knepley PetscErrorCode DMRestoreWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem)
1034aa1993deSMatthew G Knepley {
1035aa1993deSMatthew G Knepley   DMWorkLink *p,link;
1036aa1993deSMatthew G Knepley 
1037aa1993deSMatthew G Knepley   PetscFunctionBegin;
1038aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1039aa1993deSMatthew G Knepley   PetscValidPointer(mem,4);
1040aa1993deSMatthew G Knepley   for (p=&dm->workout; (link=*p); p=&link->next) {
1041aa1993deSMatthew G Knepley     if (link->mem == *(void**)mem) {
1042aa1993deSMatthew G Knepley       *p           = link->next;
1043aa1993deSMatthew G Knepley       link->next   = dm->workin;
1044aa1993deSMatthew G Knepley       dm->workin   = link;
10450298fd71SBarry Smith       *(void**)mem = NULL;
1046aa1993deSMatthew G Knepley       PetscFunctionReturn(0);
1047aa1993deSMatthew G Knepley     }
1048aa1993deSMatthew G Knepley   }
1049aa1993deSMatthew G Knepley   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out");
1050aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
1051aa1993deSMatthew G Knepley }
1052e7c4fc90SDmitry Karpeev 
1053e7c4fc90SDmitry Karpeev #undef __FUNCT__
1054435a35e8SMatthew G Knepley #define __FUNCT__ "DMSetNullSpaceConstructor"
1055435a35e8SMatthew G Knepley PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt field, MatNullSpace *nullSpace))
1056435a35e8SMatthew G Knepley {
1057435a35e8SMatthew G Knepley   PetscFunctionBegin;
1058435a35e8SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
105982f516ccSBarry Smith   if (field >= 10) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %d >= 10 fields", field);
1060435a35e8SMatthew G Knepley   dm->nullspaceConstructors[field] = nullsp;
1061435a35e8SMatthew G Knepley   PetscFunctionReturn(0);
1062435a35e8SMatthew G Knepley }
1063435a35e8SMatthew G Knepley 
1064435a35e8SMatthew G Knepley #undef __FUNCT__
10654d343eeaSMatthew G Knepley #define __FUNCT__ "DMCreateFieldIS"
10664f3b5142SJed Brown /*@C
10674d343eeaSMatthew G Knepley   DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field
10684d343eeaSMatthew G Knepley 
10694d343eeaSMatthew G Knepley   Not collective
10704d343eeaSMatthew G Knepley 
10714d343eeaSMatthew G Knepley   Input Parameter:
10724d343eeaSMatthew G Knepley . dm - the DM object
10734d343eeaSMatthew G Knepley 
10744d343eeaSMatthew G Knepley   Output Parameters:
10750298fd71SBarry Smith + numFields  - The number of fields (or NULL if not requested)
10760298fd71SBarry Smith . fieldNames - The name for each field (or NULL if not requested)
10770298fd71SBarry Smith - fields     - The global indices for each field (or NULL if not requested)
10784d343eeaSMatthew G Knepley 
10794d343eeaSMatthew G Knepley   Level: intermediate
10804d343eeaSMatthew G Knepley 
108121c9b008SJed Brown   Notes:
108221c9b008SJed Brown   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
108321c9b008SJed Brown   PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with
108421c9b008SJed Brown   PetscFree().
108521c9b008SJed Brown 
10864d343eeaSMatthew G Knepley .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
10874d343eeaSMatthew G Knepley @*/
108837d0c07bSMatthew G Knepley PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
10894d343eeaSMatthew G Knepley {
109037d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
10914d343eeaSMatthew G Knepley   PetscErrorCode ierr;
10924d343eeaSMatthew G Knepley 
10934d343eeaSMatthew G Knepley   PetscFunctionBegin;
10944d343eeaSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
109569ca1f37SDmitry Karpeev   if (numFields) {
109669ca1f37SDmitry Karpeev     PetscValidPointer(numFields,2);
109769ca1f37SDmitry Karpeev     *numFields = 0;
109869ca1f37SDmitry Karpeev   }
109937d0c07bSMatthew G Knepley   if (fieldNames) {
110037d0c07bSMatthew G Knepley     PetscValidPointer(fieldNames,3);
11010298fd71SBarry Smith     *fieldNames = NULL;
110269ca1f37SDmitry Karpeev   }
110369ca1f37SDmitry Karpeev   if (fields) {
110469ca1f37SDmitry Karpeev     PetscValidPointer(fields,4);
11050298fd71SBarry Smith     *fields = NULL;
110669ca1f37SDmitry Karpeev   }
110737d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
110837d0c07bSMatthew G Knepley   if (section) {
110937d0c07bSMatthew G Knepley     PetscInt *fieldSizes, **fieldIndices;
111037d0c07bSMatthew G Knepley     PetscInt nF, f, pStart, pEnd, p;
111137d0c07bSMatthew G Knepley 
111237d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
111337d0c07bSMatthew G Knepley     ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
1114dcca6d9dSJed Brown     ierr = PetscMalloc2(nF,&fieldSizes,nF,&fieldIndices);CHKERRQ(ierr);
111537d0c07bSMatthew G Knepley     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
111637d0c07bSMatthew G Knepley     for (f = 0; f < nF; ++f) {
111737d0c07bSMatthew G Knepley       fieldSizes[f] = 0;
111837d0c07bSMatthew G Knepley     }
111937d0c07bSMatthew G Knepley     for (p = pStart; p < pEnd; ++p) {
112037d0c07bSMatthew G Knepley       PetscInt gdof;
112137d0c07bSMatthew G Knepley 
112237d0c07bSMatthew G Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
112337d0c07bSMatthew G Knepley       if (gdof > 0) {
112437d0c07bSMatthew G Knepley         for (f = 0; f < nF; ++f) {
112537d0c07bSMatthew G Knepley           PetscInt fdof, fcdof;
112637d0c07bSMatthew G Knepley 
112737d0c07bSMatthew G Knepley           ierr           = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
112837d0c07bSMatthew G Knepley           ierr           = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
112937d0c07bSMatthew G Knepley           fieldSizes[f] += fdof-fcdof;
113037d0c07bSMatthew G Knepley         }
113137d0c07bSMatthew G Knepley       }
113237d0c07bSMatthew G Knepley     }
113337d0c07bSMatthew G Knepley     for (f = 0; f < nF; ++f) {
1134785e854fSJed Brown       ierr          = PetscMalloc1(fieldSizes[f], &fieldIndices[f]);CHKERRQ(ierr);
113537d0c07bSMatthew G Knepley       fieldSizes[f] = 0;
113637d0c07bSMatthew G Knepley     }
113737d0c07bSMatthew G Knepley     for (p = pStart; p < pEnd; ++p) {
113837d0c07bSMatthew G Knepley       PetscInt gdof, goff;
113937d0c07bSMatthew G Knepley 
114037d0c07bSMatthew G Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
114137d0c07bSMatthew G Knepley       if (gdof > 0) {
114237d0c07bSMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
114337d0c07bSMatthew G Knepley         for (f = 0; f < nF; ++f) {
114437d0c07bSMatthew G Knepley           PetscInt fdof, fcdof, fc;
114537d0c07bSMatthew G Knepley 
114637d0c07bSMatthew G Knepley           ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
114737d0c07bSMatthew G Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
114837d0c07bSMatthew G Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++fieldSizes[f]) {
114937d0c07bSMatthew G Knepley             fieldIndices[f][fieldSizes[f]] = goff++;
115037d0c07bSMatthew G Knepley           }
115137d0c07bSMatthew G Knepley         }
115237d0c07bSMatthew G Knepley       }
115337d0c07bSMatthew G Knepley     }
11548865f1eaSKarl Rupp     if (numFields) *numFields = nF;
115537d0c07bSMatthew G Knepley     if (fieldNames) {
1156785e854fSJed Brown       ierr = PetscMalloc1(nF, fieldNames);CHKERRQ(ierr);
115737d0c07bSMatthew G Knepley       for (f = 0; f < nF; ++f) {
115837d0c07bSMatthew G Knepley         const char *fieldName;
115937d0c07bSMatthew G Knepley 
116037d0c07bSMatthew G Knepley         ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr);
116137d0c07bSMatthew G Knepley         ierr = PetscStrallocpy(fieldName, (char**) &(*fieldNames)[f]);CHKERRQ(ierr);
116237d0c07bSMatthew G Knepley       }
116337d0c07bSMatthew G Knepley     }
116437d0c07bSMatthew G Knepley     if (fields) {
1165785e854fSJed Brown       ierr = PetscMalloc1(nF, fields);CHKERRQ(ierr);
116637d0c07bSMatthew G Knepley       for (f = 0; f < nF; ++f) {
116782f516ccSBarry Smith         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]);CHKERRQ(ierr);
116837d0c07bSMatthew G Knepley       }
116937d0c07bSMatthew G Knepley     }
117037d0c07bSMatthew G Knepley     ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr);
11718865f1eaSKarl Rupp   } else if (dm->ops->createfieldis) {
11728865f1eaSKarl Rupp     ierr = (*dm->ops->createfieldis)(dm, numFields, fieldNames, fields);CHKERRQ(ierr);
117369ca1f37SDmitry Karpeev   }
11744d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
11754d343eeaSMatthew G Knepley }
11764d343eeaSMatthew G Knepley 
117716621825SDmitry Karpeev 
117816621825SDmitry Karpeev #undef __FUNCT__
117916621825SDmitry Karpeev #define __FUNCT__ "DMCreateFieldDecomposition"
118016621825SDmitry Karpeev /*@C
118116621825SDmitry Karpeev   DMCreateFieldDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems
118216621825SDmitry Karpeev                           corresponding to different fields: each IS contains the global indices of the dofs of the
118316621825SDmitry Karpeev                           corresponding field. The optional list of DMs define the DM for each subproblem.
1184e7c4fc90SDmitry Karpeev                           Generalizes DMCreateFieldIS().
1185e7c4fc90SDmitry Karpeev 
1186e7c4fc90SDmitry Karpeev   Not collective
1187e7c4fc90SDmitry Karpeev 
1188e7c4fc90SDmitry Karpeev   Input Parameter:
1189e7c4fc90SDmitry Karpeev . dm - the DM object
1190e7c4fc90SDmitry Karpeev 
1191e7c4fc90SDmitry Karpeev   Output Parameters:
11920298fd71SBarry Smith + len       - The number of subproblems in the field decomposition (or NULL if not requested)
11930298fd71SBarry Smith . namelist  - The name for each field (or NULL if not requested)
11940298fd71SBarry Smith . islist    - The global indices for each field (or NULL if not requested)
11950298fd71SBarry Smith - dmlist    - The DMs for each field subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined)
1196e7c4fc90SDmitry Karpeev 
1197e7c4fc90SDmitry Karpeev   Level: intermediate
1198e7c4fc90SDmitry Karpeev 
1199e7c4fc90SDmitry Karpeev   Notes:
1200e7c4fc90SDmitry Karpeev   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
1201e7c4fc90SDmitry Karpeev   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
1202e7c4fc90SDmitry Karpeev   and all of the arrays should be freed with PetscFree().
1203e7c4fc90SDmitry Karpeev 
1204e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1205e7c4fc90SDmitry Karpeev @*/
120616621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1207e7c4fc90SDmitry Karpeev {
1208e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
1209e7c4fc90SDmitry Karpeev 
1210e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
1211e7c4fc90SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12128865f1eaSKarl Rupp   if (len) {
12138865f1eaSKarl Rupp     PetscValidPointer(len,2);
12148865f1eaSKarl Rupp     *len = 0;
12158865f1eaSKarl Rupp   }
12168865f1eaSKarl Rupp   if (namelist) {
12178865f1eaSKarl Rupp     PetscValidPointer(namelist,3);
12188865f1eaSKarl Rupp     *namelist = 0;
12198865f1eaSKarl Rupp   }
12208865f1eaSKarl Rupp   if (islist) {
12218865f1eaSKarl Rupp     PetscValidPointer(islist,4);
12228865f1eaSKarl Rupp     *islist = 0;
12238865f1eaSKarl Rupp   }
12248865f1eaSKarl Rupp   if (dmlist) {
12258865f1eaSKarl Rupp     PetscValidPointer(dmlist,5);
12268865f1eaSKarl Rupp     *dmlist = 0;
12278865f1eaSKarl Rupp   }
1228f3f0edfdSDmitry Karpeev   /*
1229f3f0edfdSDmitry Karpeev    Is it a good idea to apply the following check across all impls?
1230f3f0edfdSDmitry Karpeev    Perhaps some impls can have a well-defined decomposition before DMSetUp?
1231f3f0edfdSDmitry Karpeev    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
1232f3f0edfdSDmitry Karpeev    */
1233ce94432eSBarry Smith   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
123416621825SDmitry Karpeev   if (!dm->ops->createfielddecomposition) {
1235435a35e8SMatthew G Knepley     PetscSection section;
1236435a35e8SMatthew G Knepley     PetscInt     numFields, f;
1237435a35e8SMatthew G Knepley 
1238435a35e8SMatthew G Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1239435a35e8SMatthew G Knepley     if (section) {ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);}
1240435a35e8SMatthew G Knepley     if (section && numFields && dm->ops->createsubdm) {
1241435a35e8SMatthew G Knepley       *len = numFields;
124203dc3394SMatthew G. Knepley       if (namelist) {ierr = PetscMalloc1(numFields,namelist);CHKERRQ(ierr);}
124303dc3394SMatthew G. Knepley       if (islist)   {ierr = PetscMalloc1(numFields,islist);CHKERRQ(ierr);}
124403dc3394SMatthew G. Knepley       if (dmlist)   {ierr = PetscMalloc1(numFields,dmlist);CHKERRQ(ierr);}
1245435a35e8SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
1246435a35e8SMatthew G Knepley         const char *fieldName;
1247435a35e8SMatthew G Knepley 
124803dc3394SMatthew G. Knepley         ierr = DMCreateSubDM(dm, 1, &f, islist ? &(*islist)[f] : NULL, dmlist ? &(*dmlist)[f] : NULL);CHKERRQ(ierr);
124903dc3394SMatthew G. Knepley         if (namelist) {
1250435a35e8SMatthew G Knepley           ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr);
1251435a35e8SMatthew G Knepley           ierr = PetscStrallocpy(fieldName, (char**) &(*namelist)[f]);CHKERRQ(ierr);
1252435a35e8SMatthew G Knepley         }
125303dc3394SMatthew G. Knepley       }
1254435a35e8SMatthew G Knepley     } else {
125569ca1f37SDmitry Karpeev       ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr);
1256e7c4fc90SDmitry Karpeev       /* By default there are no DMs associated with subproblems. */
12570298fd71SBarry Smith       if (dmlist) *dmlist = NULL;
1258e7c4fc90SDmitry Karpeev     }
12598865f1eaSKarl Rupp   } else {
126016621825SDmitry Karpeev     ierr = (*dm->ops->createfielddecomposition)(dm,len,namelist,islist,dmlist);CHKERRQ(ierr);
126116621825SDmitry Karpeev   }
126216621825SDmitry Karpeev   PetscFunctionReturn(0);
126316621825SDmitry Karpeev }
126416621825SDmitry Karpeev 
126516621825SDmitry Karpeev #undef __FUNCT__
1266435a35e8SMatthew G Knepley #define __FUNCT__ "DMCreateSubDM"
1267435a35e8SMatthew G Knepley /*@C
1268435a35e8SMatthew G Knepley   DMCreateSubDM - Returns an IS and DM encapsulating a subproblem defined by the fields passed in.
1269435a35e8SMatthew G Knepley                   The fields are defined by DMCreateFieldIS().
1270435a35e8SMatthew G Knepley 
1271435a35e8SMatthew G Knepley   Not collective
1272435a35e8SMatthew G Knepley 
1273435a35e8SMatthew G Knepley   Input Parameters:
1274435a35e8SMatthew G Knepley + dm - the DM object
1275435a35e8SMatthew G Knepley . numFields - number of fields in this subproblem
12760298fd71SBarry Smith - len       - The number of subproblems in the decomposition (or NULL if not requested)
1277435a35e8SMatthew G Knepley 
1278435a35e8SMatthew G Knepley   Output Parameters:
1279435a35e8SMatthew G Knepley . is - The global indices for the subproblem
1280435a35e8SMatthew G Knepley - dm - The DM for the subproblem
1281435a35e8SMatthew G Knepley 
1282435a35e8SMatthew G Knepley   Level: intermediate
1283435a35e8SMatthew G Knepley 
1284435a35e8SMatthew G Knepley .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1285435a35e8SMatthew G Knepley @*/
1286435a35e8SMatthew G Knepley PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1287435a35e8SMatthew G Knepley {
1288435a35e8SMatthew G Knepley   PetscErrorCode ierr;
1289435a35e8SMatthew G Knepley 
1290435a35e8SMatthew G Knepley   PetscFunctionBegin;
1291435a35e8SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1292435a35e8SMatthew G Knepley   PetscValidPointer(fields,3);
12938865f1eaSKarl Rupp   if (is) PetscValidPointer(is,4);
12948865f1eaSKarl Rupp   if (subdm) PetscValidPointer(subdm,5);
1295435a35e8SMatthew G Knepley   if (dm->ops->createsubdm) {
1296435a35e8SMatthew G Knepley     ierr = (*dm->ops->createsubdm)(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
129782f516ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateSubDM implementation defined");
1298435a35e8SMatthew G Knepley   PetscFunctionReturn(0);
1299435a35e8SMatthew G Knepley }
1300435a35e8SMatthew G Knepley 
130116621825SDmitry Karpeev 
130216621825SDmitry Karpeev #undef __FUNCT__
130316621825SDmitry Karpeev #define __FUNCT__ "DMCreateDomainDecomposition"
130416621825SDmitry Karpeev /*@C
13058d4ac253SDmitry Karpeev   DMCreateDomainDecomposition - Returns lists of IS objects defining a decomposition of a problem into subproblems
13068d4ac253SDmitry Karpeev                           corresponding to restrictions to pairs nested subdomains: each IS contains the global
13078d4ac253SDmitry Karpeev                           indices of the dofs of the corresponding subdomains.  The inner subdomains conceptually
13088d4ac253SDmitry Karpeev                           define a nonoverlapping covering, while outer subdomains can overlap.
13098d4ac253SDmitry Karpeev                           The optional list of DMs define the DM for each subproblem.
131016621825SDmitry Karpeev 
131116621825SDmitry Karpeev   Not collective
131216621825SDmitry Karpeev 
131316621825SDmitry Karpeev   Input Parameter:
131416621825SDmitry Karpeev . dm - the DM object
131516621825SDmitry Karpeev 
131616621825SDmitry Karpeev   Output Parameters:
13170298fd71SBarry Smith + len         - The number of subproblems in the domain decomposition (or NULL if not requested)
13180298fd71SBarry Smith . namelist    - The name for each subdomain (or NULL if not requested)
13190298fd71SBarry Smith . innerislist - The global indices for each inner subdomain (or NULL, if not requested)
13200298fd71SBarry Smith . outerislist - The global indices for each outer subdomain (or NULL, if not requested)
13210298fd71SBarry Smith - dmlist      - The DMs for each subdomain subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined)
132216621825SDmitry Karpeev 
132316621825SDmitry Karpeev   Level: intermediate
132416621825SDmitry Karpeev 
132516621825SDmitry Karpeev   Notes:
132616621825SDmitry Karpeev   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
132716621825SDmitry Karpeev   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
132816621825SDmitry Karpeev   and all of the arrays should be freed with PetscFree().
132916621825SDmitry Karpeev 
13308d4ac253SDmitry Karpeev .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateDomainDecompositionDM(), DMCreateFieldDecomposition()
133116621825SDmitry Karpeev @*/
13328d4ac253SDmitry Karpeev PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist)
133316621825SDmitry Karpeev {
133416621825SDmitry Karpeev   PetscErrorCode      ierr;
1335be081cd6SPeter Brune   DMSubDomainHookLink link;
1336be081cd6SPeter Brune   PetscInt            i,l;
133716621825SDmitry Karpeev 
133816621825SDmitry Karpeev   PetscFunctionBegin;
133916621825SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
134014a18fd3SPeter Brune   if (len)           {PetscValidPointer(len,2);            *len         = 0;}
13410298fd71SBarry Smith   if (namelist)      {PetscValidPointer(namelist,3);       *namelist    = NULL;}
13420298fd71SBarry Smith   if (innerislist)   {PetscValidPointer(innerislist,4);    *innerislist = NULL;}
13430298fd71SBarry Smith   if (outerislist)   {PetscValidPointer(outerislist,5);    *outerislist = NULL;}
13440298fd71SBarry Smith   if (dmlist)        {PetscValidPointer(dmlist,6);         *dmlist      = NULL;}
1345f3f0edfdSDmitry Karpeev   /*
1346f3f0edfdSDmitry Karpeev    Is it a good idea to apply the following check across all impls?
1347f3f0edfdSDmitry Karpeev    Perhaps some impls can have a well-defined decomposition before DMSetUp?
1348f3f0edfdSDmitry Karpeev    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
1349f3f0edfdSDmitry Karpeev    */
1350ce94432eSBarry Smith   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
135116621825SDmitry Karpeev   if (dm->ops->createdomaindecomposition) {
1352be081cd6SPeter Brune     ierr = (*dm->ops->createdomaindecomposition)(dm,&l,namelist,innerislist,outerislist,dmlist);CHKERRQ(ierr);
135314a18fd3SPeter Brune     /* copy subdomain hooks and context over to the subdomain DMs */
135414a18fd3SPeter Brune     if (dmlist) {
1355be081cd6SPeter Brune       for (i = 0; i < l; i++) {
1356be081cd6SPeter Brune         for (link=dm->subdomainhook; link; link=link->next) {
1357be081cd6SPeter Brune           if (link->ddhook) {ierr = (*link->ddhook)(dm,(*dmlist)[i],link->ctx);CHKERRQ(ierr);}
1358be081cd6SPeter Brune         }
1359d425c654SPeter Brune         (*dmlist)[i]->ctx = dm->ctx;
1360e7c4fc90SDmitry Karpeev       }
136114a18fd3SPeter Brune     }
136214a18fd3SPeter Brune     if (len) *len = l;
136314a18fd3SPeter Brune   }
1364e30e807fSPeter Brune   PetscFunctionReturn(0);
1365e30e807fSPeter Brune }
1366e30e807fSPeter Brune 
1367e30e807fSPeter Brune 
1368e30e807fSPeter Brune #undef __FUNCT__
1369e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters"
1370e30e807fSPeter Brune /*@C
1371e30e807fSPeter Brune   DMCreateDomainDecompositionScatters - Returns scatters to the subdomain vectors from the global vector
1372e30e807fSPeter Brune 
1373e30e807fSPeter Brune   Not collective
1374e30e807fSPeter Brune 
1375e30e807fSPeter Brune   Input Parameters:
1376e30e807fSPeter Brune + dm - the DM object
1377e30e807fSPeter Brune . n  - the number of subdomain scatters
1378e30e807fSPeter Brune - subdms - the local subdomains
1379e30e807fSPeter Brune 
1380e30e807fSPeter Brune   Output Parameters:
1381e30e807fSPeter Brune + n     - the number of scatters returned
1382e30e807fSPeter Brune . iscat - scatter from global vector to nonoverlapping global vector entries on subdomain
1383e30e807fSPeter Brune . oscat - scatter from global vector to overlapping global vector entries on subdomain
1384e30e807fSPeter Brune - gscat - scatter from global vector to local vector on subdomain (fills in ghosts)
1385e30e807fSPeter Brune 
1386e30e807fSPeter Brune   Notes: This is an alternative to the iis and ois arguments in DMCreateDomainDecomposition that allow for the solution
1387e30e807fSPeter Brune   of general nonlinear problems with overlapping subdomain methods.  While merely having index sets that enable subsets
1388e30e807fSPeter Brune   of the residual equations to be created is fine for linear problems, nonlinear problems require local assembly of
1389e30e807fSPeter Brune   solution and residual data.
1390e30e807fSPeter Brune 
1391e30e807fSPeter Brune   Level: developer
1392e30e807fSPeter Brune 
1393e30e807fSPeter Brune .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1394e30e807fSPeter Brune @*/
1395e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecompositionScatters(DM dm,PetscInt n,DM *subdms,VecScatter **iscat,VecScatter **oscat,VecScatter **gscat)
1396e30e807fSPeter Brune {
1397e30e807fSPeter Brune   PetscErrorCode ierr;
1398e30e807fSPeter Brune 
1399e30e807fSPeter Brune   PetscFunctionBegin;
1400e30e807fSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1401e30e807fSPeter Brune   PetscValidPointer(subdms,3);
1402e30e807fSPeter Brune   if (dm->ops->createddscatters) {
1403e30e807fSPeter Brune     ierr = (*dm->ops->createddscatters)(dm,n,subdms,iscat,oscat,gscat);CHKERRQ(ierr);
140482f516ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateDomainDecompositionLocalScatter implementation defined");
1405e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1406e7c4fc90SDmitry Karpeev }
1407e7c4fc90SDmitry Karpeev 
1408731c8d9eSDmitry Karpeev #undef __FUNCT__
140947c6ae99SBarry Smith #define __FUNCT__ "DMRefine"
141047c6ae99SBarry Smith /*@
141147c6ae99SBarry Smith   DMRefine - Refines a DM object
141247c6ae99SBarry Smith 
141347c6ae99SBarry Smith   Collective on DM
141447c6ae99SBarry Smith 
141547c6ae99SBarry Smith   Input Parameter:
141647c6ae99SBarry Smith + dm   - the DM object
141791d95f02SJed Brown - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
141847c6ae99SBarry Smith 
141947c6ae99SBarry Smith   Output Parameter:
14200298fd71SBarry Smith . dmf - the refined DM, or NULL
1421ae0a1c52SMatthew G Knepley 
14220298fd71SBarry Smith   Note: If no refinement was done, the return value is NULL
142347c6ae99SBarry Smith 
142447c6ae99SBarry Smith   Level: developer
142547c6ae99SBarry Smith 
1426e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
142747c6ae99SBarry Smith @*/
14287087cfbeSBarry Smith PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
142947c6ae99SBarry Smith {
143047c6ae99SBarry Smith   PetscErrorCode   ierr;
1431c833c3b5SJed Brown   DMRefineHookLink link;
143247c6ae99SBarry Smith 
143347c6ae99SBarry Smith   PetscFunctionBegin;
1434732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
143547c6ae99SBarry Smith   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
14364057135bSMatthew G Knepley   if (*dmf) {
143743842a1eSJed Brown     (*dmf)->ops->creatematrix = dm->ops->creatematrix;
14388865f1eaSKarl Rupp 
14398cd211a4SJed Brown     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
14408865f1eaSKarl Rupp 
1441644e2e5bSBarry Smith     (*dmf)->ctx       = dm->ctx;
14420598a293SJed Brown     (*dmf)->leveldown = dm->leveldown;
1443656b349aSBarry Smith     (*dmf)->levelup   = dm->levelup + 1;
14448865f1eaSKarl Rupp 
1445e4b4b23bSJed Brown     ierr = DMSetMatType(*dmf,dm->mattype);CHKERRQ(ierr);
1446c833c3b5SJed Brown     for (link=dm->refinehook; link; link=link->next) {
14478865f1eaSKarl Rupp       if (link->refinehook) {
14488865f1eaSKarl Rupp         ierr = (*link->refinehook)(dm,*dmf,link->ctx);CHKERRQ(ierr);
14498865f1eaSKarl Rupp       }
1450c833c3b5SJed Brown     }
1451c833c3b5SJed Brown   }
1452c833c3b5SJed Brown   PetscFunctionReturn(0);
1453c833c3b5SJed Brown }
1454c833c3b5SJed Brown 
1455c833c3b5SJed Brown #undef __FUNCT__
1456c833c3b5SJed Brown #define __FUNCT__ "DMRefineHookAdd"
1457bb9467b5SJed Brown /*@C
1458c833c3b5SJed Brown    DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid
1459c833c3b5SJed Brown 
1460c833c3b5SJed Brown    Logically Collective
1461c833c3b5SJed Brown 
1462c833c3b5SJed Brown    Input Arguments:
1463c833c3b5SJed Brown +  coarse - nonlinear solver context on which to run a hook when restricting to a coarser level
1464c833c3b5SJed Brown .  refinehook - function to run when setting up a coarser level
1465c833c3b5SJed Brown .  interphook - function to run to update data on finer levels (once per SNESSolve())
14660298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1467c833c3b5SJed Brown 
1468c833c3b5SJed Brown    Calling sequence of refinehook:
1469c833c3b5SJed Brown $    refinehook(DM coarse,DM fine,void *ctx);
1470c833c3b5SJed Brown 
1471c833c3b5SJed Brown +  coarse - coarse level DM
1472c833c3b5SJed Brown .  fine - fine level DM to interpolate problem to
1473c833c3b5SJed Brown -  ctx - optional user-defined function context
1474c833c3b5SJed Brown 
1475c833c3b5SJed Brown    Calling sequence for interphook:
1476c833c3b5SJed Brown $    interphook(DM coarse,Mat interp,DM fine,void *ctx)
1477c833c3b5SJed Brown 
1478c833c3b5SJed Brown +  coarse - coarse level DM
1479c833c3b5SJed Brown .  interp - matrix interpolating a coarse-level solution to the finer grid
1480c833c3b5SJed Brown .  fine - fine level DM to update
1481c833c3b5SJed Brown -  ctx - optional user-defined function context
1482c833c3b5SJed Brown 
1483c833c3b5SJed Brown    Level: advanced
1484c833c3b5SJed Brown 
1485c833c3b5SJed Brown    Notes:
1486c833c3b5SJed Brown    This function is only needed if auxiliary data needs to be passed to fine grids while grid sequencing
1487c833c3b5SJed Brown 
1488c833c3b5SJed Brown    If this function is called multiple times, the hooks will be run in the order they are added.
1489c833c3b5SJed Brown 
1490bb9467b5SJed Brown    This function is currently not available from Fortran.
1491bb9467b5SJed Brown 
1492c833c3b5SJed Brown .seealso: DMCoarsenHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1493c833c3b5SJed Brown @*/
1494c833c3b5SJed Brown PetscErrorCode DMRefineHookAdd(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx)
1495c833c3b5SJed Brown {
1496c833c3b5SJed Brown   PetscErrorCode   ierr;
1497c833c3b5SJed Brown   DMRefineHookLink link,*p;
1498c833c3b5SJed Brown 
1499c833c3b5SJed Brown   PetscFunctionBegin;
1500c833c3b5SJed Brown   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
1501c833c3b5SJed Brown   for (p=&coarse->refinehook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1502c833c3b5SJed Brown   ierr             = PetscMalloc(sizeof(struct _DMRefineHookLink),&link);CHKERRQ(ierr);
1503c833c3b5SJed Brown   link->refinehook = refinehook;
1504c833c3b5SJed Brown   link->interphook = interphook;
1505c833c3b5SJed Brown   link->ctx        = ctx;
15060298fd71SBarry Smith   link->next       = NULL;
1507c833c3b5SJed Brown   *p               = link;
1508c833c3b5SJed Brown   PetscFunctionReturn(0);
1509c833c3b5SJed Brown }
1510c833c3b5SJed Brown 
1511c833c3b5SJed Brown #undef __FUNCT__
1512c833c3b5SJed Brown #define __FUNCT__ "DMInterpolate"
1513c833c3b5SJed Brown /*@
1514c833c3b5SJed Brown    DMInterpolate - interpolates user-defined problem data to a finer DM by running hooks registered by DMRefineHookAdd()
1515c833c3b5SJed Brown 
1516c833c3b5SJed Brown    Collective if any hooks are
1517c833c3b5SJed Brown 
1518c833c3b5SJed Brown    Input Arguments:
1519c833c3b5SJed Brown +  coarse - coarser DM to use as a base
1520c833c3b5SJed Brown .  restrct - interpolation matrix, apply using MatInterpolate()
1521c833c3b5SJed Brown -  fine - finer DM to update
1522c833c3b5SJed Brown 
1523c833c3b5SJed Brown    Level: developer
1524c833c3b5SJed Brown 
1525c833c3b5SJed Brown .seealso: DMRefineHookAdd(), MatInterpolate()
1526c833c3b5SJed Brown @*/
1527c833c3b5SJed Brown PetscErrorCode DMInterpolate(DM coarse,Mat interp,DM fine)
1528c833c3b5SJed Brown {
1529c833c3b5SJed Brown   PetscErrorCode   ierr;
1530c833c3b5SJed Brown   DMRefineHookLink link;
1531c833c3b5SJed Brown 
1532c833c3b5SJed Brown   PetscFunctionBegin;
1533c833c3b5SJed Brown   for (link=fine->refinehook; link; link=link->next) {
15348865f1eaSKarl Rupp     if (link->interphook) {
15358865f1eaSKarl Rupp       ierr = (*link->interphook)(coarse,interp,fine,link->ctx);CHKERRQ(ierr);
15368865f1eaSKarl Rupp     }
15374057135bSMatthew G Knepley   }
153847c6ae99SBarry Smith   PetscFunctionReturn(0);
153947c6ae99SBarry Smith }
154047c6ae99SBarry Smith 
154147c6ae99SBarry Smith #undef __FUNCT__
1542eb3f98d2SBarry Smith #define __FUNCT__ "DMGetRefineLevel"
1543eb3f98d2SBarry Smith /*@
1544eb3f98d2SBarry Smith     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
1545eb3f98d2SBarry Smith 
1546eb3f98d2SBarry Smith     Not Collective
1547eb3f98d2SBarry Smith 
1548eb3f98d2SBarry Smith     Input Parameter:
1549eb3f98d2SBarry Smith .   dm - the DM object
1550eb3f98d2SBarry Smith 
1551eb3f98d2SBarry Smith     Output Parameter:
1552eb3f98d2SBarry Smith .   level - number of refinements
1553eb3f98d2SBarry Smith 
1554eb3f98d2SBarry Smith     Level: developer
1555eb3f98d2SBarry Smith 
15566a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1557eb3f98d2SBarry Smith 
1558eb3f98d2SBarry Smith @*/
1559eb3f98d2SBarry Smith PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
1560eb3f98d2SBarry Smith {
1561eb3f98d2SBarry Smith   PetscFunctionBegin;
1562eb3f98d2SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1563eb3f98d2SBarry Smith   *level = dm->levelup;
1564eb3f98d2SBarry Smith   PetscFunctionReturn(0);
1565eb3f98d2SBarry Smith }
1566eb3f98d2SBarry Smith 
1567eb3f98d2SBarry Smith #undef __FUNCT__
1568baf369e7SPeter Brune #define __FUNCT__ "DMGlobalToLocalHookAdd"
1569bb9467b5SJed Brown /*@C
1570baf369e7SPeter Brune    DMGlobalToLocalHookAdd - adds a callback to be run when global to local is called
1571baf369e7SPeter Brune 
1572baf369e7SPeter Brune    Logically Collective
1573baf369e7SPeter Brune 
1574baf369e7SPeter Brune    Input Arguments:
1575baf369e7SPeter Brune +  dm - the DM
1576baf369e7SPeter Brune .  beginhook - function to run at the beginning of DMGlobalToLocalBegin()
1577baf369e7SPeter Brune .  endhook - function to run after DMGlobalToLocalEnd() has completed
15780298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1579baf369e7SPeter Brune 
1580baf369e7SPeter Brune    Calling sequence for beginhook:
1581baf369e7SPeter Brune $    beginhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx)
1582baf369e7SPeter Brune 
1583baf369e7SPeter Brune +  dm - global DM
1584baf369e7SPeter Brune .  g - global vector
1585baf369e7SPeter Brune .  mode - mode
1586baf369e7SPeter Brune .  l - local vector
1587baf369e7SPeter Brune -  ctx - optional user-defined function context
1588baf369e7SPeter Brune 
1589baf369e7SPeter Brune 
1590baf369e7SPeter Brune    Calling sequence for endhook:
1591ec4806b8SPeter Brune $    endhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx)
1592baf369e7SPeter Brune 
1593baf369e7SPeter Brune +  global - global DM
1594baf369e7SPeter Brune -  ctx - optional user-defined function context
1595baf369e7SPeter Brune 
1596baf369e7SPeter Brune    Level: advanced
1597baf369e7SPeter Brune 
1598baf369e7SPeter Brune .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1599baf369e7SPeter Brune @*/
1600baf369e7SPeter Brune PetscErrorCode DMGlobalToLocalHookAdd(DM dm,PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*),void *ctx)
1601baf369e7SPeter Brune {
1602baf369e7SPeter Brune   PetscErrorCode          ierr;
1603baf369e7SPeter Brune   DMGlobalToLocalHookLink link,*p;
1604baf369e7SPeter Brune 
1605baf369e7SPeter Brune   PetscFunctionBegin;
1606baf369e7SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1607baf369e7SPeter Brune   for (p=&dm->gtolhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1608baf369e7SPeter Brune   ierr            = PetscMalloc(sizeof(struct _DMGlobalToLocalHookLink),&link);CHKERRQ(ierr);
1609baf369e7SPeter Brune   link->beginhook = beginhook;
1610baf369e7SPeter Brune   link->endhook   = endhook;
1611baf369e7SPeter Brune   link->ctx       = ctx;
16120298fd71SBarry Smith   link->next      = NULL;
1613baf369e7SPeter Brune   *p              = link;
1614baf369e7SPeter Brune   PetscFunctionReturn(0);
1615baf369e7SPeter Brune }
1616baf369e7SPeter Brune 
1617baf369e7SPeter Brune #undef __FUNCT__
161847c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin"
161947c6ae99SBarry Smith /*@
162047c6ae99SBarry Smith     DMGlobalToLocalBegin - Begins updating local vectors from global vector
162147c6ae99SBarry Smith 
162247c6ae99SBarry Smith     Neighbor-wise Collective on DM
162347c6ae99SBarry Smith 
162447c6ae99SBarry Smith     Input Parameters:
162547c6ae99SBarry Smith +   dm - the DM object
162647c6ae99SBarry Smith .   g - the global vector
162747c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
162847c6ae99SBarry Smith -   l - the local vector
162947c6ae99SBarry Smith 
163047c6ae99SBarry Smith 
163147c6ae99SBarry Smith     Level: beginner
163247c6ae99SBarry Smith 
1633e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
163447c6ae99SBarry Smith 
163547c6ae99SBarry Smith @*/
16367087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
163747c6ae99SBarry Smith {
16387128ae9fSMatthew G Knepley   PetscSF                 sf;
163947c6ae99SBarry Smith   PetscErrorCode          ierr;
1640baf369e7SPeter Brune   DMGlobalToLocalHookLink link;
164147c6ae99SBarry Smith 
164247c6ae99SBarry Smith   PetscFunctionBegin;
1643171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1644baf369e7SPeter Brune   for (link=dm->gtolhook; link; link=link->next) {
16458865f1eaSKarl Rupp     if (link->beginhook) {
16468865f1eaSKarl Rupp       ierr = (*link->beginhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);
16478865f1eaSKarl Rupp     }
1648baf369e7SPeter Brune   }
16497128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
16507128ae9fSMatthew G Knepley   if (sf) {
16517128ae9fSMatthew G Knepley     PetscScalar *lArray, *gArray;
16527128ae9fSMatthew G Knepley 
165382f516ccSBarry Smith     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
16547128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
16557128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
16567128ae9fSMatthew G Knepley     ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
16577128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
16587128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
16597128ae9fSMatthew G Knepley   } else {
1660843c4018SMatthew G Knepley     ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
16617128ae9fSMatthew G Knepley   }
166247c6ae99SBarry Smith   PetscFunctionReturn(0);
166347c6ae99SBarry Smith }
166447c6ae99SBarry Smith 
166547c6ae99SBarry Smith #undef __FUNCT__
166647c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd"
166747c6ae99SBarry Smith /*@
166847c6ae99SBarry Smith     DMGlobalToLocalEnd - Ends updating local vectors from global vector
166947c6ae99SBarry Smith 
167047c6ae99SBarry Smith     Neighbor-wise Collective on DM
167147c6ae99SBarry Smith 
167247c6ae99SBarry Smith     Input Parameters:
167347c6ae99SBarry Smith +   dm - the DM object
167447c6ae99SBarry Smith .   g - the global vector
167547c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
167647c6ae99SBarry Smith -   l - the local vector
167747c6ae99SBarry Smith 
167847c6ae99SBarry Smith 
167947c6ae99SBarry Smith     Level: beginner
168047c6ae99SBarry Smith 
1681e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
168247c6ae99SBarry Smith 
168347c6ae99SBarry Smith @*/
16847087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
168547c6ae99SBarry Smith {
16867128ae9fSMatthew G Knepley   PetscSF                 sf;
168747c6ae99SBarry Smith   PetscErrorCode          ierr;
168861a3c1faSSatish Balay   PetscScalar             *lArray, *gArray;
1689baf369e7SPeter Brune   DMGlobalToLocalHookLink link;
169047c6ae99SBarry Smith 
169147c6ae99SBarry Smith   PetscFunctionBegin;
1692171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
16937128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
16947128ae9fSMatthew G Knepley   if (sf) {
169582f516ccSBarry Smith     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
16967128ae9fSMatthew G Knepley 
16977128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
16987128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
16997128ae9fSMatthew G Knepley     ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
17007128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
17017128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
17027128ae9fSMatthew G Knepley   } else {
1703843c4018SMatthew G Knepley     ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
17047128ae9fSMatthew G Knepley   }
1705baf369e7SPeter Brune   for (link=dm->gtolhook; link; link=link->next) {
1706baf369e7SPeter Brune     if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);}
1707baf369e7SPeter Brune   }
170847c6ae99SBarry Smith   PetscFunctionReturn(0);
170947c6ae99SBarry Smith }
171047c6ae99SBarry Smith 
171147c6ae99SBarry Smith #undef __FUNCT__
1712d4d07f1eSToby Isaac #define __FUNCT__ "DMLocalToGlobalHookAdd"
1713d4d07f1eSToby Isaac /*@C
1714d4d07f1eSToby Isaac    DMLocalToGlobalHookAdd - adds a callback to be run when a local to global is called
1715d4d07f1eSToby Isaac 
1716d4d07f1eSToby Isaac    Logically Collective
1717d4d07f1eSToby Isaac 
1718d4d07f1eSToby Isaac    Input Arguments:
1719d4d07f1eSToby Isaac +  dm - the DM
1720d4d07f1eSToby Isaac .  beginhook - function to run at the beginning of DMLocalToGlobalBegin()
1721d4d07f1eSToby Isaac .  endhook - function to run after DMLocalToGlobalEnd() has completed
1722d4d07f1eSToby Isaac -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1723d4d07f1eSToby Isaac 
1724d4d07f1eSToby Isaac    Calling sequence for beginhook:
1725d4d07f1eSToby Isaac $    beginhook(DM fine,Vec l,InsertMode mode,Vec g,void *ctx)
1726d4d07f1eSToby Isaac 
1727d4d07f1eSToby Isaac +  dm - global DM
1728d4d07f1eSToby Isaac .  l - local vector
1729d4d07f1eSToby Isaac .  mode - mode
1730d4d07f1eSToby Isaac .  g - global vector
1731d4d07f1eSToby Isaac -  ctx - optional user-defined function context
1732d4d07f1eSToby Isaac 
1733d4d07f1eSToby Isaac 
1734d4d07f1eSToby Isaac    Calling sequence for endhook:
1735d4d07f1eSToby Isaac $    endhook(DM fine,Vec l,InsertMode mode,Vec g,void *ctx)
1736d4d07f1eSToby Isaac 
1737d4d07f1eSToby Isaac +  global - global DM
1738d4d07f1eSToby Isaac .  l - local vector
1739d4d07f1eSToby Isaac .  mode - mode
1740d4d07f1eSToby Isaac .  g - global vector
1741d4d07f1eSToby Isaac -  ctx - optional user-defined function context
1742d4d07f1eSToby Isaac 
1743d4d07f1eSToby Isaac    Level: advanced
1744d4d07f1eSToby Isaac 
1745d4d07f1eSToby Isaac .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1746d4d07f1eSToby Isaac @*/
1747d4d07f1eSToby Isaac PetscErrorCode DMLocalToGlobalHookAdd(DM dm,PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*),void *ctx)
1748d4d07f1eSToby Isaac {
1749d4d07f1eSToby Isaac   PetscErrorCode          ierr;
1750d4d07f1eSToby Isaac   DMLocalToGlobalHookLink link,*p;
1751d4d07f1eSToby Isaac 
1752d4d07f1eSToby Isaac   PetscFunctionBegin;
1753d4d07f1eSToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1754d4d07f1eSToby Isaac   for (p=&dm->ltoghook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1755d4d07f1eSToby Isaac   ierr            = PetscMalloc(sizeof(struct _DMLocalToGlobalHookLink),&link);CHKERRQ(ierr);
1756d4d07f1eSToby Isaac   link->beginhook = beginhook;
1757d4d07f1eSToby Isaac   link->endhook   = endhook;
1758d4d07f1eSToby Isaac   link->ctx       = ctx;
1759d4d07f1eSToby Isaac   link->next      = NULL;
1760d4d07f1eSToby Isaac   *p              = link;
1761d4d07f1eSToby Isaac   PetscFunctionReturn(0);
1762d4d07f1eSToby Isaac }
1763d4d07f1eSToby Isaac 
1764d4d07f1eSToby Isaac #undef __FUNCT__
17659a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalBegin"
176647c6ae99SBarry Smith /*@
17679a42bb27SBarry Smith     DMLocalToGlobalBegin - updates global vectors from local vectors
17689a42bb27SBarry Smith 
17699a42bb27SBarry Smith     Neighbor-wise Collective on DM
17709a42bb27SBarry Smith 
17719a42bb27SBarry Smith     Input Parameters:
17729a42bb27SBarry Smith +   dm - the DM object
1773f6813fd5SJed Brown .   l - the local vector
17749a42bb27SBarry Smith .   mode - if INSERT_VALUES then no parallel communication is used, if ADD_VALUES then all ghost points from the same base point accumulate into that
17759a42bb27SBarry Smith            base point.
1776f6813fd5SJed Brown - - the global vector
17779a42bb27SBarry Smith 
17789a42bb27SBarry Smith     Notes: In the ADD_VALUES case you normally would zero the receiving vector before beginning this operation. If you would like to simply add the non-ghosted values in the local
17799a42bb27SBarry Smith            array into the global array you need to either (1) zero the ghosted locations and use ADD_VALUES or (2) use INSERT_VALUES into a work global array and then add the work
17809a42bb27SBarry Smith            global array to the final global array with VecAXPY().
17819a42bb27SBarry Smith 
17829a42bb27SBarry Smith     Level: beginner
17839a42bb27SBarry Smith 
1784e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
17859a42bb27SBarry Smith 
17869a42bb27SBarry Smith @*/
17877087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
17889a42bb27SBarry Smith {
17897128ae9fSMatthew G Knepley   PetscSF                 sf;
17909a42bb27SBarry Smith   PetscErrorCode          ierr;
1791d4d07f1eSToby Isaac   DMLocalToGlobalHookLink link;
17929a42bb27SBarry Smith 
17939a42bb27SBarry Smith   PetscFunctionBegin;
1794171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1795d4d07f1eSToby Isaac   for (link=dm->ltoghook; link; link=link->next) {
1796d4d07f1eSToby Isaac     if (link->beginhook) {
1797d4d07f1eSToby Isaac       ierr = (*link->beginhook)(dm,l,mode,g,link->ctx);CHKERRQ(ierr);
1798d4d07f1eSToby Isaac     }
1799d4d07f1eSToby Isaac   }
18007128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
18017128ae9fSMatthew G Knepley   if (sf) {
18027128ae9fSMatthew G Knepley     MPI_Op      op;
18037128ae9fSMatthew G Knepley     PetscScalar *lArray, *gArray;
18047128ae9fSMatthew G Knepley 
18057128ae9fSMatthew G Knepley     switch (mode) {
18067128ae9fSMatthew G Knepley     case INSERT_VALUES:
18077128ae9fSMatthew G Knepley     case INSERT_ALL_VALUES:
18088bfbc91cSJed Brown       op = MPIU_REPLACE; break;
18097128ae9fSMatthew G Knepley     case ADD_VALUES:
18107128ae9fSMatthew G Knepley     case ADD_ALL_VALUES:
18117128ae9fSMatthew G Knepley       op = MPI_SUM; break;
18127128ae9fSMatthew G Knepley     default:
181382f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
18147128ae9fSMatthew G Knepley     }
18157128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
18167128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
18177128ae9fSMatthew G Knepley     ierr = PetscSFReduceBegin(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr);
18187128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
18197128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
18207128ae9fSMatthew G Knepley   } else {
1821843c4018SMatthew G Knepley     ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
18227128ae9fSMatthew G Knepley   }
18239a42bb27SBarry Smith   PetscFunctionReturn(0);
18249a42bb27SBarry Smith }
18259a42bb27SBarry Smith 
18269a42bb27SBarry Smith #undef __FUNCT__
18279a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalEnd"
18289a42bb27SBarry Smith /*@
18299a42bb27SBarry Smith     DMLocalToGlobalEnd - updates global vectors from local vectors
183047c6ae99SBarry Smith 
183147c6ae99SBarry Smith     Neighbor-wise Collective on DM
183247c6ae99SBarry Smith 
183347c6ae99SBarry Smith     Input Parameters:
183447c6ae99SBarry Smith +   dm - the DM object
1835f6813fd5SJed Brown .   l - the local vector
183647c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
1837f6813fd5SJed Brown -   g - the global vector
183847c6ae99SBarry Smith 
183947c6ae99SBarry Smith 
184047c6ae99SBarry Smith     Level: beginner
184147c6ae99SBarry Smith 
1842e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
184347c6ae99SBarry Smith 
184447c6ae99SBarry Smith @*/
18457087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
184647c6ae99SBarry Smith {
18477128ae9fSMatthew G Knepley   PetscSF        sf;
184847c6ae99SBarry Smith   PetscErrorCode ierr;
1849d4d07f1eSToby Isaac   DMLocalToGlobalHookLink link;
185047c6ae99SBarry Smith 
185147c6ae99SBarry Smith   PetscFunctionBegin;
1852171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
18537128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
18547128ae9fSMatthew G Knepley   if (sf) {
18557128ae9fSMatthew G Knepley     MPI_Op      op;
18567128ae9fSMatthew G Knepley     PetscScalar *lArray, *gArray;
18577128ae9fSMatthew G Knepley 
18587128ae9fSMatthew G Knepley     switch (mode) {
18597128ae9fSMatthew G Knepley     case INSERT_VALUES:
18607128ae9fSMatthew G Knepley     case INSERT_ALL_VALUES:
18618bfbc91cSJed Brown       op = MPIU_REPLACE; break;
18627128ae9fSMatthew G Knepley     case ADD_VALUES:
18637128ae9fSMatthew G Knepley     case ADD_ALL_VALUES:
18647128ae9fSMatthew G Knepley       op = MPI_SUM; break;
18657128ae9fSMatthew G Knepley     default:
186682f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
18677128ae9fSMatthew G Knepley     }
18687128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
18697128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
18707128ae9fSMatthew G Knepley     ierr = PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr);
18717128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
18727128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
18737128ae9fSMatthew G Knepley   } else {
1874843c4018SMatthew G Knepley     ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
18757128ae9fSMatthew G Knepley   }
1876d4d07f1eSToby Isaac   for (link=dm->ltoghook; link; link=link->next) {
1877d4d07f1eSToby Isaac     if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);}
1878d4d07f1eSToby Isaac   }
187947c6ae99SBarry Smith   PetscFunctionReturn(0);
188047c6ae99SBarry Smith }
188147c6ae99SBarry Smith 
188247c6ae99SBarry Smith #undef __FUNCT__
1883f089877aSRichard Tran Mills #define __FUNCT__ "DMLocalToLocalBegin"
1884f089877aSRichard Tran Mills /*@
1885bc0a1609SRichard Tran Mills    DMLocalToLocalBegin - Maps from a local vector (including ghost points
1886bc0a1609SRichard Tran Mills    that contain irrelevant values) to another local vector where the ghost
1887d78e899eSRichard Tran Mills    points in the second are set correctly. Must be followed by DMLocalToLocalEnd().
1888f089877aSRichard Tran Mills 
1889bc0a1609SRichard Tran Mills    Neighbor-wise Collective on DM and Vec
1890f089877aSRichard Tran Mills 
1891f089877aSRichard Tran Mills    Input Parameters:
1892f089877aSRichard Tran Mills +  dm - the DM object
1893bc0a1609SRichard Tran Mills .  g - the original local vector
1894bc0a1609SRichard Tran Mills -  mode - one of INSERT_VALUES or ADD_VALUES
1895f089877aSRichard Tran Mills 
1896bc0a1609SRichard Tran Mills    Output Parameter:
1897bc0a1609SRichard Tran Mills .  l  - the local vector with correct ghost values
1898f089877aSRichard Tran Mills 
1899f089877aSRichard Tran Mills    Level: intermediate
1900f089877aSRichard Tran Mills 
1901bc0a1609SRichard Tran Mills    Notes:
1902bc0a1609SRichard Tran Mills    The local vectors used here need not be the same as those
1903bc0a1609SRichard Tran Mills    obtained from DMCreateLocalVector(), BUT they
1904bc0a1609SRichard Tran Mills    must have the same parallel data layout; they could, for example, be
1905bc0a1609SRichard Tran Mills    obtained with VecDuplicate() from the DM originating vectors.
1906bc0a1609SRichard Tran Mills 
1907bc0a1609SRichard Tran Mills .keywords: DM, local-to-local, begin
1908bc0a1609SRichard Tran Mills .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalEnd(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1909f089877aSRichard Tran Mills 
1910f089877aSRichard Tran Mills @*/
1911f089877aSRichard Tran Mills PetscErrorCode  DMLocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
1912f089877aSRichard Tran Mills {
1913f089877aSRichard Tran Mills   PetscErrorCode          ierr;
1914f089877aSRichard Tran Mills 
1915f089877aSRichard Tran Mills   PetscFunctionBegin;
1916f089877aSRichard Tran Mills   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1917f089877aSRichard Tran Mills   ierr = (*dm->ops->localtolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1918f089877aSRichard Tran Mills   PetscFunctionReturn(0);
1919f089877aSRichard Tran Mills }
1920f089877aSRichard Tran Mills 
1921f089877aSRichard Tran Mills #undef __FUNCT__
1922f089877aSRichard Tran Mills #define __FUNCT__ "DMLocalToLocalEnd"
1923f089877aSRichard Tran Mills /*@
1924bc0a1609SRichard Tran Mills    DMLocalToLocalEnd - Maps from a local vector (including ghost points
1925bc0a1609SRichard Tran Mills    that contain irrelevant values) to another local vector where the ghost
1926d78e899eSRichard Tran Mills    points in the second are set correctly. Must be preceded by DMLocalToLocalBegin().
1927f089877aSRichard Tran Mills 
1928bc0a1609SRichard Tran Mills    Neighbor-wise Collective on DM and Vec
1929f089877aSRichard Tran Mills 
1930f089877aSRichard Tran Mills    Input Parameters:
1931bc0a1609SRichard Tran Mills +  da - the DM object
1932bc0a1609SRichard Tran Mills .  g - the original local vector
1933bc0a1609SRichard Tran Mills -  mode - one of INSERT_VALUES or ADD_VALUES
1934f089877aSRichard Tran Mills 
1935bc0a1609SRichard Tran Mills    Output Parameter:
1936bc0a1609SRichard Tran Mills .  l  - the local vector with correct ghost values
1937f089877aSRichard Tran Mills 
1938f089877aSRichard Tran Mills    Level: intermediate
1939f089877aSRichard Tran Mills 
1940bc0a1609SRichard Tran Mills    Notes:
1941bc0a1609SRichard Tran Mills    The local vectors used here need not be the same as those
1942bc0a1609SRichard Tran Mills    obtained from DMCreateLocalVector(), BUT they
1943bc0a1609SRichard Tran Mills    must have the same parallel data layout; they could, for example, be
1944bc0a1609SRichard Tran Mills    obtained with VecDuplicate() from the DM originating vectors.
1945bc0a1609SRichard Tran Mills 
1946bc0a1609SRichard Tran Mills .keywords: DM, local-to-local, end
1947bc0a1609SRichard Tran Mills .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1948f089877aSRichard Tran Mills 
1949f089877aSRichard Tran Mills @*/
1950f089877aSRichard Tran Mills PetscErrorCode  DMLocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
1951f089877aSRichard Tran Mills {
1952f089877aSRichard Tran Mills   PetscErrorCode          ierr;
1953f089877aSRichard Tran Mills 
1954f089877aSRichard Tran Mills   PetscFunctionBegin;
1955f089877aSRichard Tran Mills   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1956f089877aSRichard Tran Mills   ierr = (*dm->ops->localtolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1957f089877aSRichard Tran Mills   PetscFunctionReturn(0);
1958f089877aSRichard Tran Mills }
1959f089877aSRichard Tran Mills 
1960f089877aSRichard Tran Mills 
1961f089877aSRichard Tran Mills #undef __FUNCT__
196247c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen"
196347c6ae99SBarry Smith /*@
196447c6ae99SBarry Smith     DMCoarsen - Coarsens a DM object
196547c6ae99SBarry Smith 
196647c6ae99SBarry Smith     Collective on DM
196747c6ae99SBarry Smith 
196847c6ae99SBarry Smith     Input Parameter:
196947c6ae99SBarry Smith +   dm - the DM object
197091d95f02SJed Brown -   comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
197147c6ae99SBarry Smith 
197247c6ae99SBarry Smith     Output Parameter:
197347c6ae99SBarry Smith .   dmc - the coarsened DM
197447c6ae99SBarry Smith 
197547c6ae99SBarry Smith     Level: developer
197647c6ae99SBarry Smith 
1977e727c939SJed Brown .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
197847c6ae99SBarry Smith 
197947c6ae99SBarry Smith @*/
19807087cfbeSBarry Smith PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
198147c6ae99SBarry Smith {
198247c6ae99SBarry Smith   PetscErrorCode    ierr;
1983b17ce1afSJed Brown   DMCoarsenHookLink link;
198447c6ae99SBarry Smith 
198547c6ae99SBarry Smith   PetscFunctionBegin;
1986171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
198747c6ae99SBarry Smith   ierr                      = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
198843842a1eSJed Brown   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
19898cd211a4SJed Brown   ierr                      = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1990644e2e5bSBarry Smith   (*dmc)->ctx               = dm->ctx;
19910598a293SJed Brown   (*dmc)->levelup           = dm->levelup;
1992656b349aSBarry Smith   (*dmc)->leveldown         = dm->leveldown + 1;
1993e4b4b23bSJed Brown   ierr                      = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr);
1994b17ce1afSJed Brown   for (link=dm->coarsenhook; link; link=link->next) {
1995b17ce1afSJed Brown     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1996b17ce1afSJed Brown   }
1997b17ce1afSJed Brown   PetscFunctionReturn(0);
1998b17ce1afSJed Brown }
1999b17ce1afSJed Brown 
2000b17ce1afSJed Brown #undef __FUNCT__
2001b17ce1afSJed Brown #define __FUNCT__ "DMCoarsenHookAdd"
2002bb9467b5SJed Brown /*@C
2003b17ce1afSJed Brown    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
2004b17ce1afSJed Brown 
2005b17ce1afSJed Brown    Logically Collective
2006b17ce1afSJed Brown 
2007b17ce1afSJed Brown    Input Arguments:
2008b17ce1afSJed Brown +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
2009b17ce1afSJed Brown .  coarsenhook - function to run when setting up a coarser level
2010b17ce1afSJed Brown .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
20110298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
2012b17ce1afSJed Brown 
2013b17ce1afSJed Brown    Calling sequence of coarsenhook:
2014b17ce1afSJed Brown $    coarsenhook(DM fine,DM coarse,void *ctx);
2015b17ce1afSJed Brown 
2016b17ce1afSJed Brown +  fine - fine level DM
2017b17ce1afSJed Brown .  coarse - coarse level DM to restrict problem to
2018b17ce1afSJed Brown -  ctx - optional user-defined function context
2019b17ce1afSJed Brown 
2020b17ce1afSJed Brown    Calling sequence for restricthook:
2021c833c3b5SJed Brown $    restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx)
2022b17ce1afSJed Brown 
2023b17ce1afSJed Brown +  fine - fine level DM
2024b17ce1afSJed Brown .  mrestrict - matrix restricting a fine-level solution to the coarse grid
2025c833c3b5SJed Brown .  rscale - scaling vector for restriction
2026c833c3b5SJed Brown .  inject - matrix restricting by injection
2027b17ce1afSJed Brown .  coarse - coarse level DM to update
2028b17ce1afSJed Brown -  ctx - optional user-defined function context
2029b17ce1afSJed Brown 
2030b17ce1afSJed Brown    Level: advanced
2031b17ce1afSJed Brown 
2032b17ce1afSJed Brown    Notes:
2033b17ce1afSJed Brown    This function is only needed if auxiliary data needs to be set up on coarse grids.
2034b17ce1afSJed Brown 
2035b17ce1afSJed Brown    If this function is called multiple times, the hooks will be run in the order they are added.
2036b17ce1afSJed Brown 
2037b17ce1afSJed Brown    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
2038b17ce1afSJed Brown    extract the finest level information from its context (instead of from the SNES).
2039b17ce1afSJed Brown 
2040bb9467b5SJed Brown    This function is currently not available from Fortran.
2041bb9467b5SJed Brown 
2042c833c3b5SJed Brown .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
2043b17ce1afSJed Brown @*/
2044b17ce1afSJed Brown PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
2045b17ce1afSJed Brown {
2046b17ce1afSJed Brown   PetscErrorCode    ierr;
2047b17ce1afSJed Brown   DMCoarsenHookLink link,*p;
2048b17ce1afSJed Brown 
2049b17ce1afSJed Brown   PetscFunctionBegin;
2050b17ce1afSJed Brown   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
20516bfea28cSJed Brown   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
2052b17ce1afSJed Brown   ierr               = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
2053b17ce1afSJed Brown   link->coarsenhook  = coarsenhook;
2054b17ce1afSJed Brown   link->restricthook = restricthook;
2055b17ce1afSJed Brown   link->ctx          = ctx;
20560298fd71SBarry Smith   link->next         = NULL;
2057b17ce1afSJed Brown   *p                 = link;
2058b17ce1afSJed Brown   PetscFunctionReturn(0);
2059b17ce1afSJed Brown }
2060b17ce1afSJed Brown 
2061b17ce1afSJed Brown #undef __FUNCT__
2062b17ce1afSJed Brown #define __FUNCT__ "DMRestrict"
2063b17ce1afSJed Brown /*@
2064b17ce1afSJed Brown    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
2065b17ce1afSJed Brown 
2066b17ce1afSJed Brown    Collective if any hooks are
2067b17ce1afSJed Brown 
2068b17ce1afSJed Brown    Input Arguments:
2069b17ce1afSJed Brown +  fine - finer DM to use as a base
2070b17ce1afSJed Brown .  restrct - restriction matrix, apply using MatRestrict()
2071b17ce1afSJed Brown .  inject - injection matrix, also use MatRestrict()
2072b17ce1afSJed Brown -  coarse - coarer DM to update
2073b17ce1afSJed Brown 
2074b17ce1afSJed Brown    Level: developer
2075b17ce1afSJed Brown 
2076b17ce1afSJed Brown .seealso: DMCoarsenHookAdd(), MatRestrict()
2077b17ce1afSJed Brown @*/
2078b17ce1afSJed Brown PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
2079b17ce1afSJed Brown {
2080b17ce1afSJed Brown   PetscErrorCode    ierr;
2081b17ce1afSJed Brown   DMCoarsenHookLink link;
2082b17ce1afSJed Brown 
2083b17ce1afSJed Brown   PetscFunctionBegin;
2084b17ce1afSJed Brown   for (link=fine->coarsenhook; link; link=link->next) {
20858865f1eaSKarl Rupp     if (link->restricthook) {
20868865f1eaSKarl Rupp       ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);
20878865f1eaSKarl Rupp     }
2088b17ce1afSJed Brown   }
208947c6ae99SBarry Smith   PetscFunctionReturn(0);
209047c6ae99SBarry Smith }
209147c6ae99SBarry Smith 
209247c6ae99SBarry Smith #undef __FUNCT__
2093be081cd6SPeter Brune #define __FUNCT__ "DMSubDomainHookAdd"
2094bb9467b5SJed Brown /*@C
2095be081cd6SPeter Brune    DMSubDomainHookAdd - adds a callback to be run when restricting a problem to the coarse grid
20965dbd56e3SPeter Brune 
20975dbd56e3SPeter Brune    Logically Collective
20985dbd56e3SPeter Brune 
20995dbd56e3SPeter Brune    Input Arguments:
21005dbd56e3SPeter Brune +  global - global DM
2101ec4806b8SPeter Brune .  ddhook - function to run to pass data to the decomposition DM upon its creation
21025dbd56e3SPeter Brune .  restricthook - function to run to update data on block solve (at the beginning of the block solve)
21030298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
21045dbd56e3SPeter Brune 
2105ec4806b8SPeter Brune 
2106ec4806b8SPeter Brune    Calling sequence for ddhook:
2107ec4806b8SPeter Brune $    ddhook(DM global,DM block,void *ctx)
2108ec4806b8SPeter Brune 
2109ec4806b8SPeter Brune +  global - global DM
2110ec4806b8SPeter Brune .  block  - block DM
2111ec4806b8SPeter Brune -  ctx - optional user-defined function context
2112ec4806b8SPeter Brune 
21135dbd56e3SPeter Brune    Calling sequence for restricthook:
2114ec4806b8SPeter Brune $    restricthook(DM global,VecScatter out,VecScatter in,DM block,void *ctx)
21155dbd56e3SPeter Brune 
21165dbd56e3SPeter Brune +  global - global DM
21175dbd56e3SPeter Brune .  out    - scatter to the outer (with ghost and overlap points) block vector
21185dbd56e3SPeter Brune .  in     - scatter to block vector values only owned locally
2119ec4806b8SPeter Brune .  block  - block DM
21205dbd56e3SPeter Brune -  ctx - optional user-defined function context
21215dbd56e3SPeter Brune 
21225dbd56e3SPeter Brune    Level: advanced
21235dbd56e3SPeter Brune 
21245dbd56e3SPeter Brune    Notes:
2125ec4806b8SPeter Brune    This function is only needed if auxiliary data needs to be set up on subdomain DMs.
21265dbd56e3SPeter Brune 
21275dbd56e3SPeter Brune    If this function is called multiple times, the hooks will be run in the order they are added.
21285dbd56e3SPeter Brune 
21295dbd56e3SPeter Brune    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
2130ec4806b8SPeter Brune    extract the global information from its context (instead of from the SNES).
21315dbd56e3SPeter Brune 
2132bb9467b5SJed Brown    This function is currently not available from Fortran.
2133bb9467b5SJed Brown 
21345dbd56e3SPeter Brune .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
21355dbd56e3SPeter Brune @*/
2136be081cd6SPeter Brune PetscErrorCode DMSubDomainHookAdd(DM global,PetscErrorCode (*ddhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*),void *ctx)
21375dbd56e3SPeter Brune {
21385dbd56e3SPeter Brune   PetscErrorCode      ierr;
2139be081cd6SPeter Brune   DMSubDomainHookLink link,*p;
21405dbd56e3SPeter Brune 
21415dbd56e3SPeter Brune   PetscFunctionBegin;
21425dbd56e3SPeter Brune   PetscValidHeaderSpecific(global,DM_CLASSID,1);
2143be081cd6SPeter Brune   for (p=&global->subdomainhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
2144be081cd6SPeter Brune   ierr               = PetscMalloc(sizeof(struct _DMSubDomainHookLink),&link);CHKERRQ(ierr);
21455dbd56e3SPeter Brune   link->restricthook = restricthook;
2146be081cd6SPeter Brune   link->ddhook       = ddhook;
21475dbd56e3SPeter Brune   link->ctx          = ctx;
21480298fd71SBarry Smith   link->next         = NULL;
21495dbd56e3SPeter Brune   *p                 = link;
21505dbd56e3SPeter Brune   PetscFunctionReturn(0);
21515dbd56e3SPeter Brune }
21525dbd56e3SPeter Brune 
21535dbd56e3SPeter Brune #undef __FUNCT__
2154be081cd6SPeter Brune #define __FUNCT__ "DMSubDomainRestrict"
21555dbd56e3SPeter Brune /*@
2156be081cd6SPeter Brune    DMSubDomainRestrict - restricts user-defined problem data to a block DM by running hooks registered by DMSubDomainHookAdd()
21575dbd56e3SPeter Brune 
21585dbd56e3SPeter Brune    Collective if any hooks are
21595dbd56e3SPeter Brune 
21605dbd56e3SPeter Brune    Input Arguments:
21615dbd56e3SPeter Brune +  fine - finer DM to use as a base
2162be081cd6SPeter Brune .  oscatter - scatter from domain global vector filling subdomain global vector with overlap
2163be081cd6SPeter Brune .  gscatter - scatter from domain global vector filling subdomain local vector with ghosts
21645dbd56e3SPeter Brune -  coarse - coarer DM to update
21655dbd56e3SPeter Brune 
21665dbd56e3SPeter Brune    Level: developer
21675dbd56e3SPeter Brune 
21685dbd56e3SPeter Brune .seealso: DMCoarsenHookAdd(), MatRestrict()
21695dbd56e3SPeter Brune @*/
2170be081cd6SPeter Brune PetscErrorCode DMSubDomainRestrict(DM global,VecScatter oscatter,VecScatter gscatter,DM subdm)
21715dbd56e3SPeter Brune {
21725dbd56e3SPeter Brune   PetscErrorCode      ierr;
2173be081cd6SPeter Brune   DMSubDomainHookLink link;
21745dbd56e3SPeter Brune 
21755dbd56e3SPeter Brune   PetscFunctionBegin;
2176be081cd6SPeter Brune   for (link=global->subdomainhook; link; link=link->next) {
21778865f1eaSKarl Rupp     if (link->restricthook) {
21788865f1eaSKarl Rupp       ierr = (*link->restricthook)(global,oscatter,gscatter,subdm,link->ctx);CHKERRQ(ierr);
21798865f1eaSKarl Rupp     }
21805dbd56e3SPeter Brune   }
21815dbd56e3SPeter Brune   PetscFunctionReturn(0);
21825dbd56e3SPeter Brune }
21835dbd56e3SPeter Brune 
21845dbd56e3SPeter Brune #undef __FUNCT__
21855fe1f584SPeter Brune #define __FUNCT__ "DMGetCoarsenLevel"
21865fe1f584SPeter Brune /*@
21876a7d9d85SPeter Brune     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
21885fe1f584SPeter Brune 
21895fe1f584SPeter Brune     Not Collective
21905fe1f584SPeter Brune 
21915fe1f584SPeter Brune     Input Parameter:
21925fe1f584SPeter Brune .   dm - the DM object
21935fe1f584SPeter Brune 
21945fe1f584SPeter Brune     Output Parameter:
21956a7d9d85SPeter Brune .   level - number of coarsenings
21965fe1f584SPeter Brune 
21975fe1f584SPeter Brune     Level: developer
21985fe1f584SPeter Brune 
21996a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
22005fe1f584SPeter Brune 
22015fe1f584SPeter Brune @*/
22025fe1f584SPeter Brune PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
22035fe1f584SPeter Brune {
22045fe1f584SPeter Brune   PetscFunctionBegin;
22055fe1f584SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
22065fe1f584SPeter Brune   *level = dm->leveldown;
22075fe1f584SPeter Brune   PetscFunctionReturn(0);
22085fe1f584SPeter Brune }
22095fe1f584SPeter Brune 
22105fe1f584SPeter Brune 
22115fe1f584SPeter Brune 
22125fe1f584SPeter Brune #undef __FUNCT__
221347c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy"
221447c6ae99SBarry Smith /*@C
221547c6ae99SBarry Smith     DMRefineHierarchy - Refines a DM object, all levels at once
221647c6ae99SBarry Smith 
221747c6ae99SBarry Smith     Collective on DM
221847c6ae99SBarry Smith 
221947c6ae99SBarry Smith     Input Parameter:
222047c6ae99SBarry Smith +   dm - the DM object
222147c6ae99SBarry Smith -   nlevels - the number of levels of refinement
222247c6ae99SBarry Smith 
222347c6ae99SBarry Smith     Output Parameter:
222447c6ae99SBarry Smith .   dmf - the refined DM hierarchy
222547c6ae99SBarry Smith 
222647c6ae99SBarry Smith     Level: developer
222747c6ae99SBarry Smith 
2228e727c939SJed Brown .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
222947c6ae99SBarry Smith 
223047c6ae99SBarry Smith @*/
22317087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
223247c6ae99SBarry Smith {
223347c6ae99SBarry Smith   PetscErrorCode ierr;
223447c6ae99SBarry Smith 
223547c6ae99SBarry Smith   PetscFunctionBegin;
2236171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2237ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
223847c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
223947c6ae99SBarry Smith   if (dm->ops->refinehierarchy) {
224047c6ae99SBarry Smith     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
224147c6ae99SBarry Smith   } else if (dm->ops->refine) {
224247c6ae99SBarry Smith     PetscInt i;
224347c6ae99SBarry Smith 
2244ce94432eSBarry Smith     ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr);
224547c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
2246ce94432eSBarry Smith       ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr);
224747c6ae99SBarry Smith     }
2248ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
224947c6ae99SBarry Smith   PetscFunctionReturn(0);
225047c6ae99SBarry Smith }
225147c6ae99SBarry Smith 
225247c6ae99SBarry Smith #undef __FUNCT__
225347c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy"
225447c6ae99SBarry Smith /*@C
225547c6ae99SBarry Smith     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
225647c6ae99SBarry Smith 
225747c6ae99SBarry Smith     Collective on DM
225847c6ae99SBarry Smith 
225947c6ae99SBarry Smith     Input Parameter:
226047c6ae99SBarry Smith +   dm - the DM object
226147c6ae99SBarry Smith -   nlevels - the number of levels of coarsening
226247c6ae99SBarry Smith 
226347c6ae99SBarry Smith     Output Parameter:
226447c6ae99SBarry Smith .   dmc - the coarsened DM hierarchy
226547c6ae99SBarry Smith 
226647c6ae99SBarry Smith     Level: developer
226747c6ae99SBarry Smith 
2268e727c939SJed Brown .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
226947c6ae99SBarry Smith 
227047c6ae99SBarry Smith @*/
22717087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
227247c6ae99SBarry Smith {
227347c6ae99SBarry Smith   PetscErrorCode ierr;
227447c6ae99SBarry Smith 
227547c6ae99SBarry Smith   PetscFunctionBegin;
2276171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2277ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
227847c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
227947c6ae99SBarry Smith   PetscValidPointer(dmc,3);
228047c6ae99SBarry Smith   if (dm->ops->coarsenhierarchy) {
228147c6ae99SBarry Smith     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
228247c6ae99SBarry Smith   } else if (dm->ops->coarsen) {
228347c6ae99SBarry Smith     PetscInt i;
228447c6ae99SBarry Smith 
2285ce94432eSBarry Smith     ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr);
228647c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
2287ce94432eSBarry Smith       ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr);
228847c6ae99SBarry Smith     }
2289ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
229047c6ae99SBarry Smith   PetscFunctionReturn(0);
229147c6ae99SBarry Smith }
229247c6ae99SBarry Smith 
229347c6ae99SBarry Smith #undef __FUNCT__
2294e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates"
229547c6ae99SBarry Smith /*@
2296e727c939SJed Brown    DMCreateAggregates - Gets the aggregates that map between
229747c6ae99SBarry Smith    grids associated with two DMs.
229847c6ae99SBarry Smith 
229947c6ae99SBarry Smith    Collective on DM
230047c6ae99SBarry Smith 
230147c6ae99SBarry Smith    Input Parameters:
230247c6ae99SBarry Smith +  dmc - the coarse grid DM
230347c6ae99SBarry Smith -  dmf - the fine grid DM
230447c6ae99SBarry Smith 
230547c6ae99SBarry Smith    Output Parameters:
230647c6ae99SBarry Smith .  rest - the restriction matrix (transpose of the projection matrix)
230747c6ae99SBarry Smith 
230847c6ae99SBarry Smith    Level: intermediate
230947c6ae99SBarry Smith 
231047c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid
231147c6ae99SBarry Smith 
2312e727c939SJed Brown .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
231347c6ae99SBarry Smith @*/
2314e727c939SJed Brown PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
231547c6ae99SBarry Smith {
231647c6ae99SBarry Smith   PetscErrorCode ierr;
231747c6ae99SBarry Smith 
231847c6ae99SBarry Smith   PetscFunctionBegin;
2319171400e9SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
2320171400e9SBarry Smith   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
232147c6ae99SBarry Smith   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
232247c6ae99SBarry Smith   PetscFunctionReturn(0);
232347c6ae99SBarry Smith }
232447c6ae99SBarry Smith 
232547c6ae99SBarry Smith #undef __FUNCT__
23261a266240SBarry Smith #define __FUNCT__ "DMSetApplicationContextDestroy"
23271a266240SBarry Smith /*@C
23281a266240SBarry Smith     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
23291a266240SBarry Smith 
23301a266240SBarry Smith     Not Collective
23311a266240SBarry Smith 
23321a266240SBarry Smith     Input Parameters:
23331a266240SBarry Smith +   dm - the DM object
23341a266240SBarry Smith -   destroy - the destroy function
23351a266240SBarry Smith 
23361a266240SBarry Smith     Level: intermediate
23371a266240SBarry Smith 
2338e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
23391a266240SBarry Smith 
2340f07f9ceaSJed Brown @*/
23411a266240SBarry Smith PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
23421a266240SBarry Smith {
23431a266240SBarry Smith   PetscFunctionBegin;
2344171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
23451a266240SBarry Smith   dm->ctxdestroy = destroy;
23461a266240SBarry Smith   PetscFunctionReturn(0);
23471a266240SBarry Smith }
23481a266240SBarry Smith 
23491a266240SBarry Smith #undef __FUNCT__
23501b2093e4SBarry Smith #define __FUNCT__ "DMSetApplicationContext"
2351b07ff414SBarry Smith /*@
23521b2093e4SBarry Smith     DMSetApplicationContext - Set a user context into a DM object
235347c6ae99SBarry Smith 
235447c6ae99SBarry Smith     Not Collective
235547c6ae99SBarry Smith 
235647c6ae99SBarry Smith     Input Parameters:
235747c6ae99SBarry Smith +   dm - the DM object
235847c6ae99SBarry Smith -   ctx - the user context
235947c6ae99SBarry Smith 
236047c6ae99SBarry Smith     Level: intermediate
236147c6ae99SBarry Smith 
2362e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
236347c6ae99SBarry Smith 
236447c6ae99SBarry Smith @*/
23651b2093e4SBarry Smith PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
236647c6ae99SBarry Smith {
236747c6ae99SBarry Smith   PetscFunctionBegin;
2368171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
236947c6ae99SBarry Smith   dm->ctx = ctx;
237047c6ae99SBarry Smith   PetscFunctionReturn(0);
237147c6ae99SBarry Smith }
237247c6ae99SBarry Smith 
237347c6ae99SBarry Smith #undef __FUNCT__
23741b2093e4SBarry Smith #define __FUNCT__ "DMGetApplicationContext"
237547c6ae99SBarry Smith /*@
23761b2093e4SBarry Smith     DMGetApplicationContext - Gets a user context from a DM object
237747c6ae99SBarry Smith 
237847c6ae99SBarry Smith     Not Collective
237947c6ae99SBarry Smith 
238047c6ae99SBarry Smith     Input Parameter:
238147c6ae99SBarry Smith .   dm - the DM object
238247c6ae99SBarry Smith 
238347c6ae99SBarry Smith     Output Parameter:
238447c6ae99SBarry Smith .   ctx - the user context
238547c6ae99SBarry Smith 
238647c6ae99SBarry Smith     Level: intermediate
238747c6ae99SBarry Smith 
2388e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
238947c6ae99SBarry Smith 
239047c6ae99SBarry Smith @*/
23911b2093e4SBarry Smith PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
239247c6ae99SBarry Smith {
239347c6ae99SBarry Smith   PetscFunctionBegin;
2394171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
23951b2093e4SBarry Smith   *(void**)ctx = dm->ctx;
239647c6ae99SBarry Smith   PetscFunctionReturn(0);
239747c6ae99SBarry Smith }
239847c6ae99SBarry Smith 
239947c6ae99SBarry Smith #undef __FUNCT__
240008da532bSDmitry Karpeev #define __FUNCT__ "DMSetVariableBounds"
240108da532bSDmitry Karpeev /*@C
2402df3898eeSBarry Smith     DMSetVariableBounds - sets a function to compute the lower and upper bound vectors for SNESVI.
240308da532bSDmitry Karpeev 
240408da532bSDmitry Karpeev     Logically Collective on DM
240508da532bSDmitry Karpeev 
240608da532bSDmitry Karpeev     Input Parameter:
240708da532bSDmitry Karpeev +   dm - the DM object
24080298fd71SBarry Smith -   f - the function that computes variable bounds used by SNESVI (use NULL to cancel a previous function that was set)
240908da532bSDmitry Karpeev 
241008da532bSDmitry Karpeev     Level: intermediate
241108da532bSDmitry Karpeev 
2412835c3ec7SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(),
241308da532bSDmitry Karpeev          DMSetJacobian()
241408da532bSDmitry Karpeev 
241508da532bSDmitry Karpeev @*/
241608da532bSDmitry Karpeev PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
241708da532bSDmitry Karpeev {
241808da532bSDmitry Karpeev   PetscFunctionBegin;
241908da532bSDmitry Karpeev   dm->ops->computevariablebounds = f;
242008da532bSDmitry Karpeev   PetscFunctionReturn(0);
242108da532bSDmitry Karpeev }
242208da532bSDmitry Karpeev 
242308da532bSDmitry Karpeev #undef __FUNCT__
242408da532bSDmitry Karpeev #define __FUNCT__ "DMHasVariableBounds"
242508da532bSDmitry Karpeev /*@
242608da532bSDmitry Karpeev     DMHasVariableBounds - does the DM object have a variable bounds function?
242708da532bSDmitry Karpeev 
242808da532bSDmitry Karpeev     Not Collective
242908da532bSDmitry Karpeev 
243008da532bSDmitry Karpeev     Input Parameter:
243108da532bSDmitry Karpeev .   dm - the DM object to destroy
243208da532bSDmitry Karpeev 
243308da532bSDmitry Karpeev     Output Parameter:
243408da532bSDmitry Karpeev .   flg - PETSC_TRUE if the variable bounds function exists
243508da532bSDmitry Karpeev 
243608da532bSDmitry Karpeev     Level: developer
243708da532bSDmitry Karpeev 
243874e1e8c1SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
243908da532bSDmitry Karpeev 
244008da532bSDmitry Karpeev @*/
244108da532bSDmitry Karpeev PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
244208da532bSDmitry Karpeev {
244308da532bSDmitry Karpeev   PetscFunctionBegin;
244408da532bSDmitry Karpeev   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
244508da532bSDmitry Karpeev   PetscFunctionReturn(0);
244608da532bSDmitry Karpeev }
244708da532bSDmitry Karpeev 
244808da532bSDmitry Karpeev #undef __FUNCT__
244908da532bSDmitry Karpeev #define __FUNCT__ "DMComputeVariableBounds"
245008da532bSDmitry Karpeev /*@C
245108da532bSDmitry Karpeev     DMComputeVariableBounds - compute variable bounds used by SNESVI.
245208da532bSDmitry Karpeev 
245308da532bSDmitry Karpeev     Logically Collective on DM
245408da532bSDmitry Karpeev 
245508da532bSDmitry Karpeev     Input Parameters:
245608da532bSDmitry Karpeev +   dm - the DM object to destroy
245708da532bSDmitry Karpeev -   x  - current solution at which the bounds are computed
245808da532bSDmitry Karpeev 
245908da532bSDmitry Karpeev     Output parameters:
246008da532bSDmitry Karpeev +   xl - lower bound
246108da532bSDmitry Karpeev -   xu - upper bound
246208da532bSDmitry Karpeev 
246308da532bSDmitry Karpeev     Level: intermediate
246408da532bSDmitry Karpeev 
246574e1e8c1SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
246608da532bSDmitry Karpeev 
246708da532bSDmitry Karpeev @*/
246808da532bSDmitry Karpeev PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
246908da532bSDmitry Karpeev {
247008da532bSDmitry Karpeev   PetscErrorCode ierr;
24715fd66863SKarl Rupp 
247208da532bSDmitry Karpeev   PetscFunctionBegin;
247308da532bSDmitry Karpeev   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
247408da532bSDmitry Karpeev   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
247508da532bSDmitry Karpeev   if (dm->ops->computevariablebounds) {
247608da532bSDmitry Karpeev     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu);CHKERRQ(ierr);
24778865f1eaSKarl Rupp   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "This DM is incapable of computing variable bounds.");
247808da532bSDmitry Karpeev   PetscFunctionReturn(0);
247908da532bSDmitry Karpeev }
248008da532bSDmitry Karpeev 
248108da532bSDmitry Karpeev #undef __FUNCT__
2482b0ae01b7SPeter Brune #define __FUNCT__ "DMHasColoring"
2483b0ae01b7SPeter Brune /*@
2484b0ae01b7SPeter Brune     DMHasColoring - does the DM object have a method of providing a coloring?
2485b0ae01b7SPeter Brune 
2486b0ae01b7SPeter Brune     Not Collective
2487b0ae01b7SPeter Brune 
2488b0ae01b7SPeter Brune     Input Parameter:
2489b0ae01b7SPeter Brune .   dm - the DM object
2490b0ae01b7SPeter Brune 
2491b0ae01b7SPeter Brune     Output Parameter:
2492b0ae01b7SPeter Brune .   flg - PETSC_TRUE if the DM has facilities for DMCreateColoring().
2493b0ae01b7SPeter Brune 
2494b0ae01b7SPeter Brune     Level: developer
2495b0ae01b7SPeter Brune 
2496b0ae01b7SPeter Brune .seealso DMHasFunction(), DMCreateColoring()
2497b0ae01b7SPeter Brune 
2498b0ae01b7SPeter Brune @*/
2499b0ae01b7SPeter Brune PetscErrorCode  DMHasColoring(DM dm,PetscBool  *flg)
2500b0ae01b7SPeter Brune {
2501b0ae01b7SPeter Brune   PetscFunctionBegin;
2502b0ae01b7SPeter Brune   *flg =  (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE;
2503b0ae01b7SPeter Brune   PetscFunctionReturn(0);
2504b0ae01b7SPeter Brune }
2505b0ae01b7SPeter Brune 
2506b0ae01b7SPeter Brune #undef  __FUNCT__
250708da532bSDmitry Karpeev #define __FUNCT__ "DMSetVec"
2508748fac09SDmitry Karpeev /*@C
250908da532bSDmitry Karpeev     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
251008da532bSDmitry Karpeev 
251108da532bSDmitry Karpeev     Collective on DM
251208da532bSDmitry Karpeev 
251308da532bSDmitry Karpeev     Input Parameter:
251408da532bSDmitry Karpeev +   dm - the DM object
25150298fd71SBarry Smith -   x - location to compute residual and Jacobian, if NULL is passed to those routines; will be NULL for linear problems.
251608da532bSDmitry Karpeev 
251708da532bSDmitry Karpeev     Level: developer
251808da532bSDmitry Karpeev 
251974e1e8c1SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
252008da532bSDmitry Karpeev 
252108da532bSDmitry Karpeev @*/
252208da532bSDmitry Karpeev PetscErrorCode  DMSetVec(DM dm,Vec x)
252308da532bSDmitry Karpeev {
252408da532bSDmitry Karpeev   PetscErrorCode ierr;
25255fd66863SKarl Rupp 
252608da532bSDmitry Karpeev   PetscFunctionBegin;
252708da532bSDmitry Karpeev   if (x) {
252808da532bSDmitry Karpeev     if (!dm->x) {
252908da532bSDmitry Karpeev       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
253008da532bSDmitry Karpeev     }
253108da532bSDmitry Karpeev     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
25328865f1eaSKarl Rupp   } else if (dm->x) {
253308da532bSDmitry Karpeev     ierr = VecDestroy(&dm->x);CHKERRQ(ierr);
253408da532bSDmitry Karpeev   }
253508da532bSDmitry Karpeev   PetscFunctionReturn(0);
253608da532bSDmitry Karpeev }
253708da532bSDmitry Karpeev 
25380298fd71SBarry Smith PetscFunctionList DMList              = NULL;
2539264ace61SBarry Smith PetscBool         DMRegisterAllCalled = PETSC_FALSE;
2540264ace61SBarry Smith 
2541264ace61SBarry Smith #undef __FUNCT__
2542264ace61SBarry Smith #define __FUNCT__ "DMSetType"
2543264ace61SBarry Smith /*@C
2544264ace61SBarry Smith   DMSetType - Builds a DM, for a particular DM implementation.
2545264ace61SBarry Smith 
2546264ace61SBarry Smith   Collective on DM
2547264ace61SBarry Smith 
2548264ace61SBarry Smith   Input Parameters:
2549264ace61SBarry Smith + dm     - The DM object
2550264ace61SBarry Smith - method - The name of the DM type
2551264ace61SBarry Smith 
2552264ace61SBarry Smith   Options Database Key:
2553264ace61SBarry Smith . -dm_type <type> - Sets the DM type; use -help for a list of available types
2554264ace61SBarry Smith 
2555264ace61SBarry Smith   Notes:
2556e1589f56SBarry Smith   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
2557264ace61SBarry Smith 
2558264ace61SBarry Smith   Level: intermediate
2559264ace61SBarry Smith 
2560264ace61SBarry Smith .keywords: DM, set, type
2561264ace61SBarry Smith .seealso: DMGetType(), DMCreate()
2562264ace61SBarry Smith @*/
256319fd82e9SBarry Smith PetscErrorCode  DMSetType(DM dm, DMType method)
2564264ace61SBarry Smith {
2565264ace61SBarry Smith   PetscErrorCode (*r)(DM);
2566264ace61SBarry Smith   PetscBool      match;
2567264ace61SBarry Smith   PetscErrorCode ierr;
2568264ace61SBarry Smith 
2569264ace61SBarry Smith   PetscFunctionBegin;
2570264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2571251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
2572264ace61SBarry Smith   if (match) PetscFunctionReturn(0);
2573264ace61SBarry Smith 
2574607a6623SBarry Smith   if (!DMRegisterAllCalled) {ierr = DMRegisterAll();CHKERRQ(ierr);}
25751c9cd337SJed Brown   ierr = PetscFunctionListFind(DMList,method,&r);CHKERRQ(ierr);
2576ce94432eSBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
2577264ace61SBarry Smith 
2578264ace61SBarry Smith   if (dm->ops->destroy) {
2579264ace61SBarry Smith     ierr             = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
25800298fd71SBarry Smith     dm->ops->destroy = NULL;
2581264ace61SBarry Smith   }
2582264ace61SBarry Smith   ierr = (*r)(dm);CHKERRQ(ierr);
2583264ace61SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
2584264ace61SBarry Smith   PetscFunctionReturn(0);
2585264ace61SBarry Smith }
2586264ace61SBarry Smith 
2587264ace61SBarry Smith #undef __FUNCT__
2588264ace61SBarry Smith #define __FUNCT__ "DMGetType"
2589264ace61SBarry Smith /*@C
2590264ace61SBarry Smith   DMGetType - Gets the DM type name (as a string) from the DM.
2591264ace61SBarry Smith 
2592264ace61SBarry Smith   Not Collective
2593264ace61SBarry Smith 
2594264ace61SBarry Smith   Input Parameter:
2595264ace61SBarry Smith . dm  - The DM
2596264ace61SBarry Smith 
2597264ace61SBarry Smith   Output Parameter:
2598264ace61SBarry Smith . type - The DM type name
2599264ace61SBarry Smith 
2600264ace61SBarry Smith   Level: intermediate
2601264ace61SBarry Smith 
2602264ace61SBarry Smith .keywords: DM, get, type, name
2603264ace61SBarry Smith .seealso: DMSetType(), DMCreate()
2604264ace61SBarry Smith @*/
260519fd82e9SBarry Smith PetscErrorCode  DMGetType(DM dm, DMType *type)
2606264ace61SBarry Smith {
2607264ace61SBarry Smith   PetscErrorCode ierr;
2608264ace61SBarry Smith 
2609264ace61SBarry Smith   PetscFunctionBegin;
2610264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2611c959eef4SJed Brown   PetscValidPointer(type,2);
2612264ace61SBarry Smith   if (!DMRegisterAllCalled) {
2613607a6623SBarry Smith     ierr = DMRegisterAll();CHKERRQ(ierr);
2614264ace61SBarry Smith   }
2615264ace61SBarry Smith   *type = ((PetscObject)dm)->type_name;
2616264ace61SBarry Smith   PetscFunctionReturn(0);
2617264ace61SBarry Smith }
2618264ace61SBarry Smith 
261967a56275SMatthew G Knepley #undef __FUNCT__
262067a56275SMatthew G Knepley #define __FUNCT__ "DMConvert"
262167a56275SMatthew G Knepley /*@C
262267a56275SMatthew G Knepley   DMConvert - Converts a DM to another DM, either of the same or different type.
262367a56275SMatthew G Knepley 
262467a56275SMatthew G Knepley   Collective on DM
262567a56275SMatthew G Knepley 
262667a56275SMatthew G Knepley   Input Parameters:
262767a56275SMatthew G Knepley + dm - the DM
262867a56275SMatthew G Knepley - newtype - new DM type (use "same" for the same type)
262967a56275SMatthew G Knepley 
263067a56275SMatthew G Knepley   Output Parameter:
263167a56275SMatthew G Knepley . M - pointer to new DM
263267a56275SMatthew G Knepley 
263367a56275SMatthew G Knepley   Notes:
263467a56275SMatthew G Knepley   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
263567a56275SMatthew G Knepley   the MPI communicator of the generated DM is always the same as the communicator
263667a56275SMatthew G Knepley   of the input DM.
263767a56275SMatthew G Knepley 
263867a56275SMatthew G Knepley   Level: intermediate
263967a56275SMatthew G Knepley 
264067a56275SMatthew G Knepley .seealso: DMCreate()
264167a56275SMatthew G Knepley @*/
264219fd82e9SBarry Smith PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M)
264367a56275SMatthew G Knepley {
264467a56275SMatthew G Knepley   DM             B;
264567a56275SMatthew G Knepley   char           convname[256];
264667a56275SMatthew G Knepley   PetscBool      sametype, issame;
264767a56275SMatthew G Knepley   PetscErrorCode ierr;
264867a56275SMatthew G Knepley 
264967a56275SMatthew G Knepley   PetscFunctionBegin;
265067a56275SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
265167a56275SMatthew G Knepley   PetscValidType(dm,1);
265267a56275SMatthew G Knepley   PetscValidPointer(M,3);
2653251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
265467a56275SMatthew G Knepley   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
265567a56275SMatthew G Knepley   {
26560298fd71SBarry Smith     PetscErrorCode (*conv)(DM, DMType, DM*) = NULL;
265767a56275SMatthew G Knepley 
265867a56275SMatthew G Knepley     /*
265967a56275SMatthew G Knepley        Order of precedence:
266067a56275SMatthew G Knepley        1) See if a specialized converter is known to the current DM.
266167a56275SMatthew G Knepley        2) See if a specialized converter is known to the desired DM class.
266267a56275SMatthew G Knepley        3) See if a good general converter is registered for the desired class
266367a56275SMatthew G Knepley        4) See if a good general converter is known for the current matrix.
266467a56275SMatthew G Knepley        5) Use a really basic converter.
266567a56275SMatthew G Knepley     */
266667a56275SMatthew G Knepley 
266767a56275SMatthew G Knepley     /* 1) See if a specialized converter is known to the current DM and the desired class */
266867a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
266967a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
267067a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
267167a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
267267a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
26730005d66cSJed Brown     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,&conv);CHKERRQ(ierr);
267467a56275SMatthew G Knepley     if (conv) goto foundconv;
267567a56275SMatthew G Knepley 
267667a56275SMatthew G Knepley     /* 2)  See if a specialized converter is known to the desired DM class. */
267782f516ccSBarry Smith     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &B);CHKERRQ(ierr);
267867a56275SMatthew G Knepley     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
267967a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
268067a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
268167a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
268267a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
268367a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
26840005d66cSJed Brown     ierr = PetscObjectQueryFunction((PetscObject)B,convname,&conv);CHKERRQ(ierr);
268567a56275SMatthew G Knepley     if (conv) {
2686fcfd50ebSBarry Smith       ierr = DMDestroy(&B);CHKERRQ(ierr);
268767a56275SMatthew G Knepley       goto foundconv;
268867a56275SMatthew G Knepley     }
268967a56275SMatthew G Knepley 
269067a56275SMatthew G Knepley #if 0
269167a56275SMatthew G Knepley     /* 3) See if a good general converter is registered for the desired class */
269267a56275SMatthew G Knepley     conv = B->ops->convertfrom;
2693fcfd50ebSBarry Smith     ierr = DMDestroy(&B);CHKERRQ(ierr);
269467a56275SMatthew G Knepley     if (conv) goto foundconv;
269567a56275SMatthew G Knepley 
269667a56275SMatthew G Knepley     /* 4) See if a good general converter is known for the current matrix */
269767a56275SMatthew G Knepley     if (dm->ops->convert) {
269867a56275SMatthew G Knepley       conv = dm->ops->convert;
269967a56275SMatthew G Knepley     }
270067a56275SMatthew G Knepley     if (conv) goto foundconv;
270167a56275SMatthew G Knepley #endif
270267a56275SMatthew G Knepley 
270367a56275SMatthew G Knepley     /* 5) Use a really basic converter. */
270482f516ccSBarry Smith     SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
270567a56275SMatthew G Knepley 
270667a56275SMatthew G Knepley foundconv:
270767a56275SMatthew G Knepley     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
270867a56275SMatthew G Knepley     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
270967a56275SMatthew G Knepley     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
271067a56275SMatthew G Knepley   }
271167a56275SMatthew G Knepley   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
271267a56275SMatthew G Knepley   PetscFunctionReturn(0);
271367a56275SMatthew G Knepley }
2714264ace61SBarry Smith 
2715264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
2716264ace61SBarry Smith 
2717264ace61SBarry Smith #undef __FUNCT__
2718264ace61SBarry Smith #define __FUNCT__ "DMRegister"
2719264ace61SBarry Smith /*@C
27201c84c290SBarry Smith   DMRegister -  Adds a new DM component implementation
27211c84c290SBarry Smith 
27221c84c290SBarry Smith   Not Collective
27231c84c290SBarry Smith 
27241c84c290SBarry Smith   Input Parameters:
27251c84c290SBarry Smith + name        - The name of a new user-defined creation routine
27261c84c290SBarry Smith - create_func - The creation routine itself
27271c84c290SBarry Smith 
27281c84c290SBarry Smith   Notes:
27291c84c290SBarry Smith   DMRegister() may be called multiple times to add several user-defined DMs
27301c84c290SBarry Smith 
27311c84c290SBarry Smith 
27321c84c290SBarry Smith   Sample usage:
27331c84c290SBarry Smith .vb
2734bdf89e91SBarry Smith     DMRegister("my_da", MyDMCreate);
27351c84c290SBarry Smith .ve
27361c84c290SBarry Smith 
27371c84c290SBarry Smith   Then, your DM type can be chosen with the procedural interface via
27381c84c290SBarry Smith .vb
27391c84c290SBarry Smith     DMCreate(MPI_Comm, DM *);
27401c84c290SBarry Smith     DMSetType(DM,"my_da");
27411c84c290SBarry Smith .ve
27421c84c290SBarry Smith    or at runtime via the option
27431c84c290SBarry Smith .vb
27441c84c290SBarry Smith     -da_type my_da
27451c84c290SBarry Smith .ve
2746264ace61SBarry Smith 
2747264ace61SBarry Smith   Level: advanced
27481c84c290SBarry Smith 
27491c84c290SBarry Smith .keywords: DM, register
2750bdf89e91SBarry Smith .seealso: DMRegisterAll(), DMRegisterDestroy()
27511c84c290SBarry Smith 
2752264ace61SBarry Smith @*/
2753bdf89e91SBarry Smith PetscErrorCode  DMRegister(const char sname[],PetscErrorCode (*function)(DM))
2754264ace61SBarry Smith {
2755264ace61SBarry Smith   PetscErrorCode ierr;
2756264ace61SBarry Smith 
2757264ace61SBarry Smith   PetscFunctionBegin;
2758a240a19fSJed Brown   ierr = PetscFunctionListAdd(&DMList,sname,function);CHKERRQ(ierr);
2759264ace61SBarry Smith   PetscFunctionReturn(0);
2760264ace61SBarry Smith }
2761264ace61SBarry Smith 
2762b859378eSBarry Smith #undef __FUNCT__
2763b859378eSBarry Smith #define __FUNCT__ "DMLoad"
2764b859378eSBarry Smith /*@C
276555849f57SBarry Smith   DMLoad - Loads a DM that has been stored in binary  with DMView().
2766b859378eSBarry Smith 
2767b859378eSBarry Smith   Collective on PetscViewer
2768b859378eSBarry Smith 
2769b859378eSBarry Smith   Input Parameters:
2770b859378eSBarry Smith + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2771b859378eSBarry Smith            some related function before a call to DMLoad().
2772b859378eSBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2773b859378eSBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2774b859378eSBarry Smith 
2775b859378eSBarry Smith    Level: intermediate
2776b859378eSBarry Smith 
2777b859378eSBarry Smith   Notes:
277855849f57SBarry Smith    The type is determined by the data in the file, any type set into the DM before this call is ignored.
2779b859378eSBarry Smith 
2780b859378eSBarry Smith   Notes for advanced users:
2781b859378eSBarry Smith   Most users should not need to know the details of the binary storage
2782b859378eSBarry Smith   format, since DMLoad() and DMView() completely hide these details.
2783b859378eSBarry Smith   But for anyone who's interested, the standard binary matrix storage
2784b859378eSBarry Smith   format is
2785b859378eSBarry Smith .vb
2786b859378eSBarry Smith      has not yet been determined
2787b859378eSBarry Smith .ve
2788b859378eSBarry Smith 
2789b859378eSBarry Smith .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2790b859378eSBarry Smith @*/
2791b859378eSBarry Smith PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2792b859378eSBarry Smith {
27939331c7a4SMatthew G. Knepley   PetscBool      isbinary, ishdf5;
2794b859378eSBarry Smith   PetscErrorCode ierr;
2795b859378eSBarry Smith 
2796b859378eSBarry Smith   PetscFunctionBegin;
2797b859378eSBarry Smith   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2798b859378eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
279932c0f0efSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
28009331c7a4SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
28019331c7a4SMatthew G. Knepley   if (isbinary) {
28029331c7a4SMatthew G. Knepley     PetscInt classid;
28039331c7a4SMatthew G. Knepley     char     type[256];
2804b859378eSBarry Smith 
280532c0f0efSBarry Smith     ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
28069200755eSBarry Smith     if (classid != DM_FILE_CLASSID) SETERRQ1(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not DM next in file, classid found %d",(int)classid);
280732c0f0efSBarry Smith     ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
280832c0f0efSBarry Smith     ierr = DMSetType(newdm, type);CHKERRQ(ierr);
28099331c7a4SMatthew G. Knepley     if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);}
28109331c7a4SMatthew G. Knepley   } else if (ishdf5) {
28119331c7a4SMatthew G. Knepley     if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);}
28129331c7a4SMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen() or PetscViewerHDF5Open()");
2813b859378eSBarry Smith   PetscFunctionReturn(0);
2814b859378eSBarry Smith }
2815b859378eSBarry Smith 
28167da65231SMatthew G Knepley /******************************** FEM Support **********************************/
28177da65231SMatthew G Knepley 
28187da65231SMatthew G Knepley #undef __FUNCT__
28197da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellVector"
2820a6dfd86eSKarl Rupp PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[])
2821a6dfd86eSKarl Rupp {
28221d47ebbbSSatish Balay   PetscInt       f;
28231b30c384SMatthew G Knepley   PetscErrorCode ierr;
28241b30c384SMatthew G Knepley 
28257da65231SMatthew G Knepley   PetscFunctionBegin;
282674778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
28271d47ebbbSSatish Balay   for (f = 0; f < len; ++f) {
282857622a8eSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", (double)PetscRealPart(x[f]));CHKERRQ(ierr);
28297da65231SMatthew G Knepley   }
28307da65231SMatthew G Knepley   PetscFunctionReturn(0);
28317da65231SMatthew G Knepley }
28327da65231SMatthew G Knepley 
28337da65231SMatthew G Knepley #undef __FUNCT__
28347da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellMatrix"
2835a6dfd86eSKarl Rupp PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
2836a6dfd86eSKarl Rupp {
28371b30c384SMatthew G Knepley   PetscInt       f, g;
28387da65231SMatthew G Knepley   PetscErrorCode ierr;
28397da65231SMatthew G Knepley 
28407da65231SMatthew G Knepley   PetscFunctionBegin;
284174778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
28421d47ebbbSSatish Balay   for (f = 0; f < rows; ++f) {
284374778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
28441d47ebbbSSatish Balay     for (g = 0; g < cols; ++g) {
2845e3556bceSMatthew G. Knepley       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5g", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
28467da65231SMatthew G Knepley     }
284774778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
28487da65231SMatthew G Knepley   }
28497da65231SMatthew G Knepley   PetscFunctionReturn(0);
28507da65231SMatthew G Knepley }
2851e7c4fc90SDmitry Karpeev 
2852970e74d5SMatthew G Knepley #undef __FUNCT__
2853e759306cSMatthew G. Knepley #define __FUNCT__ "DMPrintLocalVec"
28546113b454SMatthew G. Knepley PetscErrorCode DMPrintLocalVec(DM dm, const char name[], PetscReal tol, Vec X)
2855e759306cSMatthew G. Knepley {
2856e759306cSMatthew G. Knepley   PetscMPIInt    rank, numProcs;
2857e759306cSMatthew G. Knepley   PetscInt       p;
2858e759306cSMatthew G. Knepley   PetscErrorCode ierr;
2859e759306cSMatthew G. Knepley 
2860e759306cSMatthew G. Knepley   PetscFunctionBegin;
2861e759306cSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
2862e759306cSMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
2863e759306cSMatthew G. Knepley   ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s:\n", name);CHKERRQ(ierr);
2864e759306cSMatthew G. Knepley   for (p = 0; p < numProcs; ++p) {
2865e759306cSMatthew G. Knepley     if (p == rank) {
2866e759306cSMatthew G. Knepley       Vec x;
2867e759306cSMatthew G. Knepley 
2868e759306cSMatthew G. Knepley       ierr = VecDuplicate(X, &x);CHKERRQ(ierr);
2869e759306cSMatthew G. Knepley       ierr = VecCopy(X, x);CHKERRQ(ierr);
28706113b454SMatthew G. Knepley       ierr = VecChop(x, tol);CHKERRQ(ierr);
2871e759306cSMatthew G. Knepley       ierr = VecView(x, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2872e759306cSMatthew G. Knepley       ierr = VecDestroy(&x);CHKERRQ(ierr);
2873e759306cSMatthew G. Knepley       ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2874e759306cSMatthew G. Knepley     }
2875e759306cSMatthew G. Knepley     ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
2876e759306cSMatthew G. Knepley   }
2877e759306cSMatthew G. Knepley   PetscFunctionReturn(0);
2878e759306cSMatthew G. Knepley }
2879e759306cSMatthew G. Knepley 
2880e759306cSMatthew G. Knepley #undef __FUNCT__
288188ed4aceSMatthew G Knepley #define __FUNCT__ "DMGetDefaultSection"
288288ed4aceSMatthew G Knepley /*@
288388ed4aceSMatthew G Knepley   DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM.
288488ed4aceSMatthew G Knepley 
288588ed4aceSMatthew G Knepley   Input Parameter:
288688ed4aceSMatthew G Knepley . dm - The DM
288788ed4aceSMatthew G Knepley 
288888ed4aceSMatthew G Knepley   Output Parameter:
288988ed4aceSMatthew G Knepley . section - The PetscSection
289088ed4aceSMatthew G Knepley 
289188ed4aceSMatthew G Knepley   Level: intermediate
289288ed4aceSMatthew G Knepley 
289388ed4aceSMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
289488ed4aceSMatthew G Knepley 
289588ed4aceSMatthew G Knepley .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
289688ed4aceSMatthew G Knepley @*/
28970adebc6cSBarry Smith PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section)
28980adebc6cSBarry Smith {
2899fd59a867SMatthew G. Knepley   PetscErrorCode ierr;
2900fd59a867SMatthew G. Knepley 
290188ed4aceSMatthew G Knepley   PetscFunctionBegin;
290288ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
290388ed4aceSMatthew G Knepley   PetscValidPointer(section, 2);
2904fd59a867SMatthew G. Knepley   if (!dm->defaultSection && dm->ops->createdefaultsection) {ierr = (*dm->ops->createdefaultsection)(dm);CHKERRQ(ierr);}
290588ed4aceSMatthew G Knepley   *section = dm->defaultSection;
290688ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
290788ed4aceSMatthew G Knepley }
290888ed4aceSMatthew G Knepley 
290988ed4aceSMatthew G Knepley #undef __FUNCT__
291088ed4aceSMatthew G Knepley #define __FUNCT__ "DMSetDefaultSection"
291188ed4aceSMatthew G Knepley /*@
291288ed4aceSMatthew G Knepley   DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM.
291388ed4aceSMatthew G Knepley 
291488ed4aceSMatthew G Knepley   Input Parameters:
291588ed4aceSMatthew G Knepley + dm - The DM
291688ed4aceSMatthew G Knepley - section - The PetscSection
291788ed4aceSMatthew G Knepley 
291888ed4aceSMatthew G Knepley   Level: intermediate
291988ed4aceSMatthew G Knepley 
292088ed4aceSMatthew G Knepley   Note: Any existing Section will be destroyed
292188ed4aceSMatthew G Knepley 
292288ed4aceSMatthew G Knepley .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
292388ed4aceSMatthew G Knepley @*/
29240adebc6cSBarry Smith PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section)
29250adebc6cSBarry Smith {
2926af122d2aSMatthew G Knepley   PetscInt       numFields;
2927af122d2aSMatthew G Knepley   PetscInt       f;
292888ed4aceSMatthew G Knepley   PetscErrorCode ierr;
292988ed4aceSMatthew G Knepley 
293088ed4aceSMatthew G Knepley   PetscFunctionBegin;
293188ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
29321d799100SJed Brown   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
29331d799100SJed Brown   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
293488ed4aceSMatthew G Knepley   ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr);
293588ed4aceSMatthew G Knepley   dm->defaultSection = section;
2936af122d2aSMatthew G Knepley   ierr = PetscSectionGetNumFields(dm->defaultSection, &numFields);CHKERRQ(ierr);
2937af122d2aSMatthew G Knepley   if (numFields) {
2938af122d2aSMatthew G Knepley     ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr);
2939af122d2aSMatthew G Knepley     for (f = 0; f < numFields; ++f) {
29400f21e855SMatthew G. Knepley       PetscObject disc;
2941af122d2aSMatthew G Knepley       const char *name;
2942af122d2aSMatthew G Knepley 
2943af122d2aSMatthew G Knepley       ierr = PetscSectionGetFieldName(dm->defaultSection, f, &name);CHKERRQ(ierr);
29440f21e855SMatthew G. Knepley       ierr = DMGetField(dm, f, &disc);CHKERRQ(ierr);
29450f21e855SMatthew G. Knepley       ierr = PetscObjectSetName(disc, name);CHKERRQ(ierr);
2946af122d2aSMatthew G Knepley     }
2947af122d2aSMatthew G Knepley   }
29481d799100SJed Brown   /* The global section will be rebuilt in the next call to DMGetDefaultGlobalSection(). */
29491d799100SJed Brown   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
295088ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
295188ed4aceSMatthew G Knepley }
295288ed4aceSMatthew G Knepley 
295388ed4aceSMatthew G Knepley #undef __FUNCT__
2954*9435951eSToby Isaac #define __FUNCT__ "DMGetDefaultConstraints"
2955*9435951eSToby Isaac /*@
2956*9435951eSToby Isaac   DMGetDefaultConstraints - Get the PetscSection and Mat the specify the local constraint interpolation. See DMSetDefaultConstraints() for a description of the purpose of constraint interpolation.
2957*9435951eSToby Isaac 
2958*9435951eSToby Isaac   Input Parameter:
2959*9435951eSToby Isaac . dm - The DM
2960*9435951eSToby Isaac 
2961*9435951eSToby Isaac   Output Parameter:
2962*9435951eSToby Isaac + section - The PetscSection describing the range of the constraint matrix: relates rows of the constraint matrix to dofs of the default section.  Returns NULL if there are no local constraints.
2963*9435951eSToby Isaac - mat - The Mat that interpolates local constraints: its width should be the layout size of the default section.  Returns NULL if there are no local constraints.
2964*9435951eSToby Isaac 
2965*9435951eSToby Isaac   Level: advanced
2966*9435951eSToby Isaac 
2967*9435951eSToby Isaac   Note: This gets borrowed references, so the user should not destroy the PetscSection or the Mat.
2968*9435951eSToby Isaac 
2969*9435951eSToby Isaac .seealso: DMSetDefaultConstraints()
2970*9435951eSToby Isaac @*/
2971*9435951eSToby Isaac PetscErrorCode DMGetDefaultConstraints(DM dm, PetscSection *section, Mat *mat)
2972*9435951eSToby Isaac {
2973*9435951eSToby Isaac   PetscErrorCode ierr;
2974*9435951eSToby Isaac 
2975*9435951eSToby Isaac   PetscFunctionBegin;
2976*9435951eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2977*9435951eSToby Isaac   PetscValidPointer(section, 2);
2978*9435951eSToby Isaac   PetscValidPointer(mat, 3);
2979*9435951eSToby Isaac   if (!dm->defaultConstraintSection && !dm->defaultConstraintMat && dm->ops->createdefaultconstraints) {ierr = (*dm->ops->createdefaultconstraints)(dm);CHKERRQ(ierr);}
2980*9435951eSToby Isaac   *section = dm->defaultConstraintSection;
2981*9435951eSToby Isaac   *mat = dm->defaultConstraintMat;
2982*9435951eSToby Isaac   PetscFunctionReturn(0);
2983*9435951eSToby Isaac }
2984*9435951eSToby Isaac 
2985*9435951eSToby Isaac #undef __FUNCT__
2986*9435951eSToby Isaac #define __FUNCT__ "DMSetDefaultConstraints"
2987*9435951eSToby Isaac /*@
2988*9435951eSToby Isaac   DMSetDefaultConstraints - Set the PetscSection and Mat the specify the local constraint interpolation.
2989*9435951eSToby Isaac 
2990*9435951eSToby Isaac   If a constraint matrix is specified, then it is applied during DMGlobalToLocalEnd() when mode is INSERT_VALUES, INSERT_BC_VALUES, or INSERT_ALL_VALUES.  Without a constraint matrix, the local vector l returned by DMGlobalToLocalEnd() contains values that have been scattered from a global vector without modification; with a constraint matrix A, l is modified by computing c = A * l, l[s[i]] = c[i], where the scatter s is defined by the PetscSection returned by DMGetDefaultConstraintMatrix().
2991*9435951eSToby Isaac 
2992*9435951eSToby Isaac   If a constraint matrix is specified, then its adjoint is applied during DMLocalToGlobalBegin() when mode is ADD_VALUES, ADD_BC_VALUES, or ADD_ALL_VALUES.  Without a constraint matrix, the local vector l is accumulated into a global vector without modification; with a constraint matrix A, l is first modified by computing c[i] = l[s[i]], l[s[i]] = 0, l = l + A'*c, which is the adjoint of the operation described above.
2993*9435951eSToby Isaac 
2994*9435951eSToby Isaac   Input Parameters:
2995*9435951eSToby Isaac + dm - The DM
2996*9435951eSToby Isaac + section - The PetscSection describing the range of the constraint matrix: relates rows of the constraint matrix to dofs of the default section.
2997*9435951eSToby Isaac - mat - The Mat that interpolates local constraints: its width should be the layout size of the default section:  NULL indicates no constraints.
2998*9435951eSToby Isaac 
2999*9435951eSToby Isaac   Level: advanced
3000*9435951eSToby Isaac 
3001*9435951eSToby Isaac   Note: This increments the references of the PetscSection and the Mat, so they user can destroy them
3002*9435951eSToby Isaac 
3003*9435951eSToby Isaac .seealso: DMGetDefaultConstraints()
3004*9435951eSToby Isaac @*/
3005*9435951eSToby Isaac PetscErrorCode DMSetDefaultConstraints(DM dm, PetscSection section, Mat mat)
3006*9435951eSToby Isaac {
3007*9435951eSToby Isaac   PetscInt       numFields;
3008*9435951eSToby Isaac   PetscErrorCode ierr;
3009*9435951eSToby Isaac 
3010*9435951eSToby Isaac   PetscFunctionBegin;
3011*9435951eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3012*9435951eSToby Isaac   if (section) {PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);}
3013*9435951eSToby Isaac   if (mat) {PetscValidHeaderSpecific(mat,MAT_CLASSID,3);}
3014*9435951eSToby Isaac   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
3015*9435951eSToby Isaac   ierr = PetscSectionDestroy(&dm->defaultConstraintSection);CHKERRQ(ierr);
3016*9435951eSToby Isaac   dm->defaultConstraintSection = section;
3017*9435951eSToby Isaac   ierr = PetscSectionGetNumFields(dm->defaultSection, &numFields);CHKERRQ(ierr);
3018*9435951eSToby Isaac   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
3019*9435951eSToby Isaac   ierr = MatDestroy(&dm->defaultConstraintMat);CHKERRQ(ierr);
3020*9435951eSToby Isaac   dm->defaultConstraintMat = mat;
3021*9435951eSToby Isaac   PetscFunctionReturn(0);
3022*9435951eSToby Isaac }
3023*9435951eSToby Isaac 
3024*9435951eSToby Isaac #undef __FUNCT__
302588ed4aceSMatthew G Knepley #define __FUNCT__ "DMGetDefaultGlobalSection"
302688ed4aceSMatthew G Knepley /*@
302788ed4aceSMatthew G Knepley   DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM.
302888ed4aceSMatthew G Knepley 
30298b1ab98fSJed Brown   Collective on DM
30308b1ab98fSJed Brown 
303188ed4aceSMatthew G Knepley   Input Parameter:
303288ed4aceSMatthew G Knepley . dm - The DM
303388ed4aceSMatthew G Knepley 
303488ed4aceSMatthew G Knepley   Output Parameter:
303588ed4aceSMatthew G Knepley . section - The PetscSection
303688ed4aceSMatthew G Knepley 
303788ed4aceSMatthew G Knepley   Level: intermediate
303888ed4aceSMatthew G Knepley 
303988ed4aceSMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
304088ed4aceSMatthew G Knepley 
304188ed4aceSMatthew G Knepley .seealso: DMSetDefaultSection(), DMGetDefaultSection()
304288ed4aceSMatthew G Knepley @*/
30430adebc6cSBarry Smith PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section)
30440adebc6cSBarry Smith {
304588ed4aceSMatthew G Knepley   PetscErrorCode ierr;
304688ed4aceSMatthew G Knepley 
304788ed4aceSMatthew G Knepley   PetscFunctionBegin;
304888ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
304988ed4aceSMatthew G Knepley   PetscValidPointer(section, 2);
305088ed4aceSMatthew G Knepley   if (!dm->defaultGlobalSection) {
3051fd59a867SMatthew G. Knepley     PetscSection s;
3052fd59a867SMatthew G. Knepley 
3053fd59a867SMatthew G. Knepley     ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
3054fd59a867SMatthew G. Knepley     if (!s)  SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSection in order to create a global PetscSection");
3055fd59a867SMatthew G. Knepley     if (!dm->sf) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSF in order to create a global PetscSection");
3056fd59a867SMatthew G. Knepley     ierr = PetscSectionCreateGlobalSection(s, dm->sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr);
3057cf06b437SMatthew G. Knepley     ierr = PetscLayoutDestroy(&dm->map);CHKERRQ(ierr);
3058ce94432eSBarry Smith     ierr = PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm), dm->defaultGlobalSection, &dm->map);CHKERRQ(ierr);
305988ed4aceSMatthew G Knepley   }
306088ed4aceSMatthew G Knepley   *section = dm->defaultGlobalSection;
306188ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
306288ed4aceSMatthew G Knepley }
306388ed4aceSMatthew G Knepley 
306488ed4aceSMatthew G Knepley #undef __FUNCT__
3065b21d0597SMatthew G Knepley #define __FUNCT__ "DMSetDefaultGlobalSection"
3066b21d0597SMatthew G Knepley /*@
3067b21d0597SMatthew G Knepley   DMSetDefaultGlobalSection - Set the PetscSection encoding the global data layout for the DM.
3068b21d0597SMatthew G Knepley 
3069b21d0597SMatthew G Knepley   Input Parameters:
3070b21d0597SMatthew G Knepley + dm - The DM
30715080bbdbSMatthew G Knepley - section - The PetscSection, or NULL
3072b21d0597SMatthew G Knepley 
3073b21d0597SMatthew G Knepley   Level: intermediate
3074b21d0597SMatthew G Knepley 
3075b21d0597SMatthew G Knepley   Note: Any existing Section will be destroyed
3076b21d0597SMatthew G Knepley 
3077b21d0597SMatthew G Knepley .seealso: DMGetDefaultGlobalSection(), DMSetDefaultSection()
3078b21d0597SMatthew G Knepley @*/
30790adebc6cSBarry Smith PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection section)
30800adebc6cSBarry Smith {
3081b21d0597SMatthew G Knepley   PetscErrorCode ierr;
3082b21d0597SMatthew G Knepley 
3083b21d0597SMatthew G Knepley   PetscFunctionBegin;
3084b21d0597SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
30855080bbdbSMatthew G Knepley   if (section) PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
30861d799100SJed Brown   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
3087b21d0597SMatthew G Knepley   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
3088b21d0597SMatthew G Knepley   dm->defaultGlobalSection = section;
3089b21d0597SMatthew G Knepley   PetscFunctionReturn(0);
3090b21d0597SMatthew G Knepley }
3091b21d0597SMatthew G Knepley 
3092b21d0597SMatthew G Knepley #undef __FUNCT__
309388ed4aceSMatthew G Knepley #define __FUNCT__ "DMGetDefaultSF"
309488ed4aceSMatthew G Knepley /*@
309588ed4aceSMatthew G Knepley   DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set,
309688ed4aceSMatthew G Knepley   it is created from the default PetscSection layouts in the DM.
309788ed4aceSMatthew G Knepley 
309888ed4aceSMatthew G Knepley   Input Parameter:
309988ed4aceSMatthew G Knepley . dm - The DM
310088ed4aceSMatthew G Knepley 
310188ed4aceSMatthew G Knepley   Output Parameter:
310288ed4aceSMatthew G Knepley . sf - The PetscSF
310388ed4aceSMatthew G Knepley 
310488ed4aceSMatthew G Knepley   Level: intermediate
310588ed4aceSMatthew G Knepley 
310688ed4aceSMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
310788ed4aceSMatthew G Knepley 
310888ed4aceSMatthew G Knepley .seealso: DMSetDefaultSF(), DMCreateDefaultSF()
310988ed4aceSMatthew G Knepley @*/
31100adebc6cSBarry Smith PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf)
31110adebc6cSBarry Smith {
311288ed4aceSMatthew G Knepley   PetscInt       nroots;
311388ed4aceSMatthew G Knepley   PetscErrorCode ierr;
311488ed4aceSMatthew G Knepley 
311588ed4aceSMatthew G Knepley   PetscFunctionBegin;
311688ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
311788ed4aceSMatthew G Knepley   PetscValidPointer(sf, 2);
31180298fd71SBarry Smith   ierr = PetscSFGetGraph(dm->defaultSF, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
311988ed4aceSMatthew G Knepley   if (nroots < 0) {
312088ed4aceSMatthew G Knepley     PetscSection section, gSection;
312188ed4aceSMatthew G Knepley 
312288ed4aceSMatthew G Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
312331ea6d37SMatthew G Knepley     if (section) {
312488ed4aceSMatthew G Knepley       ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
312588ed4aceSMatthew G Knepley       ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr);
312631ea6d37SMatthew G Knepley     } else {
31270298fd71SBarry Smith       *sf = NULL;
312831ea6d37SMatthew G Knepley       PetscFunctionReturn(0);
312931ea6d37SMatthew G Knepley     }
313088ed4aceSMatthew G Knepley   }
313188ed4aceSMatthew G Knepley   *sf = dm->defaultSF;
313288ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
313388ed4aceSMatthew G Knepley }
313488ed4aceSMatthew G Knepley 
313588ed4aceSMatthew G Knepley #undef __FUNCT__
313688ed4aceSMatthew G Knepley #define __FUNCT__ "DMSetDefaultSF"
313788ed4aceSMatthew G Knepley /*@
313888ed4aceSMatthew G Knepley   DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM
313988ed4aceSMatthew G Knepley 
314088ed4aceSMatthew G Knepley   Input Parameters:
314188ed4aceSMatthew G Knepley + dm - The DM
314288ed4aceSMatthew G Knepley - sf - The PetscSF
314388ed4aceSMatthew G Knepley 
314488ed4aceSMatthew G Knepley   Level: intermediate
314588ed4aceSMatthew G Knepley 
314688ed4aceSMatthew G Knepley   Note: Any previous SF is destroyed
314788ed4aceSMatthew G Knepley 
314888ed4aceSMatthew G Knepley .seealso: DMGetDefaultSF(), DMCreateDefaultSF()
314988ed4aceSMatthew G Knepley @*/
31500adebc6cSBarry Smith PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf)
31510adebc6cSBarry Smith {
315288ed4aceSMatthew G Knepley   PetscErrorCode ierr;
315388ed4aceSMatthew G Knepley 
315488ed4aceSMatthew G Knepley   PetscFunctionBegin;
315588ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
315688ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
315788ed4aceSMatthew G Knepley   ierr          = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr);
315888ed4aceSMatthew G Knepley   dm->defaultSF = sf;
315988ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
316088ed4aceSMatthew G Knepley }
316188ed4aceSMatthew G Knepley 
316288ed4aceSMatthew G Knepley #undef __FUNCT__
316388ed4aceSMatthew G Knepley #define __FUNCT__ "DMCreateDefaultSF"
316488ed4aceSMatthew G Knepley /*@C
316588ed4aceSMatthew G Knepley   DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections
316688ed4aceSMatthew G Knepley   describing the data layout.
316788ed4aceSMatthew G Knepley 
316888ed4aceSMatthew G Knepley   Input Parameters:
316988ed4aceSMatthew G Knepley + dm - The DM
317088ed4aceSMatthew G Knepley . localSection - PetscSection describing the local data layout
317188ed4aceSMatthew G Knepley - globalSection - PetscSection describing the global data layout
317288ed4aceSMatthew G Knepley 
317388ed4aceSMatthew G Knepley   Level: intermediate
317488ed4aceSMatthew G Knepley 
317588ed4aceSMatthew G Knepley .seealso: DMGetDefaultSF(), DMSetDefaultSF()
317688ed4aceSMatthew G Knepley @*/
317788ed4aceSMatthew G Knepley PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection)
317888ed4aceSMatthew G Knepley {
317982f516ccSBarry Smith   MPI_Comm       comm;
318088ed4aceSMatthew G Knepley   PetscLayout    layout;
318188ed4aceSMatthew G Knepley   const PetscInt *ranges;
318288ed4aceSMatthew G Knepley   PetscInt       *local;
318388ed4aceSMatthew G Knepley   PetscSFNode    *remote;
3184ecd73843SMatthew G. Knepley   PetscInt       pStart, pEnd, p, nroots, nleaves = 0, l;
318588ed4aceSMatthew G Knepley   PetscMPIInt    size, rank;
318688ed4aceSMatthew G Knepley   PetscErrorCode ierr;
318788ed4aceSMatthew G Knepley 
318888ed4aceSMatthew G Knepley   PetscFunctionBegin;
318982f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
319088ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
319188ed4aceSMatthew G Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
319288ed4aceSMatthew G Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
319388ed4aceSMatthew G Knepley   ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr);
319488ed4aceSMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr);
319588ed4aceSMatthew G Knepley   ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
319688ed4aceSMatthew G Knepley   ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
319788ed4aceSMatthew G Knepley   ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr);
319888ed4aceSMatthew G Knepley   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
319988ed4aceSMatthew G Knepley   ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
3200ecd73843SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
32016636e97aSMatthew G Knepley     PetscInt gdof, gcdof;
320288ed4aceSMatthew G Knepley 
32036636e97aSMatthew G Knepley     ierr     = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
32046636e97aSMatthew G Knepley     ierr     = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
3205235fbf56SMatthew G. Knepley     if (gcdof > (gdof < 0 ? -(gdof+1) : gdof)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d has %d constraints > %d dof", p, gcdof, (gdof < 0 ? -(gdof+1) : gdof));
32066636e97aSMatthew G Knepley     nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
320788ed4aceSMatthew G Knepley   }
3208785e854fSJed Brown   ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr);
3209785e854fSJed Brown   ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr);
321088ed4aceSMatthew G Knepley   for (p = pStart, l = 0; p < pEnd; ++p) {
32111f588964SMatthew G Knepley     const PetscInt *cind;
32126636e97aSMatthew G Knepley     PetscInt       dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
321388ed4aceSMatthew G Knepley 
321488ed4aceSMatthew G Knepley     ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr);
321588ed4aceSMatthew G Knepley     ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr);
321688ed4aceSMatthew G Knepley     ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr);
321788ed4aceSMatthew G Knepley     ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr);
321888ed4aceSMatthew G Knepley     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
32196636e97aSMatthew G Knepley     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
322088ed4aceSMatthew G Knepley     ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
32216636e97aSMatthew G Knepley     if (!gdof) continue; /* Censored point */
32226636e97aSMatthew G Knepley     gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
32236636e97aSMatthew G Knepley     if (gsize != dof-cdof) {
3224057b4bcdSMatthew G Knepley       if (gsize != dof) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Global dof %d for point %d is neither the constrained size %d, nor the unconstrained %d", gsize, p, dof-cdof, dof);
32256636e97aSMatthew G Knepley       cdof = 0; /* Ignore constraints */
32266636e97aSMatthew G Knepley     }
322788ed4aceSMatthew G Knepley     for (d = 0, c = 0; d < dof; ++d) {
322888ed4aceSMatthew G Knepley       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
322988ed4aceSMatthew G Knepley       local[l+d-c] = off+d;
323088ed4aceSMatthew G Knepley     }
323188ed4aceSMatthew G Knepley     if (gdof < 0) {
32326636e97aSMatthew G Knepley       for (d = 0; d < gsize; ++d, ++l) {
323388ed4aceSMatthew G Knepley         PetscInt offset = -(goff+1) + d, r;
323488ed4aceSMatthew G Knepley 
323505376888SMatthew G. Knepley         ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr);
323631d3f06eSJed Brown         if (r < 0) r = -(r+2);
323705376888SMatthew G. Knepley         if ((r < 0) || (r >= size)) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d mapped to invalid process %d (%d, %d)", p, r, gdof, goff);
323888ed4aceSMatthew G Knepley         remote[l].rank  = r;
323988ed4aceSMatthew G Knepley         remote[l].index = offset - ranges[r];
324088ed4aceSMatthew G Knepley       }
324188ed4aceSMatthew G Knepley     } else {
32426636e97aSMatthew G Knepley       for (d = 0; d < gsize; ++d, ++l) {
324388ed4aceSMatthew G Knepley         remote[l].rank  = rank;
324488ed4aceSMatthew G Knepley         remote[l].index = goff+d - ranges[rank];
324588ed4aceSMatthew G Knepley       }
324688ed4aceSMatthew G Knepley     }
324788ed4aceSMatthew G Knepley   }
32486636e97aSMatthew G Knepley   if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %d != nleaves %d", l, nleaves);
324988ed4aceSMatthew G Knepley   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
325088ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
325188ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
325288ed4aceSMatthew G Knepley }
3253af122d2aSMatthew G Knepley 
3254af122d2aSMatthew G Knepley #undef __FUNCT__
3255b21d0597SMatthew G Knepley #define __FUNCT__ "DMGetPointSF"
3256b21d0597SMatthew G Knepley /*@
3257b21d0597SMatthew G Knepley   DMGetPointSF - Get the PetscSF encoding the parallel section point overlap for the DM.
3258b21d0597SMatthew G Knepley 
3259b21d0597SMatthew G Knepley   Input Parameter:
3260b21d0597SMatthew G Knepley . dm - The DM
3261b21d0597SMatthew G Knepley 
3262b21d0597SMatthew G Knepley   Output Parameter:
3263b21d0597SMatthew G Knepley . sf - The PetscSF
3264b21d0597SMatthew G Knepley 
3265b21d0597SMatthew G Knepley   Level: intermediate
3266b21d0597SMatthew G Knepley 
3267b21d0597SMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
3268b21d0597SMatthew G Knepley 
3269057b4bcdSMatthew G Knepley .seealso: DMSetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3270b21d0597SMatthew G Knepley @*/
32710adebc6cSBarry Smith PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf)
32720adebc6cSBarry Smith {
3273b21d0597SMatthew G Knepley   PetscFunctionBegin;
3274b21d0597SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3275b21d0597SMatthew G Knepley   PetscValidPointer(sf, 2);
3276b21d0597SMatthew G Knepley   *sf = dm->sf;
3277b21d0597SMatthew G Knepley   PetscFunctionReturn(0);
3278b21d0597SMatthew G Knepley }
3279b21d0597SMatthew G Knepley 
3280b21d0597SMatthew G Knepley #undef __FUNCT__
3281057b4bcdSMatthew G Knepley #define __FUNCT__ "DMSetPointSF"
3282057b4bcdSMatthew G Knepley /*@
3283057b4bcdSMatthew G Knepley   DMSetPointSF - Set the PetscSF encoding the parallel section point overlap for the DM.
3284057b4bcdSMatthew G Knepley 
3285057b4bcdSMatthew G Knepley   Input Parameters:
3286057b4bcdSMatthew G Knepley + dm - The DM
3287057b4bcdSMatthew G Knepley - sf - The PetscSF
3288057b4bcdSMatthew G Knepley 
3289057b4bcdSMatthew G Knepley   Level: intermediate
3290057b4bcdSMatthew G Knepley 
3291057b4bcdSMatthew G Knepley .seealso: DMGetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3292057b4bcdSMatthew G Knepley @*/
32930adebc6cSBarry Smith PetscErrorCode DMSetPointSF(DM dm, PetscSF sf)
32940adebc6cSBarry Smith {
3295057b4bcdSMatthew G Knepley   PetscErrorCode ierr;
3296057b4bcdSMatthew G Knepley 
3297057b4bcdSMatthew G Knepley   PetscFunctionBegin;
3298057b4bcdSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3299057b4bcdSMatthew G Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
3300057b4bcdSMatthew G Knepley   ierr   = PetscSFDestroy(&dm->sf);CHKERRQ(ierr);
3301057b4bcdSMatthew G Knepley   ierr   = PetscObjectReference((PetscObject) sf);CHKERRQ(ierr);
3302057b4bcdSMatthew G Knepley   dm->sf = sf;
3303057b4bcdSMatthew G Knepley   PetscFunctionReturn(0);
3304057b4bcdSMatthew G Knepley }
3305057b4bcdSMatthew G Knepley 
3306057b4bcdSMatthew G Knepley #undef __FUNCT__
33072764a2aaSMatthew G. Knepley #define __FUNCT__ "DMGetDS"
33082764a2aaSMatthew G. Knepley /*@
33092764a2aaSMatthew G. Knepley   DMGetDS - Get the PetscDS
33102764a2aaSMatthew G. Knepley 
33112764a2aaSMatthew G. Knepley   Input Parameter:
33122764a2aaSMatthew G. Knepley . dm - The DM
33132764a2aaSMatthew G. Knepley 
33142764a2aaSMatthew G. Knepley   Output Parameter:
33152764a2aaSMatthew G. Knepley . prob - The PetscDS
33162764a2aaSMatthew G. Knepley 
33172764a2aaSMatthew G. Knepley   Level: developer
33182764a2aaSMatthew G. Knepley 
33192764a2aaSMatthew G. Knepley .seealso: DMSetDS()
33202764a2aaSMatthew G. Knepley @*/
33212764a2aaSMatthew G. Knepley PetscErrorCode DMGetDS(DM dm, PetscDS *prob)
3322af122d2aSMatthew G Knepley {
3323af122d2aSMatthew G Knepley   PetscFunctionBegin;
3324af122d2aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
33250f21e855SMatthew G. Knepley   PetscValidPointer(prob, 2);
33260f21e855SMatthew G. Knepley   *prob = dm->prob;
33270f21e855SMatthew G. Knepley   PetscFunctionReturn(0);
33280f21e855SMatthew G. Knepley }
33290f21e855SMatthew G. Knepley 
33300f21e855SMatthew G. Knepley #undef __FUNCT__
33312764a2aaSMatthew G. Knepley #define __FUNCT__ "DMSetDS"
33322764a2aaSMatthew G. Knepley /*@
33332764a2aaSMatthew G. Knepley   DMSetDS - Set the PetscDS
33342764a2aaSMatthew G. Knepley 
33352764a2aaSMatthew G. Knepley   Input Parameters:
33362764a2aaSMatthew G. Knepley + dm - The DM
33372764a2aaSMatthew G. Knepley - prob - The PetscDS
33382764a2aaSMatthew G. Knepley 
33392764a2aaSMatthew G. Knepley   Level: developer
33402764a2aaSMatthew G. Knepley 
33412764a2aaSMatthew G. Knepley .seealso: DMGetDS()
33422764a2aaSMatthew G. Knepley @*/
33432764a2aaSMatthew G. Knepley PetscErrorCode DMSetDS(DM dm, PetscDS prob)
33440f21e855SMatthew G. Knepley {
33450f21e855SMatthew G. Knepley   PetscErrorCode ierr;
33460f21e855SMatthew G. Knepley 
33470f21e855SMatthew G. Knepley   PetscFunctionBegin;
33480f21e855SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
33492764a2aaSMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2);
33502764a2aaSMatthew G. Knepley   ierr = PetscDSDestroy(&dm->prob);CHKERRQ(ierr);
33510f21e855SMatthew G. Knepley   dm->prob = prob;
33520f21e855SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) dm->prob);CHKERRQ(ierr);
33530f21e855SMatthew G. Knepley   PetscFunctionReturn(0);
33540f21e855SMatthew G. Knepley }
33550f21e855SMatthew G. Knepley 
33560f21e855SMatthew G. Knepley #undef __FUNCT__
33570f21e855SMatthew G. Knepley #define __FUNCT__ "DMGetNumFields"
33580f21e855SMatthew G. Knepley PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields)
33590f21e855SMatthew G. Knepley {
33600f21e855SMatthew G. Knepley   PetscErrorCode ierr;
33610f21e855SMatthew G. Knepley 
33620f21e855SMatthew G. Knepley   PetscFunctionBegin;
33630f21e855SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
33642764a2aaSMatthew G. Knepley   ierr = PetscDSGetNumFields(dm->prob, numFields);CHKERRQ(ierr);
3365af122d2aSMatthew G Knepley   PetscFunctionReturn(0);
3366af122d2aSMatthew G Knepley }
3367af122d2aSMatthew G Knepley 
3368af122d2aSMatthew G Knepley #undef __FUNCT__
3369af122d2aSMatthew G Knepley #define __FUNCT__ "DMSetNumFields"
3370af122d2aSMatthew G Knepley PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields)
3371af122d2aSMatthew G Knepley {
33720f21e855SMatthew G. Knepley   PetscInt       Nf, f;
3373af122d2aSMatthew G Knepley   PetscErrorCode ierr;
3374af122d2aSMatthew G Knepley 
3375af122d2aSMatthew G Knepley   PetscFunctionBegin;
3376af122d2aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
33772764a2aaSMatthew G. Knepley   ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
33780f21e855SMatthew G. Knepley   for (f = Nf; f < numFields; ++f) {
33790f21e855SMatthew G. Knepley     PetscContainer obj;
33800f21e855SMatthew G. Knepley 
33810f21e855SMatthew G. Knepley     ierr = PetscContainerCreate(PetscObjectComm((PetscObject) dm), &obj);CHKERRQ(ierr);
33822764a2aaSMatthew G. Knepley     ierr = PetscDSSetDiscretization(dm->prob, f, (PetscObject) obj);CHKERRQ(ierr);
33830f21e855SMatthew G. Knepley     ierr = PetscContainerDestroy(&obj);CHKERRQ(ierr);
3384af122d2aSMatthew G Knepley   }
3385af122d2aSMatthew G Knepley   PetscFunctionReturn(0);
3386af122d2aSMatthew G Knepley }
3387af122d2aSMatthew G Knepley 
3388af122d2aSMatthew G Knepley #undef __FUNCT__
3389af122d2aSMatthew G Knepley #define __FUNCT__ "DMGetField"
3390c1929be8SMatthew G. Knepley /*@
3391c1929be8SMatthew G. Knepley   DMGetField - Return the discretization object for a given DM field
3392c1929be8SMatthew G. Knepley 
3393c1929be8SMatthew G. Knepley   Not collective
3394c1929be8SMatthew G. Knepley 
3395c1929be8SMatthew G. Knepley   Input Parameters:
3396c1929be8SMatthew G. Knepley + dm - The DM
3397c1929be8SMatthew G. Knepley - f  - The field number
3398c1929be8SMatthew G. Knepley 
3399c1929be8SMatthew G. Knepley   Output Parameter:
3400c1929be8SMatthew G. Knepley . field - The discretization object
3401c1929be8SMatthew G. Knepley 
3402c1929be8SMatthew G. Knepley   Level: developer
3403c1929be8SMatthew G. Knepley 
3404c1929be8SMatthew G. Knepley .seealso: DMSetField()
3405c1929be8SMatthew G. Knepley @*/
3406af122d2aSMatthew G Knepley PetscErrorCode DMGetField(DM dm, PetscInt f, PetscObject *field)
3407af122d2aSMatthew G Knepley {
34080f21e855SMatthew G. Knepley   PetscErrorCode ierr;
34090f21e855SMatthew G. Knepley 
3410af122d2aSMatthew G Knepley   PetscFunctionBegin;
3411af122d2aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
34122764a2aaSMatthew G. Knepley   ierr = PetscDSGetDiscretization(dm->prob, f, field);CHKERRQ(ierr);
3413decb47aaSMatthew G. Knepley   PetscFunctionReturn(0);
3414decb47aaSMatthew G. Knepley }
3415decb47aaSMatthew G. Knepley 
3416decb47aaSMatthew G. Knepley #undef __FUNCT__
3417decb47aaSMatthew G. Knepley #define __FUNCT__ "DMSetField"
3418c1929be8SMatthew G. Knepley /*@
3419c1929be8SMatthew G. Knepley   DMSetField - Set the discretization object for a given DM field
3420c1929be8SMatthew G. Knepley 
3421c1929be8SMatthew G. Knepley   Logically collective on DM
3422c1929be8SMatthew G. Knepley 
3423c1929be8SMatthew G. Knepley   Input Parameters:
3424c1929be8SMatthew G. Knepley + dm - The DM
3425c1929be8SMatthew G. Knepley . f  - The field number
3426c1929be8SMatthew G. Knepley - field - The discretization object
3427c1929be8SMatthew G. Knepley 
3428c1929be8SMatthew G. Knepley   Level: developer
3429c1929be8SMatthew G. Knepley 
3430c1929be8SMatthew G. Knepley .seealso: DMGetField()
3431c1929be8SMatthew G. Knepley @*/
3432decb47aaSMatthew G. Knepley PetscErrorCode DMSetField(DM dm, PetscInt f, PetscObject field)
3433decb47aaSMatthew G. Knepley {
3434decb47aaSMatthew G. Knepley   PetscErrorCode ierr;
3435decb47aaSMatthew G. Knepley 
3436decb47aaSMatthew G. Knepley   PetscFunctionBegin;
3437decb47aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
34382764a2aaSMatthew G. Knepley   ierr = PetscDSSetDiscretization(dm->prob, f, field);CHKERRQ(ierr);
3439af122d2aSMatthew G Knepley   PetscFunctionReturn(0);
3440af122d2aSMatthew G Knepley }
34416636e97aSMatthew G Knepley 
34426636e97aSMatthew G Knepley #undef __FUNCT__
3443b64e0483SPeter Brune #define __FUNCT__ "DMRestrictHook_Coordinates"
3444b64e0483SPeter Brune PetscErrorCode DMRestrictHook_Coordinates(DM dm,DM dmc,void *ctx)
3445b64e0483SPeter Brune {
3446b64e0483SPeter Brune   DM dm_coord,dmc_coord;
3447b64e0483SPeter Brune   PetscErrorCode ierr;
3448b64e0483SPeter Brune   Vec coords,ccoords;
3449b64e0483SPeter Brune   VecScatter scat;
3450b64e0483SPeter Brune   PetscFunctionBegin;
3451b64e0483SPeter Brune   ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr);
3452b64e0483SPeter Brune   ierr = DMGetCoordinateDM(dmc,&dmc_coord);CHKERRQ(ierr);
3453b64e0483SPeter Brune   ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr);
3454b64e0483SPeter Brune   ierr = DMGetCoordinates(dmc,&ccoords);CHKERRQ(ierr);
3455b64e0483SPeter Brune   if (coords && !ccoords) {
3456b64e0483SPeter Brune     ierr = DMCreateGlobalVector(dmc_coord,&ccoords);CHKERRQ(ierr);
3457b64e0483SPeter Brune     ierr = DMCreateInjection(dmc_coord,dm_coord,&scat);CHKERRQ(ierr);
3458b64e0483SPeter Brune     ierr = VecScatterBegin(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3459b64e0483SPeter Brune     ierr = VecScatterEnd(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3460b64e0483SPeter Brune     ierr = DMSetCoordinates(dmc,ccoords);CHKERRQ(ierr);
3461b64e0483SPeter Brune     ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
3462b64e0483SPeter Brune     ierr = VecDestroy(&ccoords);CHKERRQ(ierr);
3463b64e0483SPeter Brune   }
3464b64e0483SPeter Brune   PetscFunctionReturn(0);
3465b64e0483SPeter Brune }
3466b64e0483SPeter Brune 
3467b64e0483SPeter Brune #undef __FUNCT__
346803dadc2fSPeter Brune #define __FUNCT__ "DMSubDomainHook_Coordinates"
346903dadc2fSPeter Brune static PetscErrorCode DMSubDomainHook_Coordinates(DM dm,DM subdm,void *ctx)
347003dadc2fSPeter Brune {
347103dadc2fSPeter Brune   DM dm_coord,subdm_coord;
347203dadc2fSPeter Brune   PetscErrorCode ierr;
347303dadc2fSPeter Brune   Vec coords,ccoords,clcoords;
347403dadc2fSPeter Brune   VecScatter *scat_i,*scat_g;
347503dadc2fSPeter Brune   PetscFunctionBegin;
347603dadc2fSPeter Brune   ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr);
347703dadc2fSPeter Brune   ierr = DMGetCoordinateDM(subdm,&subdm_coord);CHKERRQ(ierr);
347803dadc2fSPeter Brune   ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr);
347903dadc2fSPeter Brune   ierr = DMGetCoordinates(subdm,&ccoords);CHKERRQ(ierr);
348003dadc2fSPeter Brune   if (coords && !ccoords) {
348103dadc2fSPeter Brune     ierr = DMCreateGlobalVector(subdm_coord,&ccoords);CHKERRQ(ierr);
348203dadc2fSPeter Brune     ierr = DMCreateLocalVector(subdm_coord,&clcoords);CHKERRQ(ierr);
348303dadc2fSPeter Brune     ierr = DMCreateDomainDecompositionScatters(dm_coord,1,&subdm_coord,NULL,&scat_i,&scat_g);CHKERRQ(ierr);
348403dadc2fSPeter Brune     ierr = VecScatterBegin(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
348503dadc2fSPeter Brune     ierr = VecScatterBegin(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
348603dadc2fSPeter Brune     ierr = VecScatterEnd(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
348703dadc2fSPeter Brune     ierr = VecScatterEnd(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
348803dadc2fSPeter Brune     ierr = DMSetCoordinates(subdm,ccoords);CHKERRQ(ierr);
348903dadc2fSPeter Brune     ierr = DMSetCoordinatesLocal(subdm,clcoords);CHKERRQ(ierr);
349003dadc2fSPeter Brune     ierr = VecScatterDestroy(&scat_i[0]);CHKERRQ(ierr);
349103dadc2fSPeter Brune     ierr = VecScatterDestroy(&scat_g[0]);CHKERRQ(ierr);
349203dadc2fSPeter Brune     ierr = VecDestroy(&ccoords);CHKERRQ(ierr);
349303dadc2fSPeter Brune     ierr = VecDestroy(&clcoords);CHKERRQ(ierr);
349403dadc2fSPeter Brune     ierr = PetscFree(scat_i);CHKERRQ(ierr);
349503dadc2fSPeter Brune     ierr = PetscFree(scat_g);CHKERRQ(ierr);
349603dadc2fSPeter Brune   }
349703dadc2fSPeter Brune   PetscFunctionReturn(0);
349803dadc2fSPeter Brune }
349903dadc2fSPeter Brune 
350003dadc2fSPeter Brune #undef __FUNCT__
3501c73cfb54SMatthew G. Knepley #define __FUNCT__ "DMGetDimension"
3502c73cfb54SMatthew G. Knepley /*@
3503c73cfb54SMatthew G. Knepley   DMGetDimension - Return the topological dimension of the DM
3504c73cfb54SMatthew G. Knepley 
3505c73cfb54SMatthew G. Knepley   Not collective
3506c73cfb54SMatthew G. Knepley 
3507c73cfb54SMatthew G. Knepley   Input Parameter:
3508c73cfb54SMatthew G. Knepley . dm - The DM
3509c73cfb54SMatthew G. Knepley 
3510c73cfb54SMatthew G. Knepley   Output Parameter:
3511c73cfb54SMatthew G. Knepley . dim - The topological dimension
3512c73cfb54SMatthew G. Knepley 
3513c73cfb54SMatthew G. Knepley   Level: beginner
3514c73cfb54SMatthew G. Knepley 
3515c73cfb54SMatthew G. Knepley .seealso: DMSetDimension(), DMCreate()
3516c73cfb54SMatthew G. Knepley @*/
3517c73cfb54SMatthew G. Knepley PetscErrorCode DMGetDimension(DM dm, PetscInt *dim)
3518c73cfb54SMatthew G. Knepley {
3519c73cfb54SMatthew G. Knepley   PetscFunctionBegin;
3520c73cfb54SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3521c73cfb54SMatthew G. Knepley   PetscValidPointer(dim, 2);
3522c73cfb54SMatthew G. Knepley   *dim = dm->dim;
3523c73cfb54SMatthew G. Knepley   PetscFunctionReturn(0);
3524c73cfb54SMatthew G. Knepley }
3525c73cfb54SMatthew G. Knepley 
3526c73cfb54SMatthew G. Knepley #undef __FUNCT__
3527c73cfb54SMatthew G. Knepley #define __FUNCT__ "DMSetDimension"
3528c73cfb54SMatthew G. Knepley /*@
3529c73cfb54SMatthew G. Knepley   DMSetDimension - Set the topological dimension of the DM
3530c73cfb54SMatthew G. Knepley 
3531c73cfb54SMatthew G. Knepley   Collective on dm
3532c73cfb54SMatthew G. Knepley 
3533c73cfb54SMatthew G. Knepley   Input Parameters:
3534c73cfb54SMatthew G. Knepley + dm - The DM
3535c73cfb54SMatthew G. Knepley - dim - The topological dimension
3536c73cfb54SMatthew G. Knepley 
3537c73cfb54SMatthew G. Knepley   Level: beginner
3538c73cfb54SMatthew G. Knepley 
3539c73cfb54SMatthew G. Knepley .seealso: DMGetDimension(), DMCreate()
3540c73cfb54SMatthew G. Knepley @*/
3541c73cfb54SMatthew G. Knepley PetscErrorCode DMSetDimension(DM dm, PetscInt dim)
3542c73cfb54SMatthew G. Knepley {
3543c73cfb54SMatthew G. Knepley   PetscFunctionBegin;
3544c73cfb54SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3545c73cfb54SMatthew G. Knepley   PetscValidLogicalCollectiveInt(dm, dim, 2);
3546c73cfb54SMatthew G. Knepley   dm->dim = dim;
3547c73cfb54SMatthew G. Knepley   PetscFunctionReturn(0);
3548c73cfb54SMatthew G. Knepley }
3549c73cfb54SMatthew G. Knepley 
3550793f3fe5SMatthew G. Knepley #undef __FUNCT__
3551793f3fe5SMatthew G. Knepley #define __FUNCT__ "DMGetDimPoints"
3552793f3fe5SMatthew G. Knepley /*@
3553793f3fe5SMatthew G. Knepley   DMGetDimPoints - Get the half-open interval for all points of a given dimension
3554793f3fe5SMatthew G. Knepley 
3555793f3fe5SMatthew G. Knepley   Collective on DM
3556793f3fe5SMatthew G. Knepley 
3557793f3fe5SMatthew G. Knepley   Input Parameters:
3558793f3fe5SMatthew G. Knepley + dm - the DM
3559793f3fe5SMatthew G. Knepley - dim - the dimension
3560793f3fe5SMatthew G. Knepley 
3561793f3fe5SMatthew G. Knepley   Output Parameters:
3562793f3fe5SMatthew G. Knepley + pStart - The first point of the given dimension
3563793f3fe5SMatthew G. Knepley . pEnd - The first point following points of the given dimension
3564793f3fe5SMatthew G. Knepley 
3565793f3fe5SMatthew G. Knepley   Note:
3566793f3fe5SMatthew G. Knepley   The points are vertices in the Hasse diagram encoding the topology. This is explained in
3567793f3fe5SMatthew G. Knepley   http://arxiv.org/abs/0908.4427. If not points exist of this dimension in the storage scheme,
3568793f3fe5SMatthew G. Knepley   then the interval is empty.
3569793f3fe5SMatthew G. Knepley 
3570793f3fe5SMatthew G. Knepley   Level: intermediate
3571793f3fe5SMatthew G. Knepley 
3572793f3fe5SMatthew G. Knepley .keywords: point, Hasse Diagram, dimension
3573793f3fe5SMatthew G. Knepley .seealso: DMPLEX, DMPlexGetDepthStratum(), DMPlexGetHeightStratum()
3574793f3fe5SMatthew G. Knepley @*/
3575793f3fe5SMatthew G. Knepley PetscErrorCode DMGetDimPoints(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
3576793f3fe5SMatthew G. Knepley {
3577793f3fe5SMatthew G. Knepley   PetscInt       d;
3578793f3fe5SMatthew G. Knepley   PetscErrorCode ierr;
3579793f3fe5SMatthew G. Knepley 
3580793f3fe5SMatthew G. Knepley   PetscFunctionBegin;
3581793f3fe5SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3582793f3fe5SMatthew G. Knepley   ierr = DMGetDimension(dm, &d);CHKERRQ(ierr);
3583793f3fe5SMatthew G. Knepley   if ((dim < 0) || (dim > d)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d 1", dim, d);
3584793f3fe5SMatthew G. Knepley   ierr = (*dm->ops->getdimpoints)(dm, dim, pStart, pEnd);CHKERRQ(ierr);
3585793f3fe5SMatthew G. Knepley   PetscFunctionReturn(0);
3586793f3fe5SMatthew G. Knepley }
3587793f3fe5SMatthew G. Knepley 
3588793f3fe5SMatthew G. Knepley #undef __FUNCT__
35896636e97aSMatthew G Knepley #define __FUNCT__ "DMSetCoordinates"
35906636e97aSMatthew G Knepley /*@
35916636e97aSMatthew G Knepley   DMSetCoordinates - Sets into the DM a global vector that holds the coordinates
35926636e97aSMatthew G Knepley 
35936636e97aSMatthew G Knepley   Collective on DM
35946636e97aSMatthew G Knepley 
35956636e97aSMatthew G Knepley   Input Parameters:
35966636e97aSMatthew G Knepley + dm - the DM
35976636e97aSMatthew G Knepley - c - coordinate vector
35986636e97aSMatthew G Knepley 
35996636e97aSMatthew G Knepley   Note:
36006636e97aSMatthew G Knepley   The coordinates do include those for ghost points, which are in the local vector
36016636e97aSMatthew G Knepley 
36026636e97aSMatthew G Knepley   Level: intermediate
36036636e97aSMatthew G Knepley 
36046636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
36056636e97aSMatthew G Knepley .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLoca(), DMGetCoordinateDM()
36066636e97aSMatthew G Knepley @*/
36076636e97aSMatthew G Knepley PetscErrorCode DMSetCoordinates(DM dm, Vec c)
36086636e97aSMatthew G Knepley {
36096636e97aSMatthew G Knepley   PetscErrorCode ierr;
36106636e97aSMatthew G Knepley 
36116636e97aSMatthew G Knepley   PetscFunctionBegin;
36126636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36136636e97aSMatthew G Knepley   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
36146636e97aSMatthew G Knepley   ierr            = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
36156636e97aSMatthew G Knepley   ierr            = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
36166636e97aSMatthew G Knepley   dm->coordinates = c;
36176636e97aSMatthew G Knepley   ierr            = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
3618b64e0483SPeter Brune   ierr            = DMCoarsenHookAdd(dm,DMRestrictHook_Coordinates,NULL,NULL);CHKERRQ(ierr);
361903dadc2fSPeter Brune   ierr            = DMSubDomainHookAdd(dm,DMSubDomainHook_Coordinates,NULL,NULL);CHKERRQ(ierr);
36206636e97aSMatthew G Knepley   PetscFunctionReturn(0);
36216636e97aSMatthew G Knepley }
36226636e97aSMatthew G Knepley 
36236636e97aSMatthew G Knepley #undef __FUNCT__
36246636e97aSMatthew G Knepley #define __FUNCT__ "DMSetCoordinatesLocal"
36256636e97aSMatthew G Knepley /*@
36266636e97aSMatthew G Knepley   DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates
36276636e97aSMatthew G Knepley 
36286636e97aSMatthew G Knepley   Collective on DM
36296636e97aSMatthew G Knepley 
36306636e97aSMatthew G Knepley    Input Parameters:
36316636e97aSMatthew G Knepley +  dm - the DM
36326636e97aSMatthew G Knepley -  c - coordinate vector
36336636e97aSMatthew G Knepley 
36346636e97aSMatthew G Knepley   Note:
36356636e97aSMatthew G Knepley   The coordinates of ghost points can be set using DMSetCoordinates()
36366636e97aSMatthew G Knepley   followed by DMGetCoordinatesLocal(). This is intended to enable the
36376636e97aSMatthew G Knepley   setting of ghost coordinates outside of the domain.
36386636e97aSMatthew G Knepley 
36396636e97aSMatthew G Knepley   Level: intermediate
36406636e97aSMatthew G Knepley 
36416636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
36426636e97aSMatthew G Knepley .seealso: DMGetCoordinatesLocal(), DMSetCoordinates(), DMGetCoordinates(), DMGetCoordinateDM()
36436636e97aSMatthew G Knepley @*/
36446636e97aSMatthew G Knepley PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
36456636e97aSMatthew G Knepley {
36466636e97aSMatthew G Knepley   PetscErrorCode ierr;
36476636e97aSMatthew G Knepley 
36486636e97aSMatthew G Knepley   PetscFunctionBegin;
36496636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36506636e97aSMatthew G Knepley   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
36516636e97aSMatthew G Knepley   ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
36526636e97aSMatthew G Knepley   ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
36538865f1eaSKarl Rupp 
36546636e97aSMatthew G Knepley   dm->coordinatesLocal = c;
36558865f1eaSKarl Rupp 
36566636e97aSMatthew G Knepley   ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
36576636e97aSMatthew G Knepley   PetscFunctionReturn(0);
36586636e97aSMatthew G Knepley }
36596636e97aSMatthew G Knepley 
36606636e97aSMatthew G Knepley #undef __FUNCT__
36616636e97aSMatthew G Knepley #define __FUNCT__ "DMGetCoordinates"
36626636e97aSMatthew G Knepley /*@
36636636e97aSMatthew G Knepley   DMGetCoordinates - Gets a global vector with the coordinates associated with the DM.
36646636e97aSMatthew G Knepley 
36656636e97aSMatthew G Knepley   Not Collective
36666636e97aSMatthew G Knepley 
36676636e97aSMatthew G Knepley   Input Parameter:
36686636e97aSMatthew G Knepley . dm - the DM
36696636e97aSMatthew G Knepley 
36706636e97aSMatthew G Knepley   Output Parameter:
36716636e97aSMatthew G Knepley . c - global coordinate vector
36726636e97aSMatthew G Knepley 
36736636e97aSMatthew G Knepley   Note:
36746636e97aSMatthew G Knepley   This is a borrowed reference, so the user should NOT destroy this vector
36756636e97aSMatthew G Knepley 
36766636e97aSMatthew G Knepley   Each process has only the local coordinates (does NOT have the ghost coordinates).
36776636e97aSMatthew G Knepley 
36786636e97aSMatthew G Knepley   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
36796636e97aSMatthew G Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
36806636e97aSMatthew G Knepley 
36816636e97aSMatthew G Knepley   Level: intermediate
36826636e97aSMatthew G Knepley 
36836636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
36846636e97aSMatthew G Knepley .seealso: DMSetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM()
36856636e97aSMatthew G Knepley @*/
36866636e97aSMatthew G Knepley PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
36876636e97aSMatthew G Knepley {
36886636e97aSMatthew G Knepley   PetscErrorCode ierr;
36896636e97aSMatthew G Knepley 
36906636e97aSMatthew G Knepley   PetscFunctionBegin;
36916636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36926636e97aSMatthew G Knepley   PetscValidPointer(c,2);
36931f588964SMatthew G Knepley   if (!dm->coordinates && dm->coordinatesLocal) {
36940298fd71SBarry Smith     DM cdm = NULL;
36956636e97aSMatthew G Knepley 
36966636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
36976636e97aSMatthew G Knepley     ierr = DMCreateGlobalVector(cdm, &dm->coordinates);CHKERRQ(ierr);
36986636e97aSMatthew G Knepley     ierr = PetscObjectSetName((PetscObject) dm->coordinates, "coordinates");CHKERRQ(ierr);
36996636e97aSMatthew G Knepley     ierr = DMLocalToGlobalBegin(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
37006636e97aSMatthew G Knepley     ierr = DMLocalToGlobalEnd(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
37016636e97aSMatthew G Knepley   }
37026636e97aSMatthew G Knepley   *c = dm->coordinates;
37036636e97aSMatthew G Knepley   PetscFunctionReturn(0);
37046636e97aSMatthew G Knepley }
37056636e97aSMatthew G Knepley 
37066636e97aSMatthew G Knepley #undef __FUNCT__
37076636e97aSMatthew G Knepley #define __FUNCT__ "DMGetCoordinatesLocal"
37086636e97aSMatthew G Knepley /*@
37096636e97aSMatthew G Knepley   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM.
37106636e97aSMatthew G Knepley 
37116636e97aSMatthew G Knepley   Collective on DM
37126636e97aSMatthew G Knepley 
37136636e97aSMatthew G Knepley   Input Parameter:
37146636e97aSMatthew G Knepley . dm - the DM
37156636e97aSMatthew G Knepley 
37166636e97aSMatthew G Knepley   Output Parameter:
37176636e97aSMatthew G Knepley . c - coordinate vector
37186636e97aSMatthew G Knepley 
37196636e97aSMatthew G Knepley   Note:
37206636e97aSMatthew G Knepley   This is a borrowed reference, so the user should NOT destroy this vector
37216636e97aSMatthew G Knepley 
37226636e97aSMatthew G Knepley   Each process has the local and ghost coordinates
37236636e97aSMatthew G Knepley 
37246636e97aSMatthew G Knepley   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
37256636e97aSMatthew G Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
37266636e97aSMatthew G Knepley 
37276636e97aSMatthew G Knepley   Level: intermediate
37286636e97aSMatthew G Knepley 
37296636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
37306636e97aSMatthew G Knepley .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMSetCoordinates(), DMGetCoordinateDM()
37316636e97aSMatthew G Knepley @*/
37326636e97aSMatthew G Knepley PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
37336636e97aSMatthew G Knepley {
37346636e97aSMatthew G Knepley   PetscErrorCode ierr;
37356636e97aSMatthew G Knepley 
37366636e97aSMatthew G Knepley   PetscFunctionBegin;
37376636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37386636e97aSMatthew G Knepley   PetscValidPointer(c,2);
37391f588964SMatthew G Knepley   if (!dm->coordinatesLocal && dm->coordinates) {
37400298fd71SBarry Smith     DM cdm = NULL;
37416636e97aSMatthew G Knepley 
37426636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
37436636e97aSMatthew G Knepley     ierr = DMCreateLocalVector(cdm, &dm->coordinatesLocal);CHKERRQ(ierr);
37446636e97aSMatthew G Knepley     ierr = PetscObjectSetName((PetscObject) dm->coordinatesLocal, "coordinates");CHKERRQ(ierr);
37456636e97aSMatthew G Knepley     ierr = DMGlobalToLocalBegin(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
37466636e97aSMatthew G Knepley     ierr = DMGlobalToLocalEnd(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
37476636e97aSMatthew G Knepley   }
37486636e97aSMatthew G Knepley   *c = dm->coordinatesLocal;
37496636e97aSMatthew G Knepley   PetscFunctionReturn(0);
37506636e97aSMatthew G Knepley }
37516636e97aSMatthew G Knepley 
37526636e97aSMatthew G Knepley #undef __FUNCT__
37536636e97aSMatthew G Knepley #define __FUNCT__ "DMGetCoordinateDM"
37546636e97aSMatthew G Knepley /*@
37551cfe2091SMatthew G. Knepley   DMGetCoordinateDM - Gets the DM that prescribes coordinate layout and scatters between global and local coordinates
37566636e97aSMatthew G Knepley 
37576636e97aSMatthew G Knepley   Collective on DM
37586636e97aSMatthew G Knepley 
37596636e97aSMatthew G Knepley   Input Parameter:
37606636e97aSMatthew G Knepley . dm - the DM
37616636e97aSMatthew G Knepley 
37626636e97aSMatthew G Knepley   Output Parameter:
37636636e97aSMatthew G Knepley . cdm - coordinate DM
37646636e97aSMatthew G Knepley 
37656636e97aSMatthew G Knepley   Level: intermediate
37666636e97aSMatthew G Knepley 
37676636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
37681cfe2091SMatthew G. Knepley .seealso: DMSetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
37696636e97aSMatthew G Knepley @*/
37706636e97aSMatthew G Knepley PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
37716636e97aSMatthew G Knepley {
37726636e97aSMatthew G Knepley   PetscErrorCode ierr;
37736636e97aSMatthew G Knepley 
37746636e97aSMatthew G Knepley   PetscFunctionBegin;
37756636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37766636e97aSMatthew G Knepley   PetscValidPointer(cdm,2);
37776636e97aSMatthew G Knepley   if (!dm->coordinateDM) {
377882f516ccSBarry Smith     if (!dm->ops->createcoordinatedm) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unable to create coordinates for this DM");
37796636e97aSMatthew G Knepley     ierr = (*dm->ops->createcoordinatedm)(dm, &dm->coordinateDM);CHKERRQ(ierr);
37806636e97aSMatthew G Knepley   }
37816636e97aSMatthew G Knepley   *cdm = dm->coordinateDM;
37826636e97aSMatthew G Knepley   PetscFunctionReturn(0);
37836636e97aSMatthew G Knepley }
3784e87bb0d3SMatthew G Knepley 
3785e87bb0d3SMatthew G Knepley #undef __FUNCT__
37861cfe2091SMatthew G. Knepley #define __FUNCT__ "DMSetCoordinateDM"
37871cfe2091SMatthew G. Knepley /*@
37881cfe2091SMatthew G. Knepley   DMSetCoordinateDM - Sets the DM that prescribes coordinate layout and scatters between global and local coordinates
37891cfe2091SMatthew G. Knepley 
37901cfe2091SMatthew G. Knepley   Logically Collective on DM
37911cfe2091SMatthew G. Knepley 
37921cfe2091SMatthew G. Knepley   Input Parameters:
37931cfe2091SMatthew G. Knepley + dm - the DM
37941cfe2091SMatthew G. Knepley - cdm - coordinate DM
37951cfe2091SMatthew G. Knepley 
37961cfe2091SMatthew G. Knepley   Level: intermediate
37971cfe2091SMatthew G. Knepley 
37981cfe2091SMatthew G. Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
37991cfe2091SMatthew G. Knepley .seealso: DMGetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
38001cfe2091SMatthew G. Knepley @*/
38011cfe2091SMatthew G. Knepley PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm)
38021cfe2091SMatthew G. Knepley {
38031cfe2091SMatthew G. Knepley   PetscErrorCode ierr;
38041cfe2091SMatthew G. Knepley 
38051cfe2091SMatthew G. Knepley   PetscFunctionBegin;
38061cfe2091SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
38071cfe2091SMatthew G. Knepley   PetscValidHeaderSpecific(cdm,DM_CLASSID,2);
38081cfe2091SMatthew G. Knepley   ierr = DMDestroy(&dm->coordinateDM);CHKERRQ(ierr);
38091cfe2091SMatthew G. Knepley   dm->coordinateDM = cdm;
38101cfe2091SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) dm->coordinateDM);CHKERRQ(ierr);
38111cfe2091SMatthew G. Knepley   PetscFunctionReturn(0);
38121cfe2091SMatthew G. Knepley }
38131cfe2091SMatthew G. Knepley 
38141cfe2091SMatthew G. Knepley #undef __FUNCT__
381546e270d4SMatthew G. Knepley #define __FUNCT__ "DMGetCoordinateDim"
381646e270d4SMatthew G. Knepley /*@
381746e270d4SMatthew G. Knepley   DMGetCoordinateDim - Retrieve the dimension of embedding space for coordinate values.
381846e270d4SMatthew G. Knepley 
381946e270d4SMatthew G. Knepley   Not Collective
382046e270d4SMatthew G. Knepley 
382146e270d4SMatthew G. Knepley   Input Parameter:
382246e270d4SMatthew G. Knepley . dm - The DM object
382346e270d4SMatthew G. Knepley 
382446e270d4SMatthew G. Knepley   Output Parameter:
382546e270d4SMatthew G. Knepley . dim - The embedding dimension
382646e270d4SMatthew G. Knepley 
382746e270d4SMatthew G. Knepley   Level: intermediate
382846e270d4SMatthew G. Knepley 
382946e270d4SMatthew G. Knepley .keywords: mesh, coordinates
383046e270d4SMatthew G. Knepley .seealso: DMSetCoordinateDim(), DMGetCoordinateSection(), DMGetCoordinateDM(), DMGetDefaultSection(), DMSetDefaultSection()
383146e270d4SMatthew G. Knepley @*/
383246e270d4SMatthew G. Knepley PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim)
383346e270d4SMatthew G. Knepley {
383446e270d4SMatthew G. Knepley   PetscFunctionBegin;
383546e270d4SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
383646e270d4SMatthew G. Knepley   PetscValidPointer(dim, 2);
38379a9a41abSToby Isaac   if (dm->dimEmbed == PETSC_DEFAULT) {
38389a9a41abSToby Isaac     dm->dimEmbed = dm->dim;
38399a9a41abSToby Isaac   }
384046e270d4SMatthew G. Knepley   *dim = dm->dimEmbed;
384146e270d4SMatthew G. Knepley   PetscFunctionReturn(0);
384246e270d4SMatthew G. Knepley }
384346e270d4SMatthew G. Knepley 
384446e270d4SMatthew G. Knepley #undef __FUNCT__
384546e270d4SMatthew G. Knepley #define __FUNCT__ "DMSetCoordinateDim"
384646e270d4SMatthew G. Knepley /*@
384746e270d4SMatthew G. Knepley   DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values.
384846e270d4SMatthew G. Knepley 
384946e270d4SMatthew G. Knepley   Not Collective
385046e270d4SMatthew G. Knepley 
385146e270d4SMatthew G. Knepley   Input Parameters:
385246e270d4SMatthew G. Knepley + dm  - The DM object
385346e270d4SMatthew G. Knepley - dim - The embedding dimension
385446e270d4SMatthew G. Knepley 
385546e270d4SMatthew G. Knepley   Level: intermediate
385646e270d4SMatthew G. Knepley 
385746e270d4SMatthew G. Knepley .keywords: mesh, coordinates
385846e270d4SMatthew G. Knepley .seealso: DMGetCoordinateDim(), DMSetCoordinateSection(), DMGetCoordinateSection(), DMGetDefaultSection(), DMSetDefaultSection()
385946e270d4SMatthew G. Knepley @*/
386046e270d4SMatthew G. Knepley PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim)
386146e270d4SMatthew G. Knepley {
386246e270d4SMatthew G. Knepley   PetscFunctionBegin;
386346e270d4SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
386446e270d4SMatthew G. Knepley   dm->dimEmbed = dim;
386546e270d4SMatthew G. Knepley   PetscFunctionReturn(0);
386646e270d4SMatthew G. Knepley }
386746e270d4SMatthew G. Knepley 
386846e270d4SMatthew G. Knepley #undef __FUNCT__
3869e8abe2deSMatthew G. Knepley #define __FUNCT__ "DMGetCoordinateSection"
3870e8abe2deSMatthew G. Knepley /*@
3871e8abe2deSMatthew G. Knepley   DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh.
3872e8abe2deSMatthew G. Knepley 
3873e8abe2deSMatthew G. Knepley   Not Collective
3874e8abe2deSMatthew G. Knepley 
3875e8abe2deSMatthew G. Knepley   Input Parameter:
3876e8abe2deSMatthew G. Knepley . dm - The DM object
3877e8abe2deSMatthew G. Knepley 
3878e8abe2deSMatthew G. Knepley   Output Parameter:
3879e8abe2deSMatthew G. Knepley . section - The PetscSection object
3880e8abe2deSMatthew G. Knepley 
3881e8abe2deSMatthew G. Knepley   Level: intermediate
3882e8abe2deSMatthew G. Knepley 
3883e8abe2deSMatthew G. Knepley .keywords: mesh, coordinates
3884e8abe2deSMatthew G. Knepley .seealso: DMGetCoordinateDM(), DMGetDefaultSection(), DMSetDefaultSection()
3885e8abe2deSMatthew G. Knepley @*/
3886e8abe2deSMatthew G. Knepley PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section)
3887e8abe2deSMatthew G. Knepley {
3888e8abe2deSMatthew G. Knepley   DM             cdm;
3889e8abe2deSMatthew G. Knepley   PetscErrorCode ierr;
3890e8abe2deSMatthew G. Knepley 
3891e8abe2deSMatthew G. Knepley   PetscFunctionBegin;
3892e8abe2deSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3893e8abe2deSMatthew G. Knepley   PetscValidPointer(section, 2);
3894e8abe2deSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3895e8abe2deSMatthew G. Knepley   ierr = DMGetDefaultSection(cdm, section);CHKERRQ(ierr);
3896e8abe2deSMatthew G. Knepley   PetscFunctionReturn(0);
3897e8abe2deSMatthew G. Knepley }
3898e8abe2deSMatthew G. Knepley 
3899e8abe2deSMatthew G. Knepley #undef __FUNCT__
3900e8abe2deSMatthew G. Knepley #define __FUNCT__ "DMSetCoordinateSection"
3901e8abe2deSMatthew G. Knepley /*@
3902e8abe2deSMatthew G. Knepley   DMSetCoordinateSection - Set the layout of coordinate values over the mesh.
3903e8abe2deSMatthew G. Knepley 
3904e8abe2deSMatthew G. Knepley   Not Collective
3905e8abe2deSMatthew G. Knepley 
3906e8abe2deSMatthew G. Knepley   Input Parameters:
3907e8abe2deSMatthew G. Knepley + dm      - The DM object
390846e270d4SMatthew G. Knepley . dim     - The embedding dimension, or PETSC_DETERMINE
3909e8abe2deSMatthew G. Knepley - section - The PetscSection object
3910e8abe2deSMatthew G. Knepley 
3911e8abe2deSMatthew G. Knepley   Level: intermediate
3912e8abe2deSMatthew G. Knepley 
3913e8abe2deSMatthew G. Knepley .keywords: mesh, coordinates
3914e8abe2deSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetDefaultSection(), DMSetDefaultSection()
3915e8abe2deSMatthew G. Knepley @*/
391646e270d4SMatthew G. Knepley PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section)
3917e8abe2deSMatthew G. Knepley {
3918e8abe2deSMatthew G. Knepley   DM             cdm;
3919e8abe2deSMatthew G. Knepley   PetscErrorCode ierr;
3920e8abe2deSMatthew G. Knepley 
3921e8abe2deSMatthew G. Knepley   PetscFunctionBegin;
3922e8abe2deSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
392346e270d4SMatthew G. Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,3);
3924e8abe2deSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3925e8abe2deSMatthew G. Knepley   ierr = DMSetDefaultSection(cdm, section);CHKERRQ(ierr);
392646e270d4SMatthew G. Knepley   if (dim == PETSC_DETERMINE) {
392746e270d4SMatthew G. Knepley     PetscInt d = dim;
392846e270d4SMatthew G. Knepley     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
392946e270d4SMatthew G. Knepley 
393046e270d4SMatthew G. Knepley     ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
393146e270d4SMatthew G. Knepley     ierr = DMGetDimPoints(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
393246e270d4SMatthew G. Knepley     pStart = PetscMax(vStart, pStart);
393346e270d4SMatthew G. Knepley     pEnd   = PetscMin(vEnd, pEnd);
393446e270d4SMatthew G. Knepley     for (v = pStart; v < pEnd; ++v) {
393546e270d4SMatthew G. Knepley       ierr = PetscSectionGetDof(section, v, &dd);CHKERRQ(ierr);
393646e270d4SMatthew G. Knepley       if (dd) {d = dd; break;}
393746e270d4SMatthew G. Knepley     }
393846e270d4SMatthew G. Knepley     ierr = DMSetCoordinateDim(dm, d);CHKERRQ(ierr);
393946e270d4SMatthew G. Knepley   }
3940e8abe2deSMatthew G. Knepley   PetscFunctionReturn(0);
3941e8abe2deSMatthew G. Knepley }
3942e8abe2deSMatthew G. Knepley 
3943e8abe2deSMatthew G. Knepley #undef __FUNCT__
3944c6b900c6SMatthew G. Knepley #define __FUNCT__ "DMGetPeriodicity"
3945c6b900c6SMatthew G. Knepley PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal **maxCell, const PetscReal **L)
3946c6b900c6SMatthew G. Knepley {
3947c6b900c6SMatthew G. Knepley   PetscFunctionBegin;
3948c6b900c6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3949c6b900c6SMatthew G. Knepley   if (L)       *L       = dm->L;
3950c6b900c6SMatthew G. Knepley   if (maxCell) *maxCell = dm->maxCell;
3951c6b900c6SMatthew G. Knepley   PetscFunctionReturn(0);
3952c6b900c6SMatthew G. Knepley }
3953c6b900c6SMatthew G. Knepley 
3954c6b900c6SMatthew G. Knepley #undef __FUNCT__
3955c6b900c6SMatthew G. Knepley #define __FUNCT__ "DMSetPeriodicity"
3956c6b900c6SMatthew G. Knepley PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal L[])
3957c6b900c6SMatthew G. Knepley {
3958c6b900c6SMatthew G. Knepley   Vec            coordinates;
3959c6b900c6SMatthew G. Knepley   PetscInt       dim, d;
3960c6b900c6SMatthew G. Knepley   PetscErrorCode ierr;
3961c6b900c6SMatthew G. Knepley 
3962c6b900c6SMatthew G. Knepley   PetscFunctionBegin;
3963c6b900c6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3964c6b900c6SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
3965c6b900c6SMatthew G. Knepley   ierr = VecGetBlockSize(coordinates, &dim);CHKERRQ(ierr);
3966c6b900c6SMatthew G. Knepley   ierr = PetscMalloc1(dim,&dm->L);CHKERRQ(ierr);
3967c6b900c6SMatthew G. Knepley   ierr = PetscMalloc1(dim,&dm->maxCell);CHKERRQ(ierr);
3968c6b900c6SMatthew G. Knepley   for (d = 0; d < dim; ++d) {dm->L[d] = L[d]; dm->maxCell[d] = maxCell[d];}
3969c6b900c6SMatthew G. Knepley   PetscFunctionReturn(0);
3970c6b900c6SMatthew G. Knepley }
3971c6b900c6SMatthew G. Knepley 
3972c6b900c6SMatthew G. Knepley #undef __FUNCT__
3973e87bb0d3SMatthew G Knepley #define __FUNCT__ "DMLocatePoints"
3974e87bb0d3SMatthew G Knepley /*@
3975e87bb0d3SMatthew G Knepley   DMLocatePoints - Locate the points in v in the mesh and return an IS of the containing cells
3976e87bb0d3SMatthew G Knepley 
3977e87bb0d3SMatthew G Knepley   Not collective
3978e87bb0d3SMatthew G Knepley 
3979e87bb0d3SMatthew G Knepley   Input Parameters:
3980e87bb0d3SMatthew G Knepley + dm - The DM
3981e87bb0d3SMatthew G Knepley - v - The Vec of points
3982e87bb0d3SMatthew G Knepley 
398361e3bb9bSMatthew G Knepley   Output Parameter:
3984e87bb0d3SMatthew G Knepley . cells - The local cell numbers for cells which contain the points
3985e87bb0d3SMatthew G Knepley 
3986e87bb0d3SMatthew G Knepley   Level: developer
398761e3bb9bSMatthew G Knepley 
398861e3bb9bSMatthew G Knepley .keywords: point location, mesh
398961e3bb9bSMatthew G Knepley .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
399061e3bb9bSMatthew G Knepley @*/
3991e87bb0d3SMatthew G Knepley PetscErrorCode DMLocatePoints(DM dm, Vec v, IS *cells)
3992e87bb0d3SMatthew G Knepley {
3993735aa83eSMatthew G Knepley   PetscErrorCode ierr;
3994735aa83eSMatthew G Knepley 
3995e87bb0d3SMatthew G Knepley   PetscFunctionBegin;
3996e87bb0d3SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3997e87bb0d3SMatthew G Knepley   PetscValidHeaderSpecific(v,VEC_CLASSID,2);
3998e87bb0d3SMatthew G Knepley   PetscValidPointer(cells,3);
3999735aa83eSMatthew G Knepley   if (dm->ops->locatepoints) {
4000735aa83eSMatthew G Knepley     ierr = (*dm->ops->locatepoints)(dm,v,cells);CHKERRQ(ierr);
400182f516ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Point location not available for this DM");
4002e87bb0d3SMatthew G Knepley   PetscFunctionReturn(0);
4003e87bb0d3SMatthew G Knepley }
400414f150ffSMatthew G. Knepley 
400514f150ffSMatthew G. Knepley #undef __FUNCT__
400614f150ffSMatthew G. Knepley #define __FUNCT__ "DMGetOutputDM"
4007f4d763aaSMatthew G. Knepley /*@
4008f4d763aaSMatthew G. Knepley   DMGetOutputDM - Retrieve the DM associated with the layout for output
4009f4d763aaSMatthew G. Knepley 
4010f4d763aaSMatthew G. Knepley   Input Parameter:
4011f4d763aaSMatthew G. Knepley . dm - The original DM
4012f4d763aaSMatthew G. Knepley 
4013f4d763aaSMatthew G. Knepley   Output Parameter:
4014f4d763aaSMatthew G. Knepley . odm - The DM which provides the layout for output
4015f4d763aaSMatthew G. Knepley 
4016f4d763aaSMatthew G. Knepley   Level: intermediate
4017f4d763aaSMatthew G. Knepley 
4018f4d763aaSMatthew G. Knepley .seealso: VecView(), DMGetDefaultGlobalSection()
4019f4d763aaSMatthew G. Knepley @*/
402014f150ffSMatthew G. Knepley PetscErrorCode DMGetOutputDM(DM dm, DM *odm)
402114f150ffSMatthew G. Knepley {
4022c26acbdeSMatthew G. Knepley   PetscSection   section;
4023c26acbdeSMatthew G. Knepley   PetscBool      hasConstraints;
402414f150ffSMatthew G. Knepley   PetscErrorCode ierr;
402514f150ffSMatthew G. Knepley 
402614f150ffSMatthew G. Knepley   PetscFunctionBegin;
402714f150ffSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
402814f150ffSMatthew G. Knepley   PetscValidPointer(odm,2);
4029c26acbdeSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
4030c26acbdeSMatthew G. Knepley   ierr = PetscSectionHasConstraints(section, &hasConstraints);CHKERRQ(ierr);
4031c26acbdeSMatthew G. Knepley   if (!hasConstraints) {
4032c26acbdeSMatthew G. Knepley     *odm = dm;
4033c26acbdeSMatthew G. Knepley     PetscFunctionReturn(0);
4034c26acbdeSMatthew G. Knepley   }
403514f150ffSMatthew G. Knepley   if (!dm->dmBC) {
4036c26acbdeSMatthew G. Knepley     PetscSection newSection, gsection;
403714f150ffSMatthew G. Knepley     PetscSF      sf;
403814f150ffSMatthew G. Knepley 
403914f150ffSMatthew G. Knepley     ierr = DMClone(dm, &dm->dmBC);CHKERRQ(ierr);
404014f150ffSMatthew G. Knepley     ierr = PetscSectionClone(section, &newSection);CHKERRQ(ierr);
404114f150ffSMatthew G. Knepley     ierr = DMSetDefaultSection(dm->dmBC, newSection);CHKERRQ(ierr);
404214f150ffSMatthew G. Knepley     ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);
404314f150ffSMatthew G. Knepley     ierr = DMGetPointSF(dm->dmBC, &sf);CHKERRQ(ierr);
404414f150ffSMatthew G. Knepley     ierr = PetscSectionCreateGlobalSection(section, sf, PETSC_TRUE, &gsection);CHKERRQ(ierr);
404514f150ffSMatthew G. Knepley     ierr = DMSetDefaultGlobalSection(dm->dmBC, gsection);CHKERRQ(ierr);
404614f150ffSMatthew G. Knepley     ierr = PetscSectionDestroy(&gsection);CHKERRQ(ierr);
404714f150ffSMatthew G. Knepley   }
404814f150ffSMatthew G. Knepley   *odm = dm->dmBC;
404914f150ffSMatthew G. Knepley   PetscFunctionReturn(0);
405014f150ffSMatthew G. Knepley }
4051f4d763aaSMatthew G. Knepley 
4052f4d763aaSMatthew G. Knepley #undef __FUNCT__
4053f4d763aaSMatthew G. Knepley #define __FUNCT__ "DMGetOutputSequenceNumber"
4054f4d763aaSMatthew G. Knepley /*@
4055cdb7a50dSMatthew G. Knepley   DMGetOutputSequenceNumber - Retrieve the sequence number/value for output
4056f4d763aaSMatthew G. Knepley 
4057f4d763aaSMatthew G. Knepley   Input Parameter:
4058f4d763aaSMatthew G. Knepley . dm - The original DM
4059f4d763aaSMatthew G. Knepley 
4060cdb7a50dSMatthew G. Knepley   Output Parameters:
4061cdb7a50dSMatthew G. Knepley + num - The output sequence number
4062cdb7a50dSMatthew G. Knepley - val - The output sequence value
4063f4d763aaSMatthew G. Knepley 
4064f4d763aaSMatthew G. Knepley   Level: intermediate
4065f4d763aaSMatthew G. Knepley 
4066f4d763aaSMatthew G. Knepley   Note: This is intended for output that should appear in sequence, for instance
4067f4d763aaSMatthew G. Knepley   a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system.
4068f4d763aaSMatthew G. Knepley 
4069f4d763aaSMatthew G. Knepley .seealso: VecView()
4070f4d763aaSMatthew G. Knepley @*/
4071cdb7a50dSMatthew G. Knepley PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val)
4072f4d763aaSMatthew G. Knepley {
4073f4d763aaSMatthew G. Knepley   PetscFunctionBegin;
4074f4d763aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4075cdb7a50dSMatthew G. Knepley   if (num) {PetscValidPointer(num,2); *num = dm->outputSequenceNum;}
4076cdb7a50dSMatthew G. Knepley   if (val) {PetscValidPointer(val,3);*val = dm->outputSequenceVal;}
4077f4d763aaSMatthew G. Knepley   PetscFunctionReturn(0);
4078f4d763aaSMatthew G. Knepley }
4079f4d763aaSMatthew G. Knepley 
4080f4d763aaSMatthew G. Knepley #undef __FUNCT__
4081f4d763aaSMatthew G. Knepley #define __FUNCT__ "DMSetOutputSequenceNumber"
4082f4d763aaSMatthew G. Knepley /*@
4083cdb7a50dSMatthew G. Knepley   DMSetOutputSequenceNumber - Set the sequence number/value for output
4084f4d763aaSMatthew G. Knepley 
4085f4d763aaSMatthew G. Knepley   Input Parameters:
4086f4d763aaSMatthew G. Knepley + dm - The original DM
4087cdb7a50dSMatthew G. Knepley . num - The output sequence number
4088cdb7a50dSMatthew G. Knepley - val - The output sequence value
4089f4d763aaSMatthew G. Knepley 
4090f4d763aaSMatthew G. Knepley   Level: intermediate
4091f4d763aaSMatthew G. Knepley 
4092f4d763aaSMatthew G. Knepley   Note: This is intended for output that should appear in sequence, for instance
4093f4d763aaSMatthew G. Knepley   a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system.
4094f4d763aaSMatthew G. Knepley 
4095f4d763aaSMatthew G. Knepley .seealso: VecView()
4096f4d763aaSMatthew G. Knepley @*/
4097cdb7a50dSMatthew G. Knepley PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val)
4098f4d763aaSMatthew G. Knepley {
4099f4d763aaSMatthew G. Knepley   PetscFunctionBegin;
4100f4d763aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4101f4d763aaSMatthew G. Knepley   dm->outputSequenceNum = num;
4102cdb7a50dSMatthew G. Knepley   dm->outputSequenceVal = val;
4103cdb7a50dSMatthew G. Knepley   PetscFunctionReturn(0);
4104cdb7a50dSMatthew G. Knepley }
4105cdb7a50dSMatthew G. Knepley 
4106cdb7a50dSMatthew G. Knepley #undef __FUNCT__
4107cdb7a50dSMatthew G. Knepley #define __FUNCT__ "DMOutputSequenceLoad"
4108cdb7a50dSMatthew G. Knepley /*@C
4109cdb7a50dSMatthew G. Knepley   DMOutputSequenceLoad - Retrieve the sequence value from a Viewer
4110cdb7a50dSMatthew G. Knepley 
4111cdb7a50dSMatthew G. Knepley   Input Parameters:
4112cdb7a50dSMatthew G. Knepley + dm   - The original DM
4113cdb7a50dSMatthew G. Knepley . name - The sequence name
4114cdb7a50dSMatthew G. Knepley - num  - The output sequence number
4115cdb7a50dSMatthew G. Knepley 
4116cdb7a50dSMatthew G. Knepley   Output Parameter:
4117cdb7a50dSMatthew G. Knepley . val  - The output sequence value
4118cdb7a50dSMatthew G. Knepley 
4119cdb7a50dSMatthew G. Knepley   Level: intermediate
4120cdb7a50dSMatthew G. Knepley 
4121cdb7a50dSMatthew G. Knepley   Note: This is intended for output that should appear in sequence, for instance
4122cdb7a50dSMatthew G. Knepley   a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system.
4123cdb7a50dSMatthew G. Knepley 
4124cdb7a50dSMatthew G. Knepley .seealso: DMGetOutputSequenceNumber(), DMSetOutputSequenceNumber(), VecView()
4125cdb7a50dSMatthew G. Knepley @*/
4126cdb7a50dSMatthew G. Knepley PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char *name, PetscInt num, PetscReal *val)
4127cdb7a50dSMatthew G. Knepley {
4128cdb7a50dSMatthew G. Knepley   PetscBool      ishdf5;
4129cdb7a50dSMatthew G. Knepley   PetscErrorCode ierr;
4130cdb7a50dSMatthew G. Knepley 
4131cdb7a50dSMatthew G. Knepley   PetscFunctionBegin;
4132cdb7a50dSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4133cdb7a50dSMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
4134cdb7a50dSMatthew G. Knepley   PetscValidPointer(val,4);
4135cdb7a50dSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr);
4136cdb7a50dSMatthew G. Knepley   if (ishdf5) {
4137cdb7a50dSMatthew G. Knepley #if defined(PETSC_HAVE_HDF5)
4138cdb7a50dSMatthew G. Knepley     PetscScalar value;
4139cdb7a50dSMatthew G. Knepley 
4140cdb7a50dSMatthew G. Knepley     ierr = DMSequenceLoad_HDF5(dm, name, num, &value, viewer);CHKERRQ(ierr);
41414aeb217fSMatthew G. Knepley     *val = PetscRealPart(value);
4142cdb7a50dSMatthew G. Knepley #endif
4143cdb7a50dSMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
4144f4d763aaSMatthew G. Knepley   PetscFunctionReturn(0);
4145f4d763aaSMatthew G. Knepley }
4146