xref: /petsc/src/dm/interface/dm.c (revision fc31e74df696c5bb27cacc190f3f0be2538cc7ed)
1b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/
20c312b8eSJed Brown #include <petscsf.h>
347c6ae99SBarry Smith 
4732e2eb9SMatthew G Knepley PetscClassId  DM_CLASSID;
5f089877aSRichard Tran Mills PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocalToLocal;
667a56275SMatthew G Knepley 
7bff4a2f0SMatthew G. Knepley const char *const DMBoundaryTypes[] = {"NONE","GHOSTED","MIRROR","PERIODIC","TWIST","DM_BOUNDARY_",0};
8bff4a2f0SMatthew G. Knepley 
947c6ae99SBarry Smith #undef __FUNCT__
10a4121054SBarry Smith #define __FUNCT__ "DMCreate"
11a4121054SBarry Smith /*@
12de043629SMatthew G Knepley   DMCreate - Creates an empty DM object. The type can then be set with DMSetType().
13a4121054SBarry Smith 
14a4121054SBarry Smith    If you never  call DMSetType()  it will generate an
15a4121054SBarry Smith    error when you try to use the vector.
16a4121054SBarry Smith 
17a4121054SBarry Smith   Collective on MPI_Comm
18a4121054SBarry Smith 
19a4121054SBarry Smith   Input Parameter:
20a4121054SBarry Smith . comm - The communicator for the DM object
21a4121054SBarry Smith 
22a4121054SBarry Smith   Output Parameter:
23a4121054SBarry Smith . dm - The DM object
24a4121054SBarry Smith 
25a4121054SBarry Smith   Level: beginner
26a4121054SBarry Smith 
27a4121054SBarry Smith .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE
28a4121054SBarry Smith @*/
297087cfbeSBarry Smith PetscErrorCode  DMCreate(MPI_Comm comm,DM *dm)
30a4121054SBarry Smith {
31a4121054SBarry Smith   DM             v;
32a4121054SBarry Smith   PetscErrorCode ierr;
33a4121054SBarry Smith 
34a4121054SBarry Smith   PetscFunctionBegin;
351411c6eeSJed Brown   PetscValidPointer(dm,2);
360298fd71SBarry Smith   *dm = NULL;
378b984314SMatthew G. Knepley   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
38607a6623SBarry Smith   ierr = VecInitializePackage();CHKERRQ(ierr);
39607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
40607a6623SBarry Smith   ierr = DMInitializePackage();CHKERRQ(ierr);
41a4121054SBarry Smith 
4267c2884eSBarry Smith   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
43a4121054SBarry Smith   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
441411c6eeSJed Brown 
45e7c4fc90SDmitry Karpeev 
460298fd71SBarry Smith   v->ltogmap              = NULL;
471411c6eeSJed Brown   v->bs                   = 1;
48171400e9SBarry Smith   v->coloringtype         = IS_COLORING_GLOBAL;
4988ed4aceSMatthew G Knepley   ierr                    = PetscSFCreate(comm, &v->sf);CHKERRQ(ierr);
5088ed4aceSMatthew G Knepley   ierr                    = PetscSFCreate(comm, &v->defaultSF);CHKERRQ(ierr);
510298fd71SBarry Smith   v->defaultSection       = NULL;
520298fd71SBarry Smith   v->defaultGlobalSection = NULL;
53c6b900c6SMatthew G. Knepley   v->L                    = NULL;
54c6b900c6SMatthew G. Knepley   v->maxCell              = NULL;
55435a35e8SMatthew G Knepley   {
56435a35e8SMatthew G Knepley     PetscInt i;
57435a35e8SMatthew G Knepley     for (i = 0; i < 10; ++i) {
580298fd71SBarry Smith       v->nullspaceConstructors[i] = NULL;
59435a35e8SMatthew G Knepley     }
60435a35e8SMatthew G Knepley   }
61af122d2aSMatthew G Knepley   v->numFields = 0;
620298fd71SBarry Smith   v->fields    = NULL;
6314f150ffSMatthew G. Knepley   v->dmBC      = NULL;
64f4d763aaSMatthew G. Knepley   v->outputSequenceNum = -1;
65c0dedaeaSBarry Smith   ierr = DMSetVecType(v,VECSTANDARD);CHKERRQ(ierr);
66b412c318SBarry Smith   ierr = DMSetMatType(v,MATAIJ);CHKERRQ(ierr);
671411c6eeSJed Brown   *dm = v;
68a4121054SBarry Smith   PetscFunctionReturn(0);
69a4121054SBarry Smith }
70a4121054SBarry Smith 
71a4121054SBarry Smith #undef __FUNCT__
7238221697SMatthew G. Knepley #define __FUNCT__ "DMClone"
7338221697SMatthew G. Knepley /*@
7438221697SMatthew G. Knepley   DMClone - Creates a DM object with the same topology as the original.
7538221697SMatthew G. Knepley 
7638221697SMatthew G. Knepley   Collective on MPI_Comm
7738221697SMatthew G. Knepley 
7838221697SMatthew G. Knepley   Input Parameter:
7938221697SMatthew G. Knepley . dm - The original DM object
8038221697SMatthew G. Knepley 
8138221697SMatthew G. Knepley   Output Parameter:
8238221697SMatthew G. Knepley . newdm  - The new DM object
8338221697SMatthew G. Knepley 
8438221697SMatthew G. Knepley   Level: beginner
8538221697SMatthew G. Knepley 
8638221697SMatthew G. Knepley .keywords: DM, topology, create
8738221697SMatthew G. Knepley @*/
8838221697SMatthew G. Knepley PetscErrorCode DMClone(DM dm, DM *newdm)
8938221697SMatthew G. Knepley {
9038221697SMatthew G. Knepley   PetscSF        sf;
9138221697SMatthew G. Knepley   Vec            coords;
9238221697SMatthew G. Knepley   void          *ctx;
9338221697SMatthew G. Knepley   PetscErrorCode ierr;
9438221697SMatthew G. Knepley 
9538221697SMatthew G. Knepley   PetscFunctionBegin;
9638221697SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9738221697SMatthew G. Knepley   PetscValidPointer(newdm,2);
9838221697SMatthew G. Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), newdm);CHKERRQ(ierr);
9938221697SMatthew G. Knepley   if (dm->ops->clone) {
10038221697SMatthew G. Knepley     ierr = (*dm->ops->clone)(dm, newdm);CHKERRQ(ierr);
10138221697SMatthew G. Knepley   }
102cd27c7c9SMatthew G. Knepley   (*newdm)->setupcalled = PETSC_TRUE;
10338221697SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
10438221697SMatthew G. Knepley   ierr = DMSetPointSF(*newdm, sf);CHKERRQ(ierr);
10538221697SMatthew G. Knepley   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
10638221697SMatthew G. Knepley   ierr = DMSetApplicationContext(*newdm, ctx);CHKERRQ(ierr);
10738221697SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr);
10838221697SMatthew G. Knepley   if (coords) {
10938221697SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*newdm, coords);CHKERRQ(ierr);
11038221697SMatthew G. Knepley   } else {
11138221697SMatthew G. Knepley     ierr = DMGetCoordinates(dm, &coords);CHKERRQ(ierr);
11238221697SMatthew G. Knepley     if (coords) {ierr = DMSetCoordinates(*newdm, coords);CHKERRQ(ierr);}
11338221697SMatthew G. Knepley   }
114c6b900c6SMatthew G. Knepley   if (dm->maxCell) {
115c6b900c6SMatthew G. Knepley     const PetscReal *maxCell, *L;
116c6b900c6SMatthew G. Knepley     ierr = DMGetPeriodicity(dm,     &maxCell, &L);CHKERRQ(ierr);
117c6b900c6SMatthew G. Knepley     ierr = DMSetPeriodicity(*newdm,  maxCell,  L);CHKERRQ(ierr);
118c6b900c6SMatthew G. Knepley   }
11938221697SMatthew G. Knepley   PetscFunctionReturn(0);
12038221697SMatthew G. Knepley }
12138221697SMatthew G. Knepley 
12238221697SMatthew G. Knepley #undef __FUNCT__
1239a42bb27SBarry Smith #define __FUNCT__ "DMSetVecType"
1249a42bb27SBarry Smith /*@C
125564755cdSBarry Smith        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
1269a42bb27SBarry Smith 
127aa219208SBarry Smith    Logically Collective on DMDA
1289a42bb27SBarry Smith 
1299a42bb27SBarry Smith    Input Parameter:
1309a42bb27SBarry Smith +  da - initial distributed array
1318154be41SBarry Smith .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
1329a42bb27SBarry Smith 
1339a42bb27SBarry Smith    Options Database:
134dd85299cSBarry Smith .   -dm_vec_type ctype
1359a42bb27SBarry Smith 
1369a42bb27SBarry Smith    Level: intermediate
1379a42bb27SBarry Smith 
138c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType, DMGetVecType()
1399a42bb27SBarry Smith @*/
14019fd82e9SBarry Smith PetscErrorCode  DMSetVecType(DM da,VecType ctype)
1419a42bb27SBarry Smith {
1429a42bb27SBarry Smith   PetscErrorCode ierr;
1439a42bb27SBarry Smith 
1449a42bb27SBarry Smith   PetscFunctionBegin;
1459a42bb27SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1469a42bb27SBarry Smith   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
14719fd82e9SBarry Smith   ierr = PetscStrallocpy(ctype,(char**)&da->vectype);CHKERRQ(ierr);
1489a42bb27SBarry Smith   PetscFunctionReturn(0);
1499a42bb27SBarry Smith }
1509a42bb27SBarry Smith 
1519a42bb27SBarry Smith #undef __FUNCT__
152c0dedaeaSBarry Smith #define __FUNCT__ "DMGetVecType"
153c0dedaeaSBarry Smith /*@C
154c0dedaeaSBarry Smith        DMGetVecType - Gets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
155c0dedaeaSBarry Smith 
156c0dedaeaSBarry Smith    Logically Collective on DMDA
157c0dedaeaSBarry Smith 
158c0dedaeaSBarry Smith    Input Parameter:
159c0dedaeaSBarry Smith .  da - initial distributed array
160c0dedaeaSBarry Smith 
161c0dedaeaSBarry Smith    Output Parameter:
162c0dedaeaSBarry Smith .  ctype - the vector type
163c0dedaeaSBarry Smith 
164c0dedaeaSBarry Smith    Level: intermediate
165c0dedaeaSBarry Smith 
166c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
167c0dedaeaSBarry Smith @*/
168c0dedaeaSBarry Smith PetscErrorCode  DMGetVecType(DM da,VecType *ctype)
169c0dedaeaSBarry Smith {
170c0dedaeaSBarry Smith   PetscFunctionBegin;
171c0dedaeaSBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
172c0dedaeaSBarry Smith   *ctype = da->vectype;
173c0dedaeaSBarry Smith   PetscFunctionReturn(0);
174c0dedaeaSBarry Smith }
175c0dedaeaSBarry Smith 
176c0dedaeaSBarry Smith #undef __FUNCT__
1775f1ad066SMatthew G Knepley #define __FUNCT__ "VecGetDM"
1785f1ad066SMatthew G Knepley /*@
17934f98d34SBarry Smith   VecGetDM - Gets the DM defining the data layout of the vector
1805f1ad066SMatthew G Knepley 
1815f1ad066SMatthew G Knepley   Not collective
1825f1ad066SMatthew G Knepley 
1835f1ad066SMatthew G Knepley   Input Parameter:
1845f1ad066SMatthew G Knepley . v - The Vec
1855f1ad066SMatthew G Knepley 
1865f1ad066SMatthew G Knepley   Output Parameter:
1875f1ad066SMatthew G Knepley . dm - The DM
1885f1ad066SMatthew G Knepley 
1895f1ad066SMatthew G Knepley   Level: intermediate
1905f1ad066SMatthew G Knepley 
1915f1ad066SMatthew G Knepley .seealso: VecSetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType()
1925f1ad066SMatthew G Knepley @*/
1935f1ad066SMatthew G Knepley PetscErrorCode VecGetDM(Vec v, DM *dm)
1945f1ad066SMatthew G Knepley {
1955f1ad066SMatthew G Knepley   PetscErrorCode ierr;
1965f1ad066SMatthew G Knepley 
1975f1ad066SMatthew G Knepley   PetscFunctionBegin;
1985f1ad066SMatthew G Knepley   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
1995f1ad066SMatthew G Knepley   PetscValidPointer(dm,2);
2005f1ad066SMatthew G Knepley   ierr = PetscObjectQuery((PetscObject) v, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr);
2015f1ad066SMatthew G Knepley   PetscFunctionReturn(0);
2025f1ad066SMatthew G Knepley }
2035f1ad066SMatthew G Knepley 
2045f1ad066SMatthew G Knepley #undef __FUNCT__
2055f1ad066SMatthew G Knepley #define __FUNCT__ "VecSetDM"
2065f1ad066SMatthew G Knepley /*@
2075f1ad066SMatthew G Knepley   VecSetDM - Sets the DM defining the data layout of the vector
2085f1ad066SMatthew G Knepley 
2095f1ad066SMatthew G Knepley   Not collective
2105f1ad066SMatthew G Knepley 
2115f1ad066SMatthew G Knepley   Input Parameters:
2125f1ad066SMatthew G Knepley + v - The Vec
2135f1ad066SMatthew G Knepley - dm - The DM
2145f1ad066SMatthew G Knepley 
2155f1ad066SMatthew G Knepley   Level: intermediate
2165f1ad066SMatthew G Knepley 
2175f1ad066SMatthew G Knepley .seealso: VecGetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType()
2185f1ad066SMatthew G Knepley @*/
2195f1ad066SMatthew G Knepley PetscErrorCode VecSetDM(Vec v, DM dm)
2205f1ad066SMatthew G Knepley {
2215f1ad066SMatthew G Knepley   PetscErrorCode ierr;
2225f1ad066SMatthew G Knepley 
2235f1ad066SMatthew G Knepley   PetscFunctionBegin;
2245f1ad066SMatthew G Knepley   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
225d7f50e27SLisandro Dalcin   if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2);
2265f1ad066SMatthew G Knepley   ierr = PetscObjectCompose((PetscObject) v, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr);
2275f1ad066SMatthew G Knepley   PetscFunctionReturn(0);
2285f1ad066SMatthew G Knepley }
2295f1ad066SMatthew G Knepley 
2305f1ad066SMatthew G Knepley #undef __FUNCT__
231521d9a4cSLisandro Dalcin #define __FUNCT__ "DMSetMatType"
232521d9a4cSLisandro Dalcin /*@C
233521d9a4cSLisandro Dalcin        DMSetMatType - Sets the type of matrix created with DMCreateMatrix()
234521d9a4cSLisandro Dalcin 
235521d9a4cSLisandro Dalcin    Logically Collective on DM
236521d9a4cSLisandro Dalcin 
237521d9a4cSLisandro Dalcin    Input Parameter:
238521d9a4cSLisandro Dalcin +  dm - the DM context
239521d9a4cSLisandro Dalcin .  ctype - the matrix type
240521d9a4cSLisandro Dalcin 
241521d9a4cSLisandro Dalcin    Options Database:
242521d9a4cSLisandro Dalcin .   -dm_mat_type ctype
243521d9a4cSLisandro Dalcin 
244521d9a4cSLisandro Dalcin    Level: intermediate
245521d9a4cSLisandro Dalcin 
246c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMGetMatType()
247521d9a4cSLisandro Dalcin @*/
24819fd82e9SBarry Smith PetscErrorCode  DMSetMatType(DM dm,MatType ctype)
249521d9a4cSLisandro Dalcin {
250521d9a4cSLisandro Dalcin   PetscErrorCode ierr;
25188f0584fSBarry Smith 
252521d9a4cSLisandro Dalcin   PetscFunctionBegin;
253521d9a4cSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
254521d9a4cSLisandro Dalcin   ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
25519fd82e9SBarry Smith   ierr = PetscStrallocpy(ctype,(char**)&dm->mattype);CHKERRQ(ierr);
256521d9a4cSLisandro Dalcin   PetscFunctionReturn(0);
257521d9a4cSLisandro Dalcin }
258521d9a4cSLisandro Dalcin 
259521d9a4cSLisandro Dalcin #undef __FUNCT__
260c0dedaeaSBarry Smith #define __FUNCT__ "DMGetMatType"
261c0dedaeaSBarry Smith /*@C
262c0dedaeaSBarry Smith        DMGetMatType - Gets the type of matrix created with DMCreateMatrix()
263c0dedaeaSBarry Smith 
264c0dedaeaSBarry Smith    Logically Collective on DM
265c0dedaeaSBarry Smith 
266c0dedaeaSBarry Smith    Input Parameter:
267c0dedaeaSBarry Smith .  dm - the DM context
268c0dedaeaSBarry Smith 
269c0dedaeaSBarry Smith    Output Parameter:
270c0dedaeaSBarry Smith .  ctype - the matrix type
271c0dedaeaSBarry Smith 
272c0dedaeaSBarry Smith    Options Database:
273c0dedaeaSBarry Smith .   -dm_mat_type ctype
274c0dedaeaSBarry Smith 
275c0dedaeaSBarry Smith    Level: intermediate
276c0dedaeaSBarry Smith 
277c0dedaeaSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMSetMatType()
278c0dedaeaSBarry Smith @*/
279c0dedaeaSBarry Smith PetscErrorCode  DMGetMatType(DM dm,MatType *ctype)
280c0dedaeaSBarry Smith {
281c0dedaeaSBarry Smith   PetscFunctionBegin;
282c0dedaeaSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
283c0dedaeaSBarry Smith   *ctype = dm->mattype;
284c0dedaeaSBarry Smith   PetscFunctionReturn(0);
285c0dedaeaSBarry Smith }
286c0dedaeaSBarry Smith 
287c0dedaeaSBarry Smith #undef __FUNCT__
288c688c046SMatthew G Knepley #define __FUNCT__ "MatGetDM"
289c688c046SMatthew G Knepley /*@
29034f98d34SBarry Smith   MatGetDM - Gets the DM defining the data layout of the matrix
291c688c046SMatthew G Knepley 
292c688c046SMatthew G Knepley   Not collective
293c688c046SMatthew G Knepley 
294c688c046SMatthew G Knepley   Input Parameter:
295c688c046SMatthew G Knepley . A - The Mat
296c688c046SMatthew G Knepley 
297c688c046SMatthew G Knepley   Output Parameter:
298c688c046SMatthew G Knepley . dm - The DM
299c688c046SMatthew G Knepley 
300c688c046SMatthew G Knepley   Level: intermediate
301c688c046SMatthew G Knepley 
302c688c046SMatthew G Knepley .seealso: MatSetDM(), DMCreateMatrix(), DMSetMatType()
303c688c046SMatthew G Knepley @*/
304c688c046SMatthew G Knepley PetscErrorCode MatGetDM(Mat A, DM *dm)
305c688c046SMatthew G Knepley {
306c688c046SMatthew G Knepley   PetscErrorCode ierr;
307c688c046SMatthew G Knepley 
308c688c046SMatthew G Knepley   PetscFunctionBegin;
309c688c046SMatthew G Knepley   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
310c688c046SMatthew G Knepley   PetscValidPointer(dm,2);
311c688c046SMatthew G Knepley   ierr = PetscObjectQuery((PetscObject) A, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr);
312c688c046SMatthew G Knepley   PetscFunctionReturn(0);
313c688c046SMatthew G Knepley }
314c688c046SMatthew G Knepley 
315c688c046SMatthew G Knepley #undef __FUNCT__
316c688c046SMatthew G Knepley #define __FUNCT__ "MatSetDM"
317c688c046SMatthew G Knepley /*@
318c688c046SMatthew G Knepley   MatSetDM - Sets the DM defining the data layout of the matrix
319c688c046SMatthew G Knepley 
320c688c046SMatthew G Knepley   Not collective
321c688c046SMatthew G Knepley 
322c688c046SMatthew G Knepley   Input Parameters:
323c688c046SMatthew G Knepley + A - The Mat
324c688c046SMatthew G Knepley - dm - The DM
325c688c046SMatthew G Knepley 
326c688c046SMatthew G Knepley   Level: intermediate
327c688c046SMatthew G Knepley 
328c688c046SMatthew G Knepley .seealso: MatGetDM(), DMCreateMatrix(), DMSetMatType()
329c688c046SMatthew G Knepley @*/
330c688c046SMatthew G Knepley PetscErrorCode MatSetDM(Mat A, DM dm)
331c688c046SMatthew G Knepley {
332c688c046SMatthew G Knepley   PetscErrorCode ierr;
333c688c046SMatthew G Knepley 
334c688c046SMatthew G Knepley   PetscFunctionBegin;
335c688c046SMatthew G Knepley   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
3368865f1eaSKarl Rupp   if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2);
337c688c046SMatthew G Knepley   ierr = PetscObjectCompose((PetscObject) A, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr);
338c688c046SMatthew G Knepley   PetscFunctionReturn(0);
339c688c046SMatthew G Knepley }
340c688c046SMatthew G Knepley 
341c688c046SMatthew G Knepley #undef __FUNCT__
3429a42bb27SBarry Smith #define __FUNCT__ "DMSetOptionsPrefix"
3439a42bb27SBarry Smith /*@C
3449a42bb27SBarry Smith    DMSetOptionsPrefix - Sets the prefix used for searching for all
345aa219208SBarry Smith    DMDA options in the database.
3469a42bb27SBarry Smith 
347aa219208SBarry Smith    Logically Collective on DMDA
3489a42bb27SBarry Smith 
3499a42bb27SBarry Smith    Input Parameter:
350aa219208SBarry Smith +  da - the DMDA context
3519a42bb27SBarry Smith -  prefix - the prefix to prepend to all option names
3529a42bb27SBarry Smith 
3539a42bb27SBarry Smith    Notes:
3549a42bb27SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
3559a42bb27SBarry Smith    The first character of all runtime options is AUTOMATICALLY the hyphen.
3569a42bb27SBarry Smith 
3579a42bb27SBarry Smith    Level: advanced
3589a42bb27SBarry Smith 
359aa219208SBarry Smith .keywords: DMDA, set, options, prefix, database
3609a42bb27SBarry Smith 
3619a42bb27SBarry Smith .seealso: DMSetFromOptions()
3629a42bb27SBarry Smith @*/
3637087cfbeSBarry Smith PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
3649a42bb27SBarry Smith {
3659a42bb27SBarry Smith   PetscErrorCode ierr;
3669a42bb27SBarry Smith 
3679a42bb27SBarry Smith   PetscFunctionBegin;
3689a42bb27SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3699a42bb27SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
3709a42bb27SBarry Smith   PetscFunctionReturn(0);
3719a42bb27SBarry Smith }
3729a42bb27SBarry Smith 
3739a42bb27SBarry Smith #undef __FUNCT__
37447c6ae99SBarry Smith #define __FUNCT__ "DMDestroy"
37547c6ae99SBarry Smith /*@
376aa219208SBarry Smith     DMDestroy - Destroys a vector packer or DMDA.
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith     Collective on DM
37947c6ae99SBarry Smith 
38047c6ae99SBarry Smith     Input Parameter:
38147c6ae99SBarry Smith .   dm - the DM object to destroy
38247c6ae99SBarry Smith 
38347c6ae99SBarry Smith     Level: developer
38447c6ae99SBarry Smith 
385e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith @*/
388fcfd50ebSBarry Smith PetscErrorCode  DMDestroy(DM *dm)
38947c6ae99SBarry Smith {
390af122d2aSMatthew G Knepley   PetscInt       i, cnt = 0, f;
391dfe15315SJed Brown   DMNamedVecLink nlink,nnext;
39247c6ae99SBarry Smith   PetscErrorCode ierr;
39347c6ae99SBarry Smith 
39447c6ae99SBarry Smith   PetscFunctionBegin;
3956bf464f9SBarry Smith   if (!*dm) PetscFunctionReturn(0);
3966bf464f9SBarry Smith   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
39787e657c6SBarry Smith 
398a546ed68SMatthew G. Knepley   /* I think it makes sense to dump all attached things when you are destroyed, which also eliminates circular references */
399a546ed68SMatthew G. Knepley   for (f = 0; f < (*dm)->numFields; ++f) {
400decb47aaSMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) (*dm)->fields[f], "pmat", NULL);CHKERRQ(ierr);
401decb47aaSMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) (*dm)->fields[f], "nullspace", NULL);CHKERRQ(ierr);
402decb47aaSMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) (*dm)->fields[f], "nearnullspace", NULL);CHKERRQ(ierr);
403a546ed68SMatthew G. Knepley   }
40487e657c6SBarry Smith   /* count all the circular references of DM and its contained Vecs */
405732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
4068865f1eaSKarl Rupp     if ((*dm)->localin[i])  cnt++;
4078865f1eaSKarl Rupp     if ((*dm)->globalin[i]) cnt++;
408732e2eb9SMatthew G Knepley   }
409dfe15315SJed Brown   for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++;
4102348bcf4SPeter Brune   for (nlink=(*dm)->namedlocal; nlink; nlink=nlink->next) cnt++;
4112348bcf4SPeter Brune   if ((*dm)->x) {
4122348bcf4SPeter Brune     DM obj;
4132348bcf4SPeter Brune     ierr = VecGetDM((*dm)->x, &obj);CHKERRQ(ierr);
4142348bcf4SPeter Brune     if (obj == *dm) cnt++;
4152348bcf4SPeter Brune   }
416732e2eb9SMatthew G Knepley 
4176bf464f9SBarry Smith   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
418732e2eb9SMatthew G Knepley   /*
419732e2eb9SMatthew G Knepley      Need this test because the dm references the vectors that
420732e2eb9SMatthew G Knepley      reference the dm, so destroying the dm calls destroy on the
421732e2eb9SMatthew G Knepley      vectors that cause another destroy on the dm
422732e2eb9SMatthew G Knepley   */
4236bf464f9SBarry Smith   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
4246bf464f9SBarry Smith   ((PetscObject) (*dm))->refct = 0;
425732e2eb9SMatthew G Knepley   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
4266bf464f9SBarry Smith     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
4276bf464f9SBarry Smith     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
428732e2eb9SMatthew G Knepley   }
429f490541aSPeter Brune   nnext=(*dm)->namedglobal;
4300298fd71SBarry Smith   (*dm)->namedglobal = NULL;
431f490541aSPeter Brune   for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named vectors */
4322348bcf4SPeter Brune     nnext = nlink->next;
4332348bcf4SPeter 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);
4342348bcf4SPeter Brune     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
4352348bcf4SPeter Brune     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
4362348bcf4SPeter Brune     ierr = PetscFree(nlink);CHKERRQ(ierr);
4372348bcf4SPeter Brune   }
438f490541aSPeter Brune   nnext=(*dm)->namedlocal;
4390298fd71SBarry Smith   (*dm)->namedlocal = NULL;
440f490541aSPeter Brune   for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named local vectors */
441f490541aSPeter Brune     nnext = nlink->next;
442f490541aSPeter 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);
443f490541aSPeter Brune     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
444f490541aSPeter Brune     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
445f490541aSPeter Brune     ierr = PetscFree(nlink);CHKERRQ(ierr);
446f490541aSPeter Brune   }
4472348bcf4SPeter Brune 
448b17ce1afSJed Brown   /* Destroy the list of hooks */
449c833c3b5SJed Brown   {
450c833c3b5SJed Brown     DMCoarsenHookLink link,next;
451b17ce1afSJed Brown     for (link=(*dm)->coarsenhook; link; link=next) {
452b17ce1afSJed Brown       next = link->next;
453b17ce1afSJed Brown       ierr = PetscFree(link);CHKERRQ(ierr);
454b17ce1afSJed Brown     }
4550298fd71SBarry Smith     (*dm)->coarsenhook = NULL;
456c833c3b5SJed Brown   }
457c833c3b5SJed Brown   {
458c833c3b5SJed Brown     DMRefineHookLink link,next;
459c833c3b5SJed Brown     for (link=(*dm)->refinehook; link; link=next) {
460c833c3b5SJed Brown       next = link->next;
461c833c3b5SJed Brown       ierr = PetscFree(link);CHKERRQ(ierr);
462c833c3b5SJed Brown     }
4630298fd71SBarry Smith     (*dm)->refinehook = NULL;
464c833c3b5SJed Brown   }
465be081cd6SPeter Brune   {
466be081cd6SPeter Brune     DMSubDomainHookLink link,next;
467be081cd6SPeter Brune     for (link=(*dm)->subdomainhook; link; link=next) {
468be081cd6SPeter Brune       next = link->next;
469be081cd6SPeter Brune       ierr = PetscFree(link);CHKERRQ(ierr);
470be081cd6SPeter Brune     }
4710298fd71SBarry Smith     (*dm)->subdomainhook = NULL;
472be081cd6SPeter Brune   }
473baf369e7SPeter Brune   {
474baf369e7SPeter Brune     DMGlobalToLocalHookLink link,next;
475baf369e7SPeter Brune     for (link=(*dm)->gtolhook; link; link=next) {
476baf369e7SPeter Brune       next = link->next;
477baf369e7SPeter Brune       ierr = PetscFree(link);CHKERRQ(ierr);
478baf369e7SPeter Brune     }
4790298fd71SBarry Smith     (*dm)->gtolhook = NULL;
480baf369e7SPeter Brune   }
481aa1993deSMatthew G Knepley   /* Destroy the work arrays */
482aa1993deSMatthew G Knepley   {
483aa1993deSMatthew G Knepley     DMWorkLink link,next;
484aa1993deSMatthew G Knepley     if ((*dm)->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out");
485aa1993deSMatthew G Knepley     for (link=(*dm)->workin; link; link=next) {
486aa1993deSMatthew G Knepley       next = link->next;
487aa1993deSMatthew G Knepley       ierr = PetscFree(link->mem);CHKERRQ(ierr);
488aa1993deSMatthew G Knepley       ierr = PetscFree(link);CHKERRQ(ierr);
489aa1993deSMatthew G Knepley     }
4900298fd71SBarry Smith     (*dm)->workin = NULL;
491aa1993deSMatthew G Knepley   }
492b17ce1afSJed Brown 
49352536dc3SBarry Smith   ierr = PetscObjectDestroy(&(*dm)->dmksp);CHKERRQ(ierr);
49452536dc3SBarry Smith   ierr = PetscObjectDestroy(&(*dm)->dmsnes);CHKERRQ(ierr);
49552536dc3SBarry Smith   ierr = PetscObjectDestroy(&(*dm)->dmts);CHKERRQ(ierr);
49652536dc3SBarry Smith 
4971a266240SBarry Smith   if ((*dm)->ctx && (*dm)->ctxdestroy) {
4981a266240SBarry Smith     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
4991a266240SBarry Smith   }
50087e657c6SBarry Smith   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
50171cd77b2SBarry Smith   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
5024dcab191SBarry Smith   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
5036bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
5046bf464f9SBarry Smith   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
505073dac72SJed Brown   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
50688ed4aceSMatthew G Knepley 
50788ed4aceSMatthew G Knepley   ierr = PetscSectionDestroy(&(*dm)->defaultSection);CHKERRQ(ierr);
50888ed4aceSMatthew G Knepley   ierr = PetscSectionDestroy(&(*dm)->defaultGlobalSection);CHKERRQ(ierr);
5098b1ab98fSJed Brown   ierr = PetscLayoutDestroy(&(*dm)->map);CHKERRQ(ierr);
51088ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&(*dm)->sf);CHKERRQ(ierr);
51188ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&(*dm)->defaultSF);CHKERRQ(ierr);
512af122d2aSMatthew G Knepley 
5136636e97aSMatthew G Knepley   ierr = DMDestroy(&(*dm)->coordinateDM);CHKERRQ(ierr);
5146636e97aSMatthew G Knepley   ierr = VecDestroy(&(*dm)->coordinates);CHKERRQ(ierr);
5156636e97aSMatthew G Knepley   ierr = VecDestroy(&(*dm)->coordinatesLocal);CHKERRQ(ierr);
516c6b900c6SMatthew G. Knepley   ierr = PetscFree((*dm)->maxCell);CHKERRQ(ierr);
517c6b900c6SMatthew G. Knepley   ierr = PetscFree((*dm)->L);CHKERRQ(ierr);
5186636e97aSMatthew G Knepley 
519af122d2aSMatthew G Knepley   for (f = 0; f < (*dm)->numFields; ++f) {
520decb47aaSMatthew G. Knepley     ierr = PetscObjectDestroy((PetscObject *) &(*dm)->fields[f]);CHKERRQ(ierr);
521af122d2aSMatthew G Knepley   }
522af122d2aSMatthew G Knepley   ierr = PetscFree((*dm)->fields);CHKERRQ(ierr);
52314f150ffSMatthew G. Knepley   ierr = DMDestroy(&(*dm)->dmBC);CHKERRQ(ierr);
524e04113cfSBarry Smith   /* if memory was published with SAWs then destroy it */
525e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*dm);CHKERRQ(ierr);
526732e2eb9SMatthew G Knepley 
5276bf464f9SBarry Smith   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
528435a35e8SMatthew G Knepley   /* We do not destroy (*dm)->data here so that we can reference count backend objects */
529732e2eb9SMatthew G Knepley   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
53047c6ae99SBarry Smith   PetscFunctionReturn(0);
53147c6ae99SBarry Smith }
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith #undef __FUNCT__
534d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp"
535d7bf68aeSBarry Smith /*@
536d7bf68aeSBarry Smith     DMSetUp - sets up the data structures inside a DM object
537d7bf68aeSBarry Smith 
538d7bf68aeSBarry Smith     Collective on DM
539d7bf68aeSBarry Smith 
540d7bf68aeSBarry Smith     Input Parameter:
541d7bf68aeSBarry Smith .   dm - the DM object to setup
542d7bf68aeSBarry Smith 
543d7bf68aeSBarry Smith     Level: developer
544d7bf68aeSBarry Smith 
545e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
546d7bf68aeSBarry Smith 
547d7bf68aeSBarry Smith @*/
5487087cfbeSBarry Smith PetscErrorCode  DMSetUp(DM dm)
549d7bf68aeSBarry Smith {
550d7bf68aeSBarry Smith   PetscErrorCode ierr;
551d7bf68aeSBarry Smith 
552d7bf68aeSBarry Smith   PetscFunctionBegin;
553171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5548387afaaSJed Brown   if (dm->setupcalled) PetscFunctionReturn(0);
555d7bf68aeSBarry Smith   if (dm->ops->setup) {
556d7bf68aeSBarry Smith     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
557d7bf68aeSBarry Smith   }
5588387afaaSJed Brown   dm->setupcalled = PETSC_TRUE;
559d7bf68aeSBarry Smith   PetscFunctionReturn(0);
560d7bf68aeSBarry Smith }
561d7bf68aeSBarry Smith 
562d7bf68aeSBarry Smith #undef __FUNCT__
563d7bf68aeSBarry Smith #define __FUNCT__ "DMSetFromOptions"
564d7bf68aeSBarry Smith /*@
565d7bf68aeSBarry Smith     DMSetFromOptions - sets parameters in a DM from the options database
566d7bf68aeSBarry Smith 
567d7bf68aeSBarry Smith     Collective on DM
568d7bf68aeSBarry Smith 
569d7bf68aeSBarry Smith     Input Parameter:
570d7bf68aeSBarry Smith .   dm - the DM object to set options for
571d7bf68aeSBarry Smith 
572732e2eb9SMatthew G Knepley     Options Database:
573dd85299cSBarry Smith +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
574dd85299cSBarry Smith .   -dm_vec_type <type>  type of vector to create inside DM
575171400e9SBarry Smith .   -dm_mat_type <type>  type of matrix to create inside DM
576171400e9SBarry Smith -   -dm_coloring_type <global or ghosted>
577732e2eb9SMatthew G Knepley 
578d7bf68aeSBarry Smith     Level: developer
579d7bf68aeSBarry Smith 
580e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
581d7bf68aeSBarry Smith 
582d7bf68aeSBarry Smith @*/
5837087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions(DM dm)
584d7bf68aeSBarry Smith {
5857781c08eSBarry Smith   char           typeName[256];
586ca266f36SBarry Smith   PetscBool      flg;
587d7bf68aeSBarry Smith   PetscErrorCode ierr;
588d7bf68aeSBarry Smith 
589d7bf68aeSBarry Smith   PetscFunctionBegin;
590171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5913194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
5920298fd71SBarry 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);
593a264d7a6SBarry Smith   ierr = PetscOptionsFList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
594f9ba7244SBarry Smith   if (flg) {
595f9ba7244SBarry Smith     ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
596f9ba7244SBarry Smith   }
597a264d7a6SBarry 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);
598073dac72SJed Brown   if (flg) {
599521d9a4cSLisandro Dalcin     ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr);
600073dac72SJed Brown   }
6010298fd71SBarry Smith   ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,NULL);CHKERRQ(ierr);
602f9ba7244SBarry Smith   if (dm->ops->setfromoptions) {
603f9ba7244SBarry Smith     ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
604f9ba7244SBarry Smith   }
605f9ba7244SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
606f9ba7244SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
60782fcb398SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
608d7bf68aeSBarry Smith   PetscFunctionReturn(0);
609d7bf68aeSBarry Smith }
610d7bf68aeSBarry Smith 
611d7bf68aeSBarry Smith #undef __FUNCT__
61247c6ae99SBarry Smith #define __FUNCT__ "DMView"
613fc9bc008SSatish Balay /*@C
614aa219208SBarry Smith     DMView - Views a vector packer or DMDA.
61547c6ae99SBarry Smith 
61647c6ae99SBarry Smith     Collective on DM
61747c6ae99SBarry Smith 
61847c6ae99SBarry Smith     Input Parameter:
61947c6ae99SBarry Smith +   dm - the DM object to view
62047c6ae99SBarry Smith -   v - the viewer
62147c6ae99SBarry Smith 
62247c6ae99SBarry Smith     Level: developer
62347c6ae99SBarry Smith 
624e727c939SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
62547c6ae99SBarry Smith 
62647c6ae99SBarry Smith @*/
6277087cfbeSBarry Smith PetscErrorCode  DMView(DM dm,PetscViewer v)
62847c6ae99SBarry Smith {
62947c6ae99SBarry Smith   PetscErrorCode ierr;
63032c0f0efSBarry Smith   PetscBool      isbinary;
63147c6ae99SBarry Smith 
63247c6ae99SBarry Smith   PetscFunctionBegin;
633171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6343014e516SBarry Smith   if (!v) {
635ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm),&v);CHKERRQ(ierr);
6363014e516SBarry Smith   }
63798c3331eSBarry Smith   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)dm,v);CHKERRQ(ierr);
63832c0f0efSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
63932c0f0efSBarry Smith   if (isbinary) {
64055849f57SBarry Smith     PetscInt classid = DM_FILE_CLASSID;
64132c0f0efSBarry Smith     char     type[256];
64232c0f0efSBarry Smith 
64332c0f0efSBarry Smith     ierr = PetscViewerBinaryWrite(v,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
64432c0f0efSBarry Smith     ierr = PetscStrncpy(type,((PetscObject)dm)->type_name,256);CHKERRQ(ierr);
64532c0f0efSBarry Smith     ierr = PetscViewerBinaryWrite(v,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
64632c0f0efSBarry Smith   }
6470c010503SBarry Smith   if (dm->ops->view) {
6480c010503SBarry Smith     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
64947c6ae99SBarry Smith   }
65047c6ae99SBarry Smith   PetscFunctionReturn(0);
65147c6ae99SBarry Smith }
65247c6ae99SBarry Smith 
65347c6ae99SBarry Smith #undef __FUNCT__
65447c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector"
65547c6ae99SBarry Smith /*@
656aa219208SBarry Smith     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
65747c6ae99SBarry Smith 
65847c6ae99SBarry Smith     Collective on DM
65947c6ae99SBarry Smith 
66047c6ae99SBarry Smith     Input Parameter:
66147c6ae99SBarry Smith .   dm - the DM object
66247c6ae99SBarry Smith 
66347c6ae99SBarry Smith     Output Parameter:
66447c6ae99SBarry Smith .   vec - the global vector
66547c6ae99SBarry Smith 
666073dac72SJed Brown     Level: beginner
66747c6ae99SBarry Smith 
668e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
66947c6ae99SBarry Smith 
67047c6ae99SBarry Smith @*/
6717087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
67247c6ae99SBarry Smith {
67347c6ae99SBarry Smith   PetscErrorCode ierr;
67447c6ae99SBarry Smith 
67547c6ae99SBarry Smith   PetscFunctionBegin;
676171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
67747c6ae99SBarry Smith   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
67847c6ae99SBarry Smith   PetscFunctionReturn(0);
67947c6ae99SBarry Smith }
68047c6ae99SBarry Smith 
68147c6ae99SBarry Smith #undef __FUNCT__
68247c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector"
68347c6ae99SBarry Smith /*@
684aa219208SBarry Smith     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith     Not Collective
68747c6ae99SBarry Smith 
68847c6ae99SBarry Smith     Input Parameter:
68947c6ae99SBarry Smith .   dm - the DM object
69047c6ae99SBarry Smith 
69147c6ae99SBarry Smith     Output Parameter:
69247c6ae99SBarry Smith .   vec - the local vector
69347c6ae99SBarry Smith 
694073dac72SJed Brown     Level: beginner
69547c6ae99SBarry Smith 
696e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
69747c6ae99SBarry Smith 
69847c6ae99SBarry Smith @*/
6997087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
70047c6ae99SBarry Smith {
70147c6ae99SBarry Smith   PetscErrorCode ierr;
70247c6ae99SBarry Smith 
70347c6ae99SBarry Smith   PetscFunctionBegin;
704171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
70547c6ae99SBarry Smith   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
70647c6ae99SBarry Smith   PetscFunctionReturn(0);
70747c6ae99SBarry Smith }
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith #undef __FUNCT__
7101411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping"
7111411c6eeSJed Brown /*@
7121411c6eeSJed Brown    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
7131411c6eeSJed Brown 
7141411c6eeSJed Brown    Collective on DM
7151411c6eeSJed Brown 
7161411c6eeSJed Brown    Input Parameter:
7171411c6eeSJed Brown .  dm - the DM that provides the mapping
7181411c6eeSJed Brown 
7191411c6eeSJed Brown    Output Parameter:
7201411c6eeSJed Brown .  ltog - the mapping
7211411c6eeSJed Brown 
7221411c6eeSJed Brown    Level: intermediate
7231411c6eeSJed Brown 
7241411c6eeSJed Brown    Notes:
7251411c6eeSJed Brown    This mapping can then be used by VecSetLocalToGlobalMapping() or
7261411c6eeSJed Brown    MatSetLocalToGlobalMapping().
7271411c6eeSJed Brown 
728*fc31e74dSBarry Smith .seealso: DMCreateLocalVector()
7291411c6eeSJed Brown @*/
7307087cfbeSBarry Smith PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
7311411c6eeSJed Brown {
7321411c6eeSJed Brown   PetscErrorCode ierr;
7331411c6eeSJed Brown 
7341411c6eeSJed Brown   PetscFunctionBegin;
7351411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7361411c6eeSJed Brown   PetscValidPointer(ltog,2);
7371411c6eeSJed Brown   if (!dm->ltogmap) {
73837d0c07bSMatthew G Knepley     PetscSection section, sectionGlobal;
73937d0c07bSMatthew G Knepley 
74037d0c07bSMatthew G Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
74137d0c07bSMatthew G Knepley     if (section) {
74237d0c07bSMatthew G Knepley       PetscInt *ltog;
74337d0c07bSMatthew G Knepley       PetscInt pStart, pEnd, size, p, l;
74437d0c07bSMatthew G Knepley 
74537d0c07bSMatthew G Knepley       ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
74637d0c07bSMatthew G Knepley       ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
74737d0c07bSMatthew G Knepley       ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
748785e854fSJed Brown       ierr = PetscMalloc1(size, &ltog);CHKERRQ(ierr); /* We want the local+overlap size */
74937d0c07bSMatthew G Knepley       for (p = pStart, l = 0; p < pEnd; ++p) {
75037d0c07bSMatthew G Knepley         PetscInt dof, off, c;
75137d0c07bSMatthew G Knepley 
75237d0c07bSMatthew G Knepley         /* Should probably use constrained dofs */
75337d0c07bSMatthew G Knepley         ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
75437d0c07bSMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &off);CHKERRQ(ierr);
75537d0c07bSMatthew G Knepley         for (c = 0; c < dof; ++c, ++l) {
75637d0c07bSMatthew G Knepley           ltog[l] = off+c;
75737d0c07bSMatthew G Knepley         }
75837d0c07bSMatthew G Knepley       }
759f0413b6fSBarry Smith       ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, 1,size, ltog, PETSC_OWN_POINTER, &dm->ltogmap);CHKERRQ(ierr);
7603bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)dm, (PetscObject)dm->ltogmap);CHKERRQ(ierr);
76137d0c07bSMatthew G Knepley     } else {
762184d77edSJed Brown       if (!dm->ops->getlocaltoglobalmapping) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
763184d77edSJed Brown       ierr = (*dm->ops->getlocaltoglobalmapping)(dm);CHKERRQ(ierr);
7641411c6eeSJed Brown     }
76537d0c07bSMatthew G Knepley   }
7661411c6eeSJed Brown   *ltog = dm->ltogmap;
7671411c6eeSJed Brown   PetscFunctionReturn(0);
7681411c6eeSJed Brown }
7691411c6eeSJed Brown 
7701411c6eeSJed Brown #undef __FUNCT__
7711411c6eeSJed Brown #define __FUNCT__ "DMGetBlockSize"
7721411c6eeSJed Brown /*@
7731411c6eeSJed Brown    DMGetBlockSize - Gets the inherent block size associated with a DM
7741411c6eeSJed Brown 
7751411c6eeSJed Brown    Not Collective
7761411c6eeSJed Brown 
7771411c6eeSJed Brown    Input Parameter:
7781411c6eeSJed Brown .  dm - the DM with block structure
7791411c6eeSJed Brown 
7801411c6eeSJed Brown    Output Parameter:
7811411c6eeSJed Brown .  bs - the block size, 1 implies no exploitable block structure
7821411c6eeSJed Brown 
7831411c6eeSJed Brown    Level: intermediate
7841411c6eeSJed Brown 
785*fc31e74dSBarry Smith .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMapping()
7861411c6eeSJed Brown @*/
7877087cfbeSBarry Smith PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
7881411c6eeSJed Brown {
7891411c6eeSJed Brown   PetscFunctionBegin;
7901411c6eeSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7911411c6eeSJed Brown   PetscValidPointer(bs,2);
7921411c6eeSJed 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");
7931411c6eeSJed Brown   *bs = dm->bs;
7941411c6eeSJed Brown   PetscFunctionReturn(0);
7951411c6eeSJed Brown }
7961411c6eeSJed Brown 
7971411c6eeSJed Brown #undef __FUNCT__
798e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation"
79947c6ae99SBarry Smith /*@
800e727c939SJed Brown     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith     Collective on DM
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith     Input Parameter:
80547c6ae99SBarry Smith +   dm1 - the DM object
80647c6ae99SBarry Smith -   dm2 - the second, finer DM object
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith     Output Parameter:
80947c6ae99SBarry Smith +  mat - the interpolation
81047c6ae99SBarry Smith -  vec - the scaling (optional)
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith     Level: developer
81347c6ae99SBarry Smith 
81485afcc9aSBarry 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
81585afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
816d52bd9f3SBarry Smith 
8171f588964SMatthew G Knepley         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMGetCoordinateDM()) to interpolate the mesh coordinate vectors
818d52bd9f3SBarry Smith         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
81985afcc9aSBarry Smith 
82085afcc9aSBarry Smith 
821e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
82247c6ae99SBarry Smith 
82347c6ae99SBarry Smith @*/
824e727c939SJed Brown PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
82547c6ae99SBarry Smith {
82647c6ae99SBarry Smith   PetscErrorCode ierr;
82747c6ae99SBarry Smith 
82847c6ae99SBarry Smith   PetscFunctionBegin;
829171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
830171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
83125296bd5SBarry Smith   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
83247c6ae99SBarry Smith   PetscFunctionReturn(0);
83347c6ae99SBarry Smith }
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith #undef __FUNCT__
836e727c939SJed Brown #define __FUNCT__ "DMCreateInjection"
83747c6ae99SBarry Smith /*@
838e727c939SJed Brown     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith     Collective on DM
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith     Input Parameter:
84347c6ae99SBarry Smith +   dm1 - the DM object
84447c6ae99SBarry Smith -   dm2 - the second, finer DM object
84547c6ae99SBarry Smith 
84647c6ae99SBarry Smith     Output Parameter:
84747c6ae99SBarry Smith .   ctx - the injection
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith     Level: developer
85047c6ae99SBarry Smith 
85185afcc9aSBarry 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
85285afcc9aSBarry Smith         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
85385afcc9aSBarry Smith 
854e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith @*/
857e727c939SJed Brown PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
85847c6ae99SBarry Smith {
85947c6ae99SBarry Smith   PetscErrorCode ierr;
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith   PetscFunctionBegin;
862171400e9SBarry Smith   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
863171400e9SBarry Smith   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
86447c6ae99SBarry Smith   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
86547c6ae99SBarry Smith   PetscFunctionReturn(0);
86647c6ae99SBarry Smith }
86747c6ae99SBarry Smith 
86847c6ae99SBarry Smith #undef __FUNCT__
869e727c939SJed Brown #define __FUNCT__ "DMCreateColoring"
870b412c318SBarry Smith /*@
871b412c318SBarry Smith     DMCreateColoring - Gets coloring for a DM
87247c6ae99SBarry Smith 
87347c6ae99SBarry Smith     Collective on DM
87447c6ae99SBarry Smith 
87547c6ae99SBarry Smith     Input Parameter:
87647c6ae99SBarry Smith +   dm - the DM object
877b412c318SBarry Smith -   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
87847c6ae99SBarry Smith 
87947c6ae99SBarry Smith     Output Parameter:
88047c6ae99SBarry Smith .   coloring - the coloring
88147c6ae99SBarry Smith 
88247c6ae99SBarry Smith     Level: developer
88347c6ae99SBarry Smith 
884b412c318SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix(), DMSetMatType()
88547c6ae99SBarry Smith 
886aab9d709SJed Brown @*/
887b412c318SBarry Smith PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,ISColoring *coloring)
88847c6ae99SBarry Smith {
88947c6ae99SBarry Smith   PetscErrorCode ierr;
89047c6ae99SBarry Smith 
89147c6ae99SBarry Smith   PetscFunctionBegin;
892171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
893ce94432eSBarry Smith   if (!dm->ops->getcoloring) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No coloring for this type of DM yet");
894b412c318SBarry Smith   ierr = (*dm->ops->getcoloring)(dm,ctype,coloring);CHKERRQ(ierr);
89547c6ae99SBarry Smith   PetscFunctionReturn(0);
89647c6ae99SBarry Smith }
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith #undef __FUNCT__
899950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix"
900b412c318SBarry Smith /*@
901950540a4SJed Brown     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
90247c6ae99SBarry Smith 
90347c6ae99SBarry Smith     Collective on DM
90447c6ae99SBarry Smith 
90547c6ae99SBarry Smith     Input Parameter:
906b412c318SBarry Smith .   dm - the DM object
90747c6ae99SBarry Smith 
90847c6ae99SBarry Smith     Output Parameter:
90947c6ae99SBarry Smith .   mat - the empty Jacobian
91047c6ae99SBarry Smith 
911073dac72SJed Brown     Level: beginner
91247c6ae99SBarry Smith 
91394013140SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
91494013140SBarry Smith        do not need to do it yourself.
91594013140SBarry Smith 
91694013140SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
917aa219208SBarry Smith        the nonzero pattern call DMDASetMatPreallocateOnly()
91894013140SBarry Smith 
91994013140SBarry 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
92094013140SBarry Smith        internally by PETSc.
92194013140SBarry Smith 
92294013140SBarry Smith        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
923aa219208SBarry Smith        the indices for the global numbering for DMDAs which is complicated.
92494013140SBarry Smith 
925b412c318SBarry Smith .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMSetMatType()
92647c6ae99SBarry Smith 
927aab9d709SJed Brown @*/
928b412c318SBarry Smith PetscErrorCode  DMCreateMatrix(DM dm,Mat *mat)
92947c6ae99SBarry Smith {
93047c6ae99SBarry Smith   PetscErrorCode ierr;
93147c6ae99SBarry Smith 
93247c6ae99SBarry Smith   PetscFunctionBegin;
933171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
934607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
935c7b7c8a4SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
936c7b7c8a4SJed Brown   PetscValidPointer(mat,3);
937b412c318SBarry Smith   ierr = (*dm->ops->creatematrix)(dm,mat);CHKERRQ(ierr);
93847c6ae99SBarry Smith   PetscFunctionReturn(0);
93947c6ae99SBarry Smith }
94047c6ae99SBarry Smith 
94147c6ae99SBarry Smith #undef __FUNCT__
942732e2eb9SMatthew G Knepley #define __FUNCT__ "DMSetMatrixPreallocateOnly"
943732e2eb9SMatthew G Knepley /*@
944950540a4SJed Brown   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
945732e2eb9SMatthew G Knepley     preallocated but the nonzero structure and zero values will not be set.
946732e2eb9SMatthew G Knepley 
947732e2eb9SMatthew G Knepley   Logically Collective on DMDA
948732e2eb9SMatthew G Knepley 
949732e2eb9SMatthew G Knepley   Input Parameter:
950732e2eb9SMatthew G Knepley + dm - the DM
951732e2eb9SMatthew G Knepley - only - PETSC_TRUE if only want preallocation
952732e2eb9SMatthew G Knepley 
953732e2eb9SMatthew G Knepley   Level: developer
954950540a4SJed Brown .seealso DMCreateMatrix()
955732e2eb9SMatthew G Knepley @*/
956732e2eb9SMatthew G Knepley PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
957732e2eb9SMatthew G Knepley {
958732e2eb9SMatthew G Knepley   PetscFunctionBegin;
959732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
960732e2eb9SMatthew G Knepley   dm->prealloc_only = only;
961732e2eb9SMatthew G Knepley   PetscFunctionReturn(0);
962732e2eb9SMatthew G Knepley }
963732e2eb9SMatthew G Knepley 
964732e2eb9SMatthew G Knepley #undef __FUNCT__
965a89ea682SMatthew G Knepley #define __FUNCT__ "DMGetWorkArray"
966a89ea682SMatthew G Knepley /*@C
967aa1993deSMatthew G Knepley   DMGetWorkArray - Gets a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray()
968a89ea682SMatthew G Knepley 
969a89ea682SMatthew G Knepley   Not Collective
970a89ea682SMatthew G Knepley 
971a89ea682SMatthew G Knepley   Input Parameters:
972a89ea682SMatthew G Knepley + dm - the DM object
973aa1993deSMatthew G Knepley . count - The minium size
974aa1993deSMatthew G Knepley - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT)
975a89ea682SMatthew G Knepley 
976a89ea682SMatthew G Knepley   Output Parameter:
977a89ea682SMatthew G Knepley . array - the work array
978a89ea682SMatthew G Knepley 
979a89ea682SMatthew G Knepley   Level: developer
980a89ea682SMatthew G Knepley 
981a89ea682SMatthew G Knepley .seealso DMDestroy(), DMCreate()
982a89ea682SMatthew G Knepley @*/
983aa1993deSMatthew G Knepley PetscErrorCode DMGetWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem)
984a89ea682SMatthew G Knepley {
985a89ea682SMatthew G Knepley   PetscErrorCode ierr;
986aa1993deSMatthew G Knepley   DMWorkLink     link;
987aa1993deSMatthew G Knepley   size_t         size;
988a89ea682SMatthew G Knepley 
989a89ea682SMatthew G Knepley   PetscFunctionBegin;
990a89ea682SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
991aa1993deSMatthew G Knepley   PetscValidPointer(mem,4);
992aa1993deSMatthew G Knepley   if (dm->workin) {
993aa1993deSMatthew G Knepley     link       = dm->workin;
994aa1993deSMatthew G Knepley     dm->workin = dm->workin->next;
995aa1993deSMatthew G Knepley   } else {
996b00a9115SJed Brown     ierr = PetscNewLog(dm,&link);CHKERRQ(ierr);
997a89ea682SMatthew G Knepley   }
998aa1993deSMatthew G Knepley   ierr = PetscDataTypeGetSize(dtype,&size);CHKERRQ(ierr);
999aa1993deSMatthew G Knepley   if (size*count > link->bytes) {
1000aa1993deSMatthew G Knepley     ierr        = PetscFree(link->mem);CHKERRQ(ierr);
1001aa1993deSMatthew G Knepley     ierr        = PetscMalloc(size*count,&link->mem);CHKERRQ(ierr);
1002aa1993deSMatthew G Knepley     link->bytes = size*count;
1003aa1993deSMatthew G Knepley   }
1004aa1993deSMatthew G Knepley   link->next   = dm->workout;
1005aa1993deSMatthew G Knepley   dm->workout  = link;
1006aa1993deSMatthew G Knepley   *(void**)mem = link->mem;
1007a89ea682SMatthew G Knepley   PetscFunctionReturn(0);
1008a89ea682SMatthew G Knepley }
1009a89ea682SMatthew G Knepley 
1010aa1993deSMatthew G Knepley #undef __FUNCT__
1011aa1993deSMatthew G Knepley #define __FUNCT__ "DMRestoreWorkArray"
1012aa1993deSMatthew G Knepley /*@C
1013aa1993deSMatthew G Knepley   DMRestoreWorkArray - Restores a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray()
1014aa1993deSMatthew G Knepley 
1015aa1993deSMatthew G Knepley   Not Collective
1016aa1993deSMatthew G Knepley 
1017aa1993deSMatthew G Knepley   Input Parameters:
1018aa1993deSMatthew G Knepley + dm - the DM object
1019aa1993deSMatthew G Knepley . count - The minium size
1020aa1993deSMatthew G Knepley - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT)
1021aa1993deSMatthew G Knepley 
1022aa1993deSMatthew G Knepley   Output Parameter:
1023aa1993deSMatthew G Knepley . array - the work array
1024aa1993deSMatthew G Knepley 
1025aa1993deSMatthew G Knepley   Level: developer
1026aa1993deSMatthew G Knepley 
1027aa1993deSMatthew G Knepley .seealso DMDestroy(), DMCreate()
1028aa1993deSMatthew G Knepley @*/
1029aa1993deSMatthew G Knepley PetscErrorCode DMRestoreWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem)
1030aa1993deSMatthew G Knepley {
1031aa1993deSMatthew G Knepley   DMWorkLink *p,link;
1032aa1993deSMatthew G Knepley 
1033aa1993deSMatthew G Knepley   PetscFunctionBegin;
1034aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1035aa1993deSMatthew G Knepley   PetscValidPointer(mem,4);
1036aa1993deSMatthew G Knepley   for (p=&dm->workout; (link=*p); p=&link->next) {
1037aa1993deSMatthew G Knepley     if (link->mem == *(void**)mem) {
1038aa1993deSMatthew G Knepley       *p           = link->next;
1039aa1993deSMatthew G Knepley       link->next   = dm->workin;
1040aa1993deSMatthew G Knepley       dm->workin   = link;
10410298fd71SBarry Smith       *(void**)mem = NULL;
1042aa1993deSMatthew G Knepley       PetscFunctionReturn(0);
1043aa1993deSMatthew G Knepley     }
1044aa1993deSMatthew G Knepley   }
1045aa1993deSMatthew G Knepley   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out");
1046aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
1047aa1993deSMatthew G Knepley }
1048e7c4fc90SDmitry Karpeev 
1049e7c4fc90SDmitry Karpeev #undef __FUNCT__
1050435a35e8SMatthew G Knepley #define __FUNCT__ "DMSetNullSpaceConstructor"
1051435a35e8SMatthew G Knepley PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt field, MatNullSpace *nullSpace))
1052435a35e8SMatthew G Knepley {
1053435a35e8SMatthew G Knepley   PetscFunctionBegin;
1054435a35e8SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
105582f516ccSBarry Smith   if (field >= 10) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %d >= 10 fields", field);
1056435a35e8SMatthew G Knepley   dm->nullspaceConstructors[field] = nullsp;
1057435a35e8SMatthew G Knepley   PetscFunctionReturn(0);
1058435a35e8SMatthew G Knepley }
1059435a35e8SMatthew G Knepley 
1060435a35e8SMatthew G Knepley #undef __FUNCT__
10614d343eeaSMatthew G Knepley #define __FUNCT__ "DMCreateFieldIS"
10624f3b5142SJed Brown /*@C
10634d343eeaSMatthew G Knepley   DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field
10644d343eeaSMatthew G Knepley 
10654d343eeaSMatthew G Knepley   Not collective
10664d343eeaSMatthew G Knepley 
10674d343eeaSMatthew G Knepley   Input Parameter:
10684d343eeaSMatthew G Knepley . dm - the DM object
10694d343eeaSMatthew G Knepley 
10704d343eeaSMatthew G Knepley   Output Parameters:
10710298fd71SBarry Smith + numFields  - The number of fields (or NULL if not requested)
10720298fd71SBarry Smith . fieldNames - The name for each field (or NULL if not requested)
10730298fd71SBarry Smith - fields     - The global indices for each field (or NULL if not requested)
10744d343eeaSMatthew G Knepley 
10754d343eeaSMatthew G Knepley   Level: intermediate
10764d343eeaSMatthew G Knepley 
107721c9b008SJed Brown   Notes:
107821c9b008SJed Brown   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
107921c9b008SJed Brown   PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with
108021c9b008SJed Brown   PetscFree().
108121c9b008SJed Brown 
10824d343eeaSMatthew G Knepley .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
10834d343eeaSMatthew G Knepley @*/
108437d0c07bSMatthew G Knepley PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
10854d343eeaSMatthew G Knepley {
108637d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
10874d343eeaSMatthew G Knepley   PetscErrorCode ierr;
10884d343eeaSMatthew G Knepley 
10894d343eeaSMatthew G Knepley   PetscFunctionBegin;
10904d343eeaSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
109169ca1f37SDmitry Karpeev   if (numFields) {
109269ca1f37SDmitry Karpeev     PetscValidPointer(numFields,2);
109369ca1f37SDmitry Karpeev     *numFields = 0;
109469ca1f37SDmitry Karpeev   }
109537d0c07bSMatthew G Knepley   if (fieldNames) {
109637d0c07bSMatthew G Knepley     PetscValidPointer(fieldNames,3);
10970298fd71SBarry Smith     *fieldNames = NULL;
109869ca1f37SDmitry Karpeev   }
109969ca1f37SDmitry Karpeev   if (fields) {
110069ca1f37SDmitry Karpeev     PetscValidPointer(fields,4);
11010298fd71SBarry Smith     *fields = NULL;
110269ca1f37SDmitry Karpeev   }
110337d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
110437d0c07bSMatthew G Knepley   if (section) {
110537d0c07bSMatthew G Knepley     PetscInt *fieldSizes, **fieldIndices;
110637d0c07bSMatthew G Knepley     PetscInt nF, f, pStart, pEnd, p;
110737d0c07bSMatthew G Knepley 
110837d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
110937d0c07bSMatthew G Knepley     ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
1110dcca6d9dSJed Brown     ierr = PetscMalloc2(nF,&fieldSizes,nF,&fieldIndices);CHKERRQ(ierr);
111137d0c07bSMatthew G Knepley     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
111237d0c07bSMatthew G Knepley     for (f = 0; f < nF; ++f) {
111337d0c07bSMatthew G Knepley       fieldSizes[f] = 0;
111437d0c07bSMatthew G Knepley     }
111537d0c07bSMatthew G Knepley     for (p = pStart; p < pEnd; ++p) {
111637d0c07bSMatthew G Knepley       PetscInt gdof;
111737d0c07bSMatthew G Knepley 
111837d0c07bSMatthew G Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
111937d0c07bSMatthew G Knepley       if (gdof > 0) {
112037d0c07bSMatthew G Knepley         for (f = 0; f < nF; ++f) {
112137d0c07bSMatthew G Knepley           PetscInt fdof, fcdof;
112237d0c07bSMatthew G Knepley 
112337d0c07bSMatthew G Knepley           ierr           = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
112437d0c07bSMatthew G Knepley           ierr           = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
112537d0c07bSMatthew G Knepley           fieldSizes[f] += fdof-fcdof;
112637d0c07bSMatthew G Knepley         }
112737d0c07bSMatthew G Knepley       }
112837d0c07bSMatthew G Knepley     }
112937d0c07bSMatthew G Knepley     for (f = 0; f < nF; ++f) {
1130785e854fSJed Brown       ierr          = PetscMalloc1(fieldSizes[f], &fieldIndices[f]);CHKERRQ(ierr);
113137d0c07bSMatthew G Knepley       fieldSizes[f] = 0;
113237d0c07bSMatthew G Knepley     }
113337d0c07bSMatthew G Knepley     for (p = pStart; p < pEnd; ++p) {
113437d0c07bSMatthew G Knepley       PetscInt gdof, goff;
113537d0c07bSMatthew G Knepley 
113637d0c07bSMatthew G Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
113737d0c07bSMatthew G Knepley       if (gdof > 0) {
113837d0c07bSMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
113937d0c07bSMatthew G Knepley         for (f = 0; f < nF; ++f) {
114037d0c07bSMatthew G Knepley           PetscInt fdof, fcdof, fc;
114137d0c07bSMatthew G Knepley 
114237d0c07bSMatthew G Knepley           ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
114337d0c07bSMatthew G Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
114437d0c07bSMatthew G Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++fieldSizes[f]) {
114537d0c07bSMatthew G Knepley             fieldIndices[f][fieldSizes[f]] = goff++;
114637d0c07bSMatthew G Knepley           }
114737d0c07bSMatthew G Knepley         }
114837d0c07bSMatthew G Knepley       }
114937d0c07bSMatthew G Knepley     }
11508865f1eaSKarl Rupp     if (numFields) *numFields = nF;
115137d0c07bSMatthew G Knepley     if (fieldNames) {
1152785e854fSJed Brown       ierr = PetscMalloc1(nF, fieldNames);CHKERRQ(ierr);
115337d0c07bSMatthew G Knepley       for (f = 0; f < nF; ++f) {
115437d0c07bSMatthew G Knepley         const char *fieldName;
115537d0c07bSMatthew G Knepley 
115637d0c07bSMatthew G Knepley         ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr);
115737d0c07bSMatthew G Knepley         ierr = PetscStrallocpy(fieldName, (char**) &(*fieldNames)[f]);CHKERRQ(ierr);
115837d0c07bSMatthew G Knepley       }
115937d0c07bSMatthew G Knepley     }
116037d0c07bSMatthew G Knepley     if (fields) {
1161785e854fSJed Brown       ierr = PetscMalloc1(nF, fields);CHKERRQ(ierr);
116237d0c07bSMatthew G Knepley       for (f = 0; f < nF; ++f) {
116382f516ccSBarry Smith         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]);CHKERRQ(ierr);
116437d0c07bSMatthew G Knepley       }
116537d0c07bSMatthew G Knepley     }
116637d0c07bSMatthew G Knepley     ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr);
11678865f1eaSKarl Rupp   } else if (dm->ops->createfieldis) {
11688865f1eaSKarl Rupp     ierr = (*dm->ops->createfieldis)(dm, numFields, fieldNames, fields);CHKERRQ(ierr);
116969ca1f37SDmitry Karpeev   }
11704d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
11714d343eeaSMatthew G Knepley }
11724d343eeaSMatthew G Knepley 
117316621825SDmitry Karpeev 
117416621825SDmitry Karpeev #undef __FUNCT__
117516621825SDmitry Karpeev #define __FUNCT__ "DMCreateFieldDecomposition"
117616621825SDmitry Karpeev /*@C
117716621825SDmitry Karpeev   DMCreateFieldDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems
117816621825SDmitry Karpeev                           corresponding to different fields: each IS contains the global indices of the dofs of the
117916621825SDmitry Karpeev                           corresponding field. The optional list of DMs define the DM for each subproblem.
1180e7c4fc90SDmitry Karpeev                           Generalizes DMCreateFieldIS().
1181e7c4fc90SDmitry Karpeev 
1182e7c4fc90SDmitry Karpeev   Not collective
1183e7c4fc90SDmitry Karpeev 
1184e7c4fc90SDmitry Karpeev   Input Parameter:
1185e7c4fc90SDmitry Karpeev . dm - the DM object
1186e7c4fc90SDmitry Karpeev 
1187e7c4fc90SDmitry Karpeev   Output Parameters:
11880298fd71SBarry Smith + len       - The number of subproblems in the field decomposition (or NULL if not requested)
11890298fd71SBarry Smith . namelist  - The name for each field (or NULL if not requested)
11900298fd71SBarry Smith . islist    - The global indices for each field (or NULL if not requested)
11910298fd71SBarry Smith - dmlist    - The DMs for each field subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined)
1192e7c4fc90SDmitry Karpeev 
1193e7c4fc90SDmitry Karpeev   Level: intermediate
1194e7c4fc90SDmitry Karpeev 
1195e7c4fc90SDmitry Karpeev   Notes:
1196e7c4fc90SDmitry Karpeev   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
1197e7c4fc90SDmitry Karpeev   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
1198e7c4fc90SDmitry Karpeev   and all of the arrays should be freed with PetscFree().
1199e7c4fc90SDmitry Karpeev 
1200e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1201e7c4fc90SDmitry Karpeev @*/
120216621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1203e7c4fc90SDmitry Karpeev {
1204e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
1205e7c4fc90SDmitry Karpeev 
1206e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
1207e7c4fc90SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12088865f1eaSKarl Rupp   if (len) {
12098865f1eaSKarl Rupp     PetscValidPointer(len,2);
12108865f1eaSKarl Rupp     *len = 0;
12118865f1eaSKarl Rupp   }
12128865f1eaSKarl Rupp   if (namelist) {
12138865f1eaSKarl Rupp     PetscValidPointer(namelist,3);
12148865f1eaSKarl Rupp     *namelist = 0;
12158865f1eaSKarl Rupp   }
12168865f1eaSKarl Rupp   if (islist) {
12178865f1eaSKarl Rupp     PetscValidPointer(islist,4);
12188865f1eaSKarl Rupp     *islist = 0;
12198865f1eaSKarl Rupp   }
12208865f1eaSKarl Rupp   if (dmlist) {
12218865f1eaSKarl Rupp     PetscValidPointer(dmlist,5);
12228865f1eaSKarl Rupp     *dmlist = 0;
12238865f1eaSKarl Rupp   }
1224f3f0edfdSDmitry Karpeev   /*
1225f3f0edfdSDmitry Karpeev    Is it a good idea to apply the following check across all impls?
1226f3f0edfdSDmitry Karpeev    Perhaps some impls can have a well-defined decomposition before DMSetUp?
1227f3f0edfdSDmitry Karpeev    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
1228f3f0edfdSDmitry Karpeev    */
1229ce94432eSBarry Smith   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
123016621825SDmitry Karpeev   if (!dm->ops->createfielddecomposition) {
1231435a35e8SMatthew G Knepley     PetscSection section;
1232435a35e8SMatthew G Knepley     PetscInt     numFields, f;
1233435a35e8SMatthew G Knepley 
1234435a35e8SMatthew G Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1235435a35e8SMatthew G Knepley     if (section) {ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);}
1236435a35e8SMatthew G Knepley     if (section && numFields && dm->ops->createsubdm) {
1237435a35e8SMatthew G Knepley       *len = numFields;
123803dc3394SMatthew G. Knepley       if (namelist) {ierr = PetscMalloc1(numFields,namelist);CHKERRQ(ierr);}
123903dc3394SMatthew G. Knepley       if (islist)   {ierr = PetscMalloc1(numFields,islist);CHKERRQ(ierr);}
124003dc3394SMatthew G. Knepley       if (dmlist)   {ierr = PetscMalloc1(numFields,dmlist);CHKERRQ(ierr);}
1241435a35e8SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
1242435a35e8SMatthew G Knepley         const char *fieldName;
1243435a35e8SMatthew G Knepley 
124403dc3394SMatthew G. Knepley         ierr = DMCreateSubDM(dm, 1, &f, islist ? &(*islist)[f] : NULL, dmlist ? &(*dmlist)[f] : NULL);CHKERRQ(ierr);
124503dc3394SMatthew G. Knepley         if (namelist) {
1246435a35e8SMatthew G Knepley           ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr);
1247435a35e8SMatthew G Knepley           ierr = PetscStrallocpy(fieldName, (char**) &(*namelist)[f]);CHKERRQ(ierr);
1248435a35e8SMatthew G Knepley         }
124903dc3394SMatthew G. Knepley       }
1250435a35e8SMatthew G Knepley     } else {
125169ca1f37SDmitry Karpeev       ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr);
1252e7c4fc90SDmitry Karpeev       /* By default there are no DMs associated with subproblems. */
12530298fd71SBarry Smith       if (dmlist) *dmlist = NULL;
1254e7c4fc90SDmitry Karpeev     }
12558865f1eaSKarl Rupp   } else {
125616621825SDmitry Karpeev     ierr = (*dm->ops->createfielddecomposition)(dm,len,namelist,islist,dmlist);CHKERRQ(ierr);
125716621825SDmitry Karpeev   }
125816621825SDmitry Karpeev   PetscFunctionReturn(0);
125916621825SDmitry Karpeev }
126016621825SDmitry Karpeev 
126116621825SDmitry Karpeev #undef __FUNCT__
1262435a35e8SMatthew G Knepley #define __FUNCT__ "DMCreateSubDM"
1263435a35e8SMatthew G Knepley /*@C
1264435a35e8SMatthew G Knepley   DMCreateSubDM - Returns an IS and DM encapsulating a subproblem defined by the fields passed in.
1265435a35e8SMatthew G Knepley                   The fields are defined by DMCreateFieldIS().
1266435a35e8SMatthew G Knepley 
1267435a35e8SMatthew G Knepley   Not collective
1268435a35e8SMatthew G Knepley 
1269435a35e8SMatthew G Knepley   Input Parameters:
1270435a35e8SMatthew G Knepley + dm - the DM object
1271435a35e8SMatthew G Knepley . numFields - number of fields in this subproblem
12720298fd71SBarry Smith - len       - The number of subproblems in the decomposition (or NULL if not requested)
1273435a35e8SMatthew G Knepley 
1274435a35e8SMatthew G Knepley   Output Parameters:
1275435a35e8SMatthew G Knepley . is - The global indices for the subproblem
1276435a35e8SMatthew G Knepley - dm - The DM for the subproblem
1277435a35e8SMatthew G Knepley 
1278435a35e8SMatthew G Knepley   Level: intermediate
1279435a35e8SMatthew G Knepley 
1280435a35e8SMatthew G Knepley .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1281435a35e8SMatthew G Knepley @*/
1282435a35e8SMatthew G Knepley PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1283435a35e8SMatthew G Knepley {
1284435a35e8SMatthew G Knepley   PetscErrorCode ierr;
1285435a35e8SMatthew G Knepley 
1286435a35e8SMatthew G Knepley   PetscFunctionBegin;
1287435a35e8SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1288435a35e8SMatthew G Knepley   PetscValidPointer(fields,3);
12898865f1eaSKarl Rupp   if (is) PetscValidPointer(is,4);
12908865f1eaSKarl Rupp   if (subdm) PetscValidPointer(subdm,5);
1291435a35e8SMatthew G Knepley   if (dm->ops->createsubdm) {
1292435a35e8SMatthew G Knepley     ierr = (*dm->ops->createsubdm)(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
129382f516ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateSubDM implementation defined");
1294435a35e8SMatthew G Knepley   PetscFunctionReturn(0);
1295435a35e8SMatthew G Knepley }
1296435a35e8SMatthew G Knepley 
129716621825SDmitry Karpeev 
129816621825SDmitry Karpeev #undef __FUNCT__
129916621825SDmitry Karpeev #define __FUNCT__ "DMCreateDomainDecomposition"
130016621825SDmitry Karpeev /*@C
13018d4ac253SDmitry Karpeev   DMCreateDomainDecomposition - Returns lists of IS objects defining a decomposition of a problem into subproblems
13028d4ac253SDmitry Karpeev                           corresponding to restrictions to pairs nested subdomains: each IS contains the global
13038d4ac253SDmitry Karpeev                           indices of the dofs of the corresponding subdomains.  The inner subdomains conceptually
13048d4ac253SDmitry Karpeev                           define a nonoverlapping covering, while outer subdomains can overlap.
13058d4ac253SDmitry Karpeev                           The optional list of DMs define the DM for each subproblem.
130616621825SDmitry Karpeev 
130716621825SDmitry Karpeev   Not collective
130816621825SDmitry Karpeev 
130916621825SDmitry Karpeev   Input Parameter:
131016621825SDmitry Karpeev . dm - the DM object
131116621825SDmitry Karpeev 
131216621825SDmitry Karpeev   Output Parameters:
13130298fd71SBarry Smith + len         - The number of subproblems in the domain decomposition (or NULL if not requested)
13140298fd71SBarry Smith . namelist    - The name for each subdomain (or NULL if not requested)
13150298fd71SBarry Smith . innerislist - The global indices for each inner subdomain (or NULL, if not requested)
13160298fd71SBarry Smith . outerislist - The global indices for each outer subdomain (or NULL, if not requested)
13170298fd71SBarry Smith - dmlist      - The DMs for each subdomain subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined)
131816621825SDmitry Karpeev 
131916621825SDmitry Karpeev   Level: intermediate
132016621825SDmitry Karpeev 
132116621825SDmitry Karpeev   Notes:
132216621825SDmitry Karpeev   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
132316621825SDmitry Karpeev   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
132416621825SDmitry Karpeev   and all of the arrays should be freed with PetscFree().
132516621825SDmitry Karpeev 
13268d4ac253SDmitry Karpeev .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateDomainDecompositionDM(), DMCreateFieldDecomposition()
132716621825SDmitry Karpeev @*/
13288d4ac253SDmitry Karpeev PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist)
132916621825SDmitry Karpeev {
133016621825SDmitry Karpeev   PetscErrorCode      ierr;
1331be081cd6SPeter Brune   DMSubDomainHookLink link;
1332be081cd6SPeter Brune   PetscInt            i,l;
133316621825SDmitry Karpeev 
133416621825SDmitry Karpeev   PetscFunctionBegin;
133516621825SDmitry Karpeev   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
133614a18fd3SPeter Brune   if (len)           {PetscValidPointer(len,2);            *len         = 0;}
13370298fd71SBarry Smith   if (namelist)      {PetscValidPointer(namelist,3);       *namelist    = NULL;}
13380298fd71SBarry Smith   if (innerislist)   {PetscValidPointer(innerislist,4);    *innerislist = NULL;}
13390298fd71SBarry Smith   if (outerislist)   {PetscValidPointer(outerislist,5);    *outerislist = NULL;}
13400298fd71SBarry Smith   if (dmlist)        {PetscValidPointer(dmlist,6);         *dmlist      = NULL;}
1341f3f0edfdSDmitry Karpeev   /*
1342f3f0edfdSDmitry Karpeev    Is it a good idea to apply the following check across all impls?
1343f3f0edfdSDmitry Karpeev    Perhaps some impls can have a well-defined decomposition before DMSetUp?
1344f3f0edfdSDmitry Karpeev    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
1345f3f0edfdSDmitry Karpeev    */
1346ce94432eSBarry Smith   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
134716621825SDmitry Karpeev   if (dm->ops->createdomaindecomposition) {
1348be081cd6SPeter Brune     ierr = (*dm->ops->createdomaindecomposition)(dm,&l,namelist,innerislist,outerislist,dmlist);CHKERRQ(ierr);
134914a18fd3SPeter Brune     /* copy subdomain hooks and context over to the subdomain DMs */
135014a18fd3SPeter Brune     if (dmlist) {
1351be081cd6SPeter Brune       for (i = 0; i < l; i++) {
1352be081cd6SPeter Brune         for (link=dm->subdomainhook; link; link=link->next) {
1353be081cd6SPeter Brune           if (link->ddhook) {ierr = (*link->ddhook)(dm,(*dmlist)[i],link->ctx);CHKERRQ(ierr);}
1354be081cd6SPeter Brune         }
1355d425c654SPeter Brune         (*dmlist)[i]->ctx = dm->ctx;
1356e7c4fc90SDmitry Karpeev       }
135714a18fd3SPeter Brune     }
135814a18fd3SPeter Brune     if (len) *len = l;
135914a18fd3SPeter Brune   }
1360e30e807fSPeter Brune   PetscFunctionReturn(0);
1361e30e807fSPeter Brune }
1362e30e807fSPeter Brune 
1363e30e807fSPeter Brune 
1364e30e807fSPeter Brune #undef __FUNCT__
1365e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters"
1366e30e807fSPeter Brune /*@C
1367e30e807fSPeter Brune   DMCreateDomainDecompositionScatters - Returns scatters to the subdomain vectors from the global vector
1368e30e807fSPeter Brune 
1369e30e807fSPeter Brune   Not collective
1370e30e807fSPeter Brune 
1371e30e807fSPeter Brune   Input Parameters:
1372e30e807fSPeter Brune + dm - the DM object
1373e30e807fSPeter Brune . n  - the number of subdomain scatters
1374e30e807fSPeter Brune - subdms - the local subdomains
1375e30e807fSPeter Brune 
1376e30e807fSPeter Brune   Output Parameters:
1377e30e807fSPeter Brune + n     - the number of scatters returned
1378e30e807fSPeter Brune . iscat - scatter from global vector to nonoverlapping global vector entries on subdomain
1379e30e807fSPeter Brune . oscat - scatter from global vector to overlapping global vector entries on subdomain
1380e30e807fSPeter Brune - gscat - scatter from global vector to local vector on subdomain (fills in ghosts)
1381e30e807fSPeter Brune 
1382e30e807fSPeter Brune   Notes: This is an alternative to the iis and ois arguments in DMCreateDomainDecomposition that allow for the solution
1383e30e807fSPeter Brune   of general nonlinear problems with overlapping subdomain methods.  While merely having index sets that enable subsets
1384e30e807fSPeter Brune   of the residual equations to be created is fine for linear problems, nonlinear problems require local assembly of
1385e30e807fSPeter Brune   solution and residual data.
1386e30e807fSPeter Brune 
1387e30e807fSPeter Brune   Level: developer
1388e30e807fSPeter Brune 
1389e30e807fSPeter Brune .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1390e30e807fSPeter Brune @*/
1391e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecompositionScatters(DM dm,PetscInt n,DM *subdms,VecScatter **iscat,VecScatter **oscat,VecScatter **gscat)
1392e30e807fSPeter Brune {
1393e30e807fSPeter Brune   PetscErrorCode ierr;
1394e30e807fSPeter Brune 
1395e30e807fSPeter Brune   PetscFunctionBegin;
1396e30e807fSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1397e30e807fSPeter Brune   PetscValidPointer(subdms,3);
1398e30e807fSPeter Brune   if (dm->ops->createddscatters) {
1399e30e807fSPeter Brune     ierr = (*dm->ops->createddscatters)(dm,n,subdms,iscat,oscat,gscat);CHKERRQ(ierr);
140082f516ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateDomainDecompositionLocalScatter implementation defined");
1401e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1402e7c4fc90SDmitry Karpeev }
1403e7c4fc90SDmitry Karpeev 
1404731c8d9eSDmitry Karpeev #undef __FUNCT__
140547c6ae99SBarry Smith #define __FUNCT__ "DMRefine"
140647c6ae99SBarry Smith /*@
140747c6ae99SBarry Smith   DMRefine - Refines a DM object
140847c6ae99SBarry Smith 
140947c6ae99SBarry Smith   Collective on DM
141047c6ae99SBarry Smith 
141147c6ae99SBarry Smith   Input Parameter:
141247c6ae99SBarry Smith + dm   - the DM object
141391d95f02SJed Brown - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
141447c6ae99SBarry Smith 
141547c6ae99SBarry Smith   Output Parameter:
14160298fd71SBarry Smith . dmf - the refined DM, or NULL
1417ae0a1c52SMatthew G Knepley 
14180298fd71SBarry Smith   Note: If no refinement was done, the return value is NULL
141947c6ae99SBarry Smith 
142047c6ae99SBarry Smith   Level: developer
142147c6ae99SBarry Smith 
1422e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
142347c6ae99SBarry Smith @*/
14247087cfbeSBarry Smith PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
142547c6ae99SBarry Smith {
142647c6ae99SBarry Smith   PetscErrorCode   ierr;
1427c833c3b5SJed Brown   DMRefineHookLink link;
142847c6ae99SBarry Smith 
142947c6ae99SBarry Smith   PetscFunctionBegin;
1430732e2eb9SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
143147c6ae99SBarry Smith   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
14324057135bSMatthew G Knepley   if (*dmf) {
143343842a1eSJed Brown     (*dmf)->ops->creatematrix = dm->ops->creatematrix;
14348865f1eaSKarl Rupp 
14358cd211a4SJed Brown     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
14368865f1eaSKarl Rupp 
1437644e2e5bSBarry Smith     (*dmf)->ctx       = dm->ctx;
14380598a293SJed Brown     (*dmf)->leveldown = dm->leveldown;
1439656b349aSBarry Smith     (*dmf)->levelup   = dm->levelup + 1;
14408865f1eaSKarl Rupp 
1441e4b4b23bSJed Brown     ierr = DMSetMatType(*dmf,dm->mattype);CHKERRQ(ierr);
1442c833c3b5SJed Brown     for (link=dm->refinehook; link; link=link->next) {
14438865f1eaSKarl Rupp       if (link->refinehook) {
14448865f1eaSKarl Rupp         ierr = (*link->refinehook)(dm,*dmf,link->ctx);CHKERRQ(ierr);
14458865f1eaSKarl Rupp       }
1446c833c3b5SJed Brown     }
1447c833c3b5SJed Brown   }
1448c833c3b5SJed Brown   PetscFunctionReturn(0);
1449c833c3b5SJed Brown }
1450c833c3b5SJed Brown 
1451c833c3b5SJed Brown #undef __FUNCT__
1452c833c3b5SJed Brown #define __FUNCT__ "DMRefineHookAdd"
1453bb9467b5SJed Brown /*@C
1454c833c3b5SJed Brown    DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid
1455c833c3b5SJed Brown 
1456c833c3b5SJed Brown    Logically Collective
1457c833c3b5SJed Brown 
1458c833c3b5SJed Brown    Input Arguments:
1459c833c3b5SJed Brown +  coarse - nonlinear solver context on which to run a hook when restricting to a coarser level
1460c833c3b5SJed Brown .  refinehook - function to run when setting up a coarser level
1461c833c3b5SJed Brown .  interphook - function to run to update data on finer levels (once per SNESSolve())
14620298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1463c833c3b5SJed Brown 
1464c833c3b5SJed Brown    Calling sequence of refinehook:
1465c833c3b5SJed Brown $    refinehook(DM coarse,DM fine,void *ctx);
1466c833c3b5SJed Brown 
1467c833c3b5SJed Brown +  coarse - coarse level DM
1468c833c3b5SJed Brown .  fine - fine level DM to interpolate problem to
1469c833c3b5SJed Brown -  ctx - optional user-defined function context
1470c833c3b5SJed Brown 
1471c833c3b5SJed Brown    Calling sequence for interphook:
1472c833c3b5SJed Brown $    interphook(DM coarse,Mat interp,DM fine,void *ctx)
1473c833c3b5SJed Brown 
1474c833c3b5SJed Brown +  coarse - coarse level DM
1475c833c3b5SJed Brown .  interp - matrix interpolating a coarse-level solution to the finer grid
1476c833c3b5SJed Brown .  fine - fine level DM to update
1477c833c3b5SJed Brown -  ctx - optional user-defined function context
1478c833c3b5SJed Brown 
1479c833c3b5SJed Brown    Level: advanced
1480c833c3b5SJed Brown 
1481c833c3b5SJed Brown    Notes:
1482c833c3b5SJed Brown    This function is only needed if auxiliary data needs to be passed to fine grids while grid sequencing
1483c833c3b5SJed Brown 
1484c833c3b5SJed Brown    If this function is called multiple times, the hooks will be run in the order they are added.
1485c833c3b5SJed Brown 
1486bb9467b5SJed Brown    This function is currently not available from Fortran.
1487bb9467b5SJed Brown 
1488c833c3b5SJed Brown .seealso: DMCoarsenHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1489c833c3b5SJed Brown @*/
1490c833c3b5SJed Brown PetscErrorCode DMRefineHookAdd(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx)
1491c833c3b5SJed Brown {
1492c833c3b5SJed Brown   PetscErrorCode   ierr;
1493c833c3b5SJed Brown   DMRefineHookLink link,*p;
1494c833c3b5SJed Brown 
1495c833c3b5SJed Brown   PetscFunctionBegin;
1496c833c3b5SJed Brown   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
1497c833c3b5SJed Brown   for (p=&coarse->refinehook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1498c833c3b5SJed Brown   ierr             = PetscMalloc(sizeof(struct _DMRefineHookLink),&link);CHKERRQ(ierr);
1499c833c3b5SJed Brown   link->refinehook = refinehook;
1500c833c3b5SJed Brown   link->interphook = interphook;
1501c833c3b5SJed Brown   link->ctx        = ctx;
15020298fd71SBarry Smith   link->next       = NULL;
1503c833c3b5SJed Brown   *p               = link;
1504c833c3b5SJed Brown   PetscFunctionReturn(0);
1505c833c3b5SJed Brown }
1506c833c3b5SJed Brown 
1507c833c3b5SJed Brown #undef __FUNCT__
1508c833c3b5SJed Brown #define __FUNCT__ "DMInterpolate"
1509c833c3b5SJed Brown /*@
1510c833c3b5SJed Brown    DMInterpolate - interpolates user-defined problem data to a finer DM by running hooks registered by DMRefineHookAdd()
1511c833c3b5SJed Brown 
1512c833c3b5SJed Brown    Collective if any hooks are
1513c833c3b5SJed Brown 
1514c833c3b5SJed Brown    Input Arguments:
1515c833c3b5SJed Brown +  coarse - coarser DM to use as a base
1516c833c3b5SJed Brown .  restrct - interpolation matrix, apply using MatInterpolate()
1517c833c3b5SJed Brown -  fine - finer DM to update
1518c833c3b5SJed Brown 
1519c833c3b5SJed Brown    Level: developer
1520c833c3b5SJed Brown 
1521c833c3b5SJed Brown .seealso: DMRefineHookAdd(), MatInterpolate()
1522c833c3b5SJed Brown @*/
1523c833c3b5SJed Brown PetscErrorCode DMInterpolate(DM coarse,Mat interp,DM fine)
1524c833c3b5SJed Brown {
1525c833c3b5SJed Brown   PetscErrorCode   ierr;
1526c833c3b5SJed Brown   DMRefineHookLink link;
1527c833c3b5SJed Brown 
1528c833c3b5SJed Brown   PetscFunctionBegin;
1529c833c3b5SJed Brown   for (link=fine->refinehook; link; link=link->next) {
15308865f1eaSKarl Rupp     if (link->interphook) {
15318865f1eaSKarl Rupp       ierr = (*link->interphook)(coarse,interp,fine,link->ctx);CHKERRQ(ierr);
15328865f1eaSKarl Rupp     }
15334057135bSMatthew G Knepley   }
153447c6ae99SBarry Smith   PetscFunctionReturn(0);
153547c6ae99SBarry Smith }
153647c6ae99SBarry Smith 
153747c6ae99SBarry Smith #undef __FUNCT__
1538eb3f98d2SBarry Smith #define __FUNCT__ "DMGetRefineLevel"
1539eb3f98d2SBarry Smith /*@
1540eb3f98d2SBarry Smith     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
1541eb3f98d2SBarry Smith 
1542eb3f98d2SBarry Smith     Not Collective
1543eb3f98d2SBarry Smith 
1544eb3f98d2SBarry Smith     Input Parameter:
1545eb3f98d2SBarry Smith .   dm - the DM object
1546eb3f98d2SBarry Smith 
1547eb3f98d2SBarry Smith     Output Parameter:
1548eb3f98d2SBarry Smith .   level - number of refinements
1549eb3f98d2SBarry Smith 
1550eb3f98d2SBarry Smith     Level: developer
1551eb3f98d2SBarry Smith 
15526a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1553eb3f98d2SBarry Smith 
1554eb3f98d2SBarry Smith @*/
1555eb3f98d2SBarry Smith PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
1556eb3f98d2SBarry Smith {
1557eb3f98d2SBarry Smith   PetscFunctionBegin;
1558eb3f98d2SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1559eb3f98d2SBarry Smith   *level = dm->levelup;
1560eb3f98d2SBarry Smith   PetscFunctionReturn(0);
1561eb3f98d2SBarry Smith }
1562eb3f98d2SBarry Smith 
1563eb3f98d2SBarry Smith #undef __FUNCT__
1564baf369e7SPeter Brune #define __FUNCT__ "DMGlobalToLocalHookAdd"
1565bb9467b5SJed Brown /*@C
1566baf369e7SPeter Brune    DMGlobalToLocalHookAdd - adds a callback to be run when global to local is called
1567baf369e7SPeter Brune 
1568baf369e7SPeter Brune    Logically Collective
1569baf369e7SPeter Brune 
1570baf369e7SPeter Brune    Input Arguments:
1571baf369e7SPeter Brune +  dm - the DM
1572baf369e7SPeter Brune .  beginhook - function to run at the beginning of DMGlobalToLocalBegin()
1573baf369e7SPeter Brune .  endhook - function to run after DMGlobalToLocalEnd() has completed
15740298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1575baf369e7SPeter Brune 
1576baf369e7SPeter Brune    Calling sequence for beginhook:
1577baf369e7SPeter Brune $    beginhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx)
1578baf369e7SPeter Brune 
1579baf369e7SPeter Brune +  dm - global DM
1580baf369e7SPeter Brune .  g - global vector
1581baf369e7SPeter Brune .  mode - mode
1582baf369e7SPeter Brune .  l - local vector
1583baf369e7SPeter Brune -  ctx - optional user-defined function context
1584baf369e7SPeter Brune 
1585baf369e7SPeter Brune 
1586baf369e7SPeter Brune    Calling sequence for endhook:
1587ec4806b8SPeter Brune $    endhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx)
1588baf369e7SPeter Brune 
1589baf369e7SPeter Brune +  global - global DM
1590baf369e7SPeter Brune -  ctx - optional user-defined function context
1591baf369e7SPeter Brune 
1592baf369e7SPeter Brune    Level: advanced
1593baf369e7SPeter Brune 
1594baf369e7SPeter Brune .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1595baf369e7SPeter Brune @*/
1596baf369e7SPeter Brune PetscErrorCode DMGlobalToLocalHookAdd(DM dm,PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*),void *ctx)
1597baf369e7SPeter Brune {
1598baf369e7SPeter Brune   PetscErrorCode          ierr;
1599baf369e7SPeter Brune   DMGlobalToLocalHookLink link,*p;
1600baf369e7SPeter Brune 
1601baf369e7SPeter Brune   PetscFunctionBegin;
1602baf369e7SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1603baf369e7SPeter Brune   for (p=&dm->gtolhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1604baf369e7SPeter Brune   ierr            = PetscMalloc(sizeof(struct _DMGlobalToLocalHookLink),&link);CHKERRQ(ierr);
1605baf369e7SPeter Brune   link->beginhook = beginhook;
1606baf369e7SPeter Brune   link->endhook   = endhook;
1607baf369e7SPeter Brune   link->ctx       = ctx;
16080298fd71SBarry Smith   link->next      = NULL;
1609baf369e7SPeter Brune   *p              = link;
1610baf369e7SPeter Brune   PetscFunctionReturn(0);
1611baf369e7SPeter Brune }
1612baf369e7SPeter Brune 
1613baf369e7SPeter Brune #undef __FUNCT__
161447c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin"
161547c6ae99SBarry Smith /*@
161647c6ae99SBarry Smith     DMGlobalToLocalBegin - Begins updating local vectors from global vector
161747c6ae99SBarry Smith 
161847c6ae99SBarry Smith     Neighbor-wise Collective on DM
161947c6ae99SBarry Smith 
162047c6ae99SBarry Smith     Input Parameters:
162147c6ae99SBarry Smith +   dm - the DM object
162247c6ae99SBarry Smith .   g - the global vector
162347c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
162447c6ae99SBarry Smith -   l - the local vector
162547c6ae99SBarry Smith 
162647c6ae99SBarry Smith 
162747c6ae99SBarry Smith     Level: beginner
162847c6ae99SBarry Smith 
1629e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
163047c6ae99SBarry Smith 
163147c6ae99SBarry Smith @*/
16327087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
163347c6ae99SBarry Smith {
16347128ae9fSMatthew G Knepley   PetscSF                 sf;
163547c6ae99SBarry Smith   PetscErrorCode          ierr;
1636baf369e7SPeter Brune   DMGlobalToLocalHookLink link;
163747c6ae99SBarry Smith 
163847c6ae99SBarry Smith   PetscFunctionBegin;
1639171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1640baf369e7SPeter Brune   for (link=dm->gtolhook; link; link=link->next) {
16418865f1eaSKarl Rupp     if (link->beginhook) {
16428865f1eaSKarl Rupp       ierr = (*link->beginhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);
16438865f1eaSKarl Rupp     }
1644baf369e7SPeter Brune   }
16457128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
16467128ae9fSMatthew G Knepley   if (sf) {
16477128ae9fSMatthew G Knepley     PetscScalar *lArray, *gArray;
16487128ae9fSMatthew G Knepley 
164982f516ccSBarry Smith     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
16507128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
16517128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
16527128ae9fSMatthew G Knepley     ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
16537128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
16547128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
16557128ae9fSMatthew G Knepley   } else {
1656843c4018SMatthew G Knepley     ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
16577128ae9fSMatthew G Knepley   }
165847c6ae99SBarry Smith   PetscFunctionReturn(0);
165947c6ae99SBarry Smith }
166047c6ae99SBarry Smith 
166147c6ae99SBarry Smith #undef __FUNCT__
166247c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd"
166347c6ae99SBarry Smith /*@
166447c6ae99SBarry Smith     DMGlobalToLocalEnd - Ends updating local vectors from global vector
166547c6ae99SBarry Smith 
166647c6ae99SBarry Smith     Neighbor-wise Collective on DM
166747c6ae99SBarry Smith 
166847c6ae99SBarry Smith     Input Parameters:
166947c6ae99SBarry Smith +   dm - the DM object
167047c6ae99SBarry Smith .   g - the global vector
167147c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
167247c6ae99SBarry Smith -   l - the local vector
167347c6ae99SBarry Smith 
167447c6ae99SBarry Smith 
167547c6ae99SBarry Smith     Level: beginner
167647c6ae99SBarry Smith 
1677e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
167847c6ae99SBarry Smith 
167947c6ae99SBarry Smith @*/
16807087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
168147c6ae99SBarry Smith {
16827128ae9fSMatthew G Knepley   PetscSF                 sf;
168347c6ae99SBarry Smith   PetscErrorCode          ierr;
168461a3c1faSSatish Balay   PetscScalar             *lArray, *gArray;
1685baf369e7SPeter Brune   DMGlobalToLocalHookLink link;
168647c6ae99SBarry Smith 
168747c6ae99SBarry Smith   PetscFunctionBegin;
1688171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
16897128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
16907128ae9fSMatthew G Knepley   if (sf) {
169182f516ccSBarry Smith     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
16927128ae9fSMatthew G Knepley 
16937128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
16947128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
16957128ae9fSMatthew G Knepley     ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
16967128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
16977128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
16987128ae9fSMatthew G Knepley   } else {
1699843c4018SMatthew G Knepley     ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
17007128ae9fSMatthew G Knepley   }
1701baf369e7SPeter Brune   for (link=dm->gtolhook; link; link=link->next) {
1702baf369e7SPeter Brune     if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);}
1703baf369e7SPeter Brune   }
170447c6ae99SBarry Smith   PetscFunctionReturn(0);
170547c6ae99SBarry Smith }
170647c6ae99SBarry Smith 
170747c6ae99SBarry Smith #undef __FUNCT__
17089a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalBegin"
170947c6ae99SBarry Smith /*@
17109a42bb27SBarry Smith     DMLocalToGlobalBegin - updates global vectors from local vectors
17119a42bb27SBarry Smith 
17129a42bb27SBarry Smith     Neighbor-wise Collective on DM
17139a42bb27SBarry Smith 
17149a42bb27SBarry Smith     Input Parameters:
17159a42bb27SBarry Smith +   dm - the DM object
1716f6813fd5SJed Brown .   l - the local vector
17179a42bb27SBarry 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
17189a42bb27SBarry Smith            base point.
1719f6813fd5SJed Brown - - the global vector
17209a42bb27SBarry Smith 
17219a42bb27SBarry 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
17229a42bb27SBarry 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
17239a42bb27SBarry Smith            global array to the final global array with VecAXPY().
17249a42bb27SBarry Smith 
17259a42bb27SBarry Smith     Level: beginner
17269a42bb27SBarry Smith 
1727e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
17289a42bb27SBarry Smith 
17299a42bb27SBarry Smith @*/
17307087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
17319a42bb27SBarry Smith {
17327128ae9fSMatthew G Knepley   PetscSF        sf;
17339a42bb27SBarry Smith   PetscErrorCode ierr;
17349a42bb27SBarry Smith 
17359a42bb27SBarry Smith   PetscFunctionBegin;
1736171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
17377128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
17387128ae9fSMatthew G Knepley   if (sf) {
17397128ae9fSMatthew G Knepley     MPI_Op      op;
17407128ae9fSMatthew G Knepley     PetscScalar *lArray, *gArray;
17417128ae9fSMatthew G Knepley 
17427128ae9fSMatthew G Knepley     switch (mode) {
17437128ae9fSMatthew G Knepley     case INSERT_VALUES:
17447128ae9fSMatthew G Knepley     case INSERT_ALL_VALUES:
17458bfbc91cSJed Brown       op = MPIU_REPLACE; break;
17467128ae9fSMatthew G Knepley     case ADD_VALUES:
17477128ae9fSMatthew G Knepley     case ADD_ALL_VALUES:
17487128ae9fSMatthew G Knepley       op = MPI_SUM; break;
17497128ae9fSMatthew G Knepley     default:
175082f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
17517128ae9fSMatthew G Knepley     }
17527128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
17537128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
17547128ae9fSMatthew G Knepley     ierr = PetscSFReduceBegin(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr);
17557128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
17567128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
17577128ae9fSMatthew G Knepley   } else {
1758843c4018SMatthew G Knepley     ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
17597128ae9fSMatthew G Knepley   }
17609a42bb27SBarry Smith   PetscFunctionReturn(0);
17619a42bb27SBarry Smith }
17629a42bb27SBarry Smith 
17639a42bb27SBarry Smith #undef __FUNCT__
17649a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalEnd"
17659a42bb27SBarry Smith /*@
17669a42bb27SBarry Smith     DMLocalToGlobalEnd - updates global vectors from local vectors
176747c6ae99SBarry Smith 
176847c6ae99SBarry Smith     Neighbor-wise Collective on DM
176947c6ae99SBarry Smith 
177047c6ae99SBarry Smith     Input Parameters:
177147c6ae99SBarry Smith +   dm - the DM object
1772f6813fd5SJed Brown .   l - the local vector
177347c6ae99SBarry Smith .   mode - INSERT_VALUES or ADD_VALUES
1774f6813fd5SJed Brown -   g - the global vector
177547c6ae99SBarry Smith 
177647c6ae99SBarry Smith 
177747c6ae99SBarry Smith     Level: beginner
177847c6ae99SBarry Smith 
1779e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
178047c6ae99SBarry Smith 
178147c6ae99SBarry Smith @*/
17827087cfbeSBarry Smith PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
178347c6ae99SBarry Smith {
17847128ae9fSMatthew G Knepley   PetscSF        sf;
178547c6ae99SBarry Smith   PetscErrorCode ierr;
178647c6ae99SBarry Smith 
178747c6ae99SBarry Smith   PetscFunctionBegin;
1788171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
17897128ae9fSMatthew G Knepley   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
17907128ae9fSMatthew G Knepley   if (sf) {
17917128ae9fSMatthew G Knepley     MPI_Op      op;
17927128ae9fSMatthew G Knepley     PetscScalar *lArray, *gArray;
17937128ae9fSMatthew G Knepley 
17947128ae9fSMatthew G Knepley     switch (mode) {
17957128ae9fSMatthew G Knepley     case INSERT_VALUES:
17967128ae9fSMatthew G Knepley     case INSERT_ALL_VALUES:
17978bfbc91cSJed Brown       op = MPIU_REPLACE; break;
17987128ae9fSMatthew G Knepley     case ADD_VALUES:
17997128ae9fSMatthew G Knepley     case ADD_ALL_VALUES:
18007128ae9fSMatthew G Knepley       op = MPI_SUM; break;
18017128ae9fSMatthew G Knepley     default:
180282f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
18037128ae9fSMatthew G Knepley     }
18047128ae9fSMatthew G Knepley     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
18057128ae9fSMatthew G Knepley     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
18067128ae9fSMatthew G Knepley     ierr = PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr);
18077128ae9fSMatthew G Knepley     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
18087128ae9fSMatthew G Knepley     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
18097128ae9fSMatthew G Knepley   } else {
1810843c4018SMatthew G Knepley     ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
18117128ae9fSMatthew G Knepley   }
181247c6ae99SBarry Smith   PetscFunctionReturn(0);
181347c6ae99SBarry Smith }
181447c6ae99SBarry Smith 
181547c6ae99SBarry Smith #undef __FUNCT__
1816f089877aSRichard Tran Mills #define __FUNCT__ "DMLocalToLocalBegin"
1817f089877aSRichard Tran Mills /*@
1818bc0a1609SRichard Tran Mills    DMLocalToLocalBegin - Maps from a local vector (including ghost points
1819bc0a1609SRichard Tran Mills    that contain irrelevant values) to another local vector where the ghost
1820d78e899eSRichard Tran Mills    points in the second are set correctly. Must be followed by DMLocalToLocalEnd().
1821f089877aSRichard Tran Mills 
1822bc0a1609SRichard Tran Mills    Neighbor-wise Collective on DM and Vec
1823f089877aSRichard Tran Mills 
1824f089877aSRichard Tran Mills    Input Parameters:
1825f089877aSRichard Tran Mills +  dm - the DM object
1826bc0a1609SRichard Tran Mills .  g - the original local vector
1827bc0a1609SRichard Tran Mills -  mode - one of INSERT_VALUES or ADD_VALUES
1828f089877aSRichard Tran Mills 
1829bc0a1609SRichard Tran Mills    Output Parameter:
1830bc0a1609SRichard Tran Mills .  l  - the local vector with correct ghost values
1831f089877aSRichard Tran Mills 
1832f089877aSRichard Tran Mills    Level: intermediate
1833f089877aSRichard Tran Mills 
1834bc0a1609SRichard Tran Mills    Notes:
1835bc0a1609SRichard Tran Mills    The local vectors used here need not be the same as those
1836bc0a1609SRichard Tran Mills    obtained from DMCreateLocalVector(), BUT they
1837bc0a1609SRichard Tran Mills    must have the same parallel data layout; they could, for example, be
1838bc0a1609SRichard Tran Mills    obtained with VecDuplicate() from the DM originating vectors.
1839bc0a1609SRichard Tran Mills 
1840bc0a1609SRichard Tran Mills .keywords: DM, local-to-local, begin
1841bc0a1609SRichard Tran Mills .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalEnd(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1842f089877aSRichard Tran Mills 
1843f089877aSRichard Tran Mills @*/
1844f089877aSRichard Tran Mills PetscErrorCode  DMLocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
1845f089877aSRichard Tran Mills {
1846f089877aSRichard Tran Mills   PetscErrorCode          ierr;
1847f089877aSRichard Tran Mills 
1848f089877aSRichard Tran Mills   PetscFunctionBegin;
1849f089877aSRichard Tran Mills   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1850f089877aSRichard Tran Mills   ierr = (*dm->ops->localtolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1851f089877aSRichard Tran Mills   PetscFunctionReturn(0);
1852f089877aSRichard Tran Mills }
1853f089877aSRichard Tran Mills 
1854f089877aSRichard Tran Mills #undef __FUNCT__
1855f089877aSRichard Tran Mills #define __FUNCT__ "DMLocalToLocalEnd"
1856f089877aSRichard Tran Mills /*@
1857bc0a1609SRichard Tran Mills    DMLocalToLocalEnd - Maps from a local vector (including ghost points
1858bc0a1609SRichard Tran Mills    that contain irrelevant values) to another local vector where the ghost
1859d78e899eSRichard Tran Mills    points in the second are set correctly. Must be preceded by DMLocalToLocalBegin().
1860f089877aSRichard Tran Mills 
1861bc0a1609SRichard Tran Mills    Neighbor-wise Collective on DM and Vec
1862f089877aSRichard Tran Mills 
1863f089877aSRichard Tran Mills    Input Parameters:
1864bc0a1609SRichard Tran Mills +  da - the DM object
1865bc0a1609SRichard Tran Mills .  g - the original local vector
1866bc0a1609SRichard Tran Mills -  mode - one of INSERT_VALUES or ADD_VALUES
1867f089877aSRichard Tran Mills 
1868bc0a1609SRichard Tran Mills    Output Parameter:
1869bc0a1609SRichard Tran Mills .  l  - the local vector with correct ghost values
1870f089877aSRichard Tran Mills 
1871f089877aSRichard Tran Mills    Level: intermediate
1872f089877aSRichard Tran Mills 
1873bc0a1609SRichard Tran Mills    Notes:
1874bc0a1609SRichard Tran Mills    The local vectors used here need not be the same as those
1875bc0a1609SRichard Tran Mills    obtained from DMCreateLocalVector(), BUT they
1876bc0a1609SRichard Tran Mills    must have the same parallel data layout; they could, for example, be
1877bc0a1609SRichard Tran Mills    obtained with VecDuplicate() from the DM originating vectors.
1878bc0a1609SRichard Tran Mills 
1879bc0a1609SRichard Tran Mills .keywords: DM, local-to-local, end
1880bc0a1609SRichard Tran Mills .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1881f089877aSRichard Tran Mills 
1882f089877aSRichard Tran Mills @*/
1883f089877aSRichard Tran Mills PetscErrorCode  DMLocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
1884f089877aSRichard Tran Mills {
1885f089877aSRichard Tran Mills   PetscErrorCode          ierr;
1886f089877aSRichard Tran Mills 
1887f089877aSRichard Tran Mills   PetscFunctionBegin;
1888f089877aSRichard Tran Mills   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1889f089877aSRichard Tran Mills   ierr = (*dm->ops->localtolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1890f089877aSRichard Tran Mills   PetscFunctionReturn(0);
1891f089877aSRichard Tran Mills }
1892f089877aSRichard Tran Mills 
1893f089877aSRichard Tran Mills 
1894f089877aSRichard Tran Mills #undef __FUNCT__
189547c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen"
189647c6ae99SBarry Smith /*@
189747c6ae99SBarry Smith     DMCoarsen - Coarsens a DM object
189847c6ae99SBarry Smith 
189947c6ae99SBarry Smith     Collective on DM
190047c6ae99SBarry Smith 
190147c6ae99SBarry Smith     Input Parameter:
190247c6ae99SBarry Smith +   dm - the DM object
190391d95f02SJed Brown -   comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
190447c6ae99SBarry Smith 
190547c6ae99SBarry Smith     Output Parameter:
190647c6ae99SBarry Smith .   dmc - the coarsened DM
190747c6ae99SBarry Smith 
190847c6ae99SBarry Smith     Level: developer
190947c6ae99SBarry Smith 
1910e727c939SJed Brown .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
191147c6ae99SBarry Smith 
191247c6ae99SBarry Smith @*/
19137087cfbeSBarry Smith PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
191447c6ae99SBarry Smith {
191547c6ae99SBarry Smith   PetscErrorCode    ierr;
1916b17ce1afSJed Brown   DMCoarsenHookLink link;
191747c6ae99SBarry Smith 
191847c6ae99SBarry Smith   PetscFunctionBegin;
1919171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
192047c6ae99SBarry Smith   ierr                      = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
192143842a1eSJed Brown   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
19228cd211a4SJed Brown   ierr                      = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1923644e2e5bSBarry Smith   (*dmc)->ctx               = dm->ctx;
19240598a293SJed Brown   (*dmc)->levelup           = dm->levelup;
1925656b349aSBarry Smith   (*dmc)->leveldown         = dm->leveldown + 1;
1926e4b4b23bSJed Brown   ierr                      = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr);
1927b17ce1afSJed Brown   for (link=dm->coarsenhook; link; link=link->next) {
1928b17ce1afSJed Brown     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1929b17ce1afSJed Brown   }
1930b17ce1afSJed Brown   PetscFunctionReturn(0);
1931b17ce1afSJed Brown }
1932b17ce1afSJed Brown 
1933b17ce1afSJed Brown #undef __FUNCT__
1934b17ce1afSJed Brown #define __FUNCT__ "DMCoarsenHookAdd"
1935bb9467b5SJed Brown /*@C
1936b17ce1afSJed Brown    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1937b17ce1afSJed Brown 
1938b17ce1afSJed Brown    Logically Collective
1939b17ce1afSJed Brown 
1940b17ce1afSJed Brown    Input Arguments:
1941b17ce1afSJed Brown +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1942b17ce1afSJed Brown .  coarsenhook - function to run when setting up a coarser level
1943b17ce1afSJed Brown .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
19440298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1945b17ce1afSJed Brown 
1946b17ce1afSJed Brown    Calling sequence of coarsenhook:
1947b17ce1afSJed Brown $    coarsenhook(DM fine,DM coarse,void *ctx);
1948b17ce1afSJed Brown 
1949b17ce1afSJed Brown +  fine - fine level DM
1950b17ce1afSJed Brown .  coarse - coarse level DM to restrict problem to
1951b17ce1afSJed Brown -  ctx - optional user-defined function context
1952b17ce1afSJed Brown 
1953b17ce1afSJed Brown    Calling sequence for restricthook:
1954c833c3b5SJed Brown $    restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx)
1955b17ce1afSJed Brown 
1956b17ce1afSJed Brown +  fine - fine level DM
1957b17ce1afSJed Brown .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1958c833c3b5SJed Brown .  rscale - scaling vector for restriction
1959c833c3b5SJed Brown .  inject - matrix restricting by injection
1960b17ce1afSJed Brown .  coarse - coarse level DM to update
1961b17ce1afSJed Brown -  ctx - optional user-defined function context
1962b17ce1afSJed Brown 
1963b17ce1afSJed Brown    Level: advanced
1964b17ce1afSJed Brown 
1965b17ce1afSJed Brown    Notes:
1966b17ce1afSJed Brown    This function is only needed if auxiliary data needs to be set up on coarse grids.
1967b17ce1afSJed Brown 
1968b17ce1afSJed Brown    If this function is called multiple times, the hooks will be run in the order they are added.
1969b17ce1afSJed Brown 
1970b17ce1afSJed Brown    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1971b17ce1afSJed Brown    extract the finest level information from its context (instead of from the SNES).
1972b17ce1afSJed Brown 
1973bb9467b5SJed Brown    This function is currently not available from Fortran.
1974bb9467b5SJed Brown 
1975c833c3b5SJed Brown .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1976b17ce1afSJed Brown @*/
1977b17ce1afSJed Brown PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1978b17ce1afSJed Brown {
1979b17ce1afSJed Brown   PetscErrorCode    ierr;
1980b17ce1afSJed Brown   DMCoarsenHookLink link,*p;
1981b17ce1afSJed Brown 
1982b17ce1afSJed Brown   PetscFunctionBegin;
1983b17ce1afSJed Brown   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
19846bfea28cSJed Brown   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1985b17ce1afSJed Brown   ierr               = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1986b17ce1afSJed Brown   link->coarsenhook  = coarsenhook;
1987b17ce1afSJed Brown   link->restricthook = restricthook;
1988b17ce1afSJed Brown   link->ctx          = ctx;
19890298fd71SBarry Smith   link->next         = NULL;
1990b17ce1afSJed Brown   *p                 = link;
1991b17ce1afSJed Brown   PetscFunctionReturn(0);
1992b17ce1afSJed Brown }
1993b17ce1afSJed Brown 
1994b17ce1afSJed Brown #undef __FUNCT__
1995b17ce1afSJed Brown #define __FUNCT__ "DMRestrict"
1996b17ce1afSJed Brown /*@
1997b17ce1afSJed Brown    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1998b17ce1afSJed Brown 
1999b17ce1afSJed Brown    Collective if any hooks are
2000b17ce1afSJed Brown 
2001b17ce1afSJed Brown    Input Arguments:
2002b17ce1afSJed Brown +  fine - finer DM to use as a base
2003b17ce1afSJed Brown .  restrct - restriction matrix, apply using MatRestrict()
2004b17ce1afSJed Brown .  inject - injection matrix, also use MatRestrict()
2005b17ce1afSJed Brown -  coarse - coarer DM to update
2006b17ce1afSJed Brown 
2007b17ce1afSJed Brown    Level: developer
2008b17ce1afSJed Brown 
2009b17ce1afSJed Brown .seealso: DMCoarsenHookAdd(), MatRestrict()
2010b17ce1afSJed Brown @*/
2011b17ce1afSJed Brown PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
2012b17ce1afSJed Brown {
2013b17ce1afSJed Brown   PetscErrorCode    ierr;
2014b17ce1afSJed Brown   DMCoarsenHookLink link;
2015b17ce1afSJed Brown 
2016b17ce1afSJed Brown   PetscFunctionBegin;
2017b17ce1afSJed Brown   for (link=fine->coarsenhook; link; link=link->next) {
20188865f1eaSKarl Rupp     if (link->restricthook) {
20198865f1eaSKarl Rupp       ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);
20208865f1eaSKarl Rupp     }
2021b17ce1afSJed Brown   }
202247c6ae99SBarry Smith   PetscFunctionReturn(0);
202347c6ae99SBarry Smith }
202447c6ae99SBarry Smith 
202547c6ae99SBarry Smith #undef __FUNCT__
2026be081cd6SPeter Brune #define __FUNCT__ "DMSubDomainHookAdd"
2027bb9467b5SJed Brown /*@C
2028be081cd6SPeter Brune    DMSubDomainHookAdd - adds a callback to be run when restricting a problem to the coarse grid
20295dbd56e3SPeter Brune 
20305dbd56e3SPeter Brune    Logically Collective
20315dbd56e3SPeter Brune 
20325dbd56e3SPeter Brune    Input Arguments:
20335dbd56e3SPeter Brune +  global - global DM
2034ec4806b8SPeter Brune .  ddhook - function to run to pass data to the decomposition DM upon its creation
20355dbd56e3SPeter Brune .  restricthook - function to run to update data on block solve (at the beginning of the block solve)
20360298fd71SBarry Smith -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
20375dbd56e3SPeter Brune 
2038ec4806b8SPeter Brune 
2039ec4806b8SPeter Brune    Calling sequence for ddhook:
2040ec4806b8SPeter Brune $    ddhook(DM global,DM block,void *ctx)
2041ec4806b8SPeter Brune 
2042ec4806b8SPeter Brune +  global - global DM
2043ec4806b8SPeter Brune .  block  - block DM
2044ec4806b8SPeter Brune -  ctx - optional user-defined function context
2045ec4806b8SPeter Brune 
20465dbd56e3SPeter Brune    Calling sequence for restricthook:
2047ec4806b8SPeter Brune $    restricthook(DM global,VecScatter out,VecScatter in,DM block,void *ctx)
20485dbd56e3SPeter Brune 
20495dbd56e3SPeter Brune +  global - global DM
20505dbd56e3SPeter Brune .  out    - scatter to the outer (with ghost and overlap points) block vector
20515dbd56e3SPeter Brune .  in     - scatter to block vector values only owned locally
2052ec4806b8SPeter Brune .  block  - block DM
20535dbd56e3SPeter Brune -  ctx - optional user-defined function context
20545dbd56e3SPeter Brune 
20555dbd56e3SPeter Brune    Level: advanced
20565dbd56e3SPeter Brune 
20575dbd56e3SPeter Brune    Notes:
2058ec4806b8SPeter Brune    This function is only needed if auxiliary data needs to be set up on subdomain DMs.
20595dbd56e3SPeter Brune 
20605dbd56e3SPeter Brune    If this function is called multiple times, the hooks will be run in the order they are added.
20615dbd56e3SPeter Brune 
20625dbd56e3SPeter Brune    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
2063ec4806b8SPeter Brune    extract the global information from its context (instead of from the SNES).
20645dbd56e3SPeter Brune 
2065bb9467b5SJed Brown    This function is currently not available from Fortran.
2066bb9467b5SJed Brown 
20675dbd56e3SPeter Brune .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
20685dbd56e3SPeter Brune @*/
2069be081cd6SPeter Brune PetscErrorCode DMSubDomainHookAdd(DM global,PetscErrorCode (*ddhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*),void *ctx)
20705dbd56e3SPeter Brune {
20715dbd56e3SPeter Brune   PetscErrorCode      ierr;
2072be081cd6SPeter Brune   DMSubDomainHookLink link,*p;
20735dbd56e3SPeter Brune 
20745dbd56e3SPeter Brune   PetscFunctionBegin;
20755dbd56e3SPeter Brune   PetscValidHeaderSpecific(global,DM_CLASSID,1);
2076be081cd6SPeter Brune   for (p=&global->subdomainhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
2077be081cd6SPeter Brune   ierr               = PetscMalloc(sizeof(struct _DMSubDomainHookLink),&link);CHKERRQ(ierr);
20785dbd56e3SPeter Brune   link->restricthook = restricthook;
2079be081cd6SPeter Brune   link->ddhook       = ddhook;
20805dbd56e3SPeter Brune   link->ctx          = ctx;
20810298fd71SBarry Smith   link->next         = NULL;
20825dbd56e3SPeter Brune   *p                 = link;
20835dbd56e3SPeter Brune   PetscFunctionReturn(0);
20845dbd56e3SPeter Brune }
20855dbd56e3SPeter Brune 
20865dbd56e3SPeter Brune #undef __FUNCT__
2087be081cd6SPeter Brune #define __FUNCT__ "DMSubDomainRestrict"
20885dbd56e3SPeter Brune /*@
2089be081cd6SPeter Brune    DMSubDomainRestrict - restricts user-defined problem data to a block DM by running hooks registered by DMSubDomainHookAdd()
20905dbd56e3SPeter Brune 
20915dbd56e3SPeter Brune    Collective if any hooks are
20925dbd56e3SPeter Brune 
20935dbd56e3SPeter Brune    Input Arguments:
20945dbd56e3SPeter Brune +  fine - finer DM to use as a base
2095be081cd6SPeter Brune .  oscatter - scatter from domain global vector filling subdomain global vector with overlap
2096be081cd6SPeter Brune .  gscatter - scatter from domain global vector filling subdomain local vector with ghosts
20975dbd56e3SPeter Brune -  coarse - coarer DM to update
20985dbd56e3SPeter Brune 
20995dbd56e3SPeter Brune    Level: developer
21005dbd56e3SPeter Brune 
21015dbd56e3SPeter Brune .seealso: DMCoarsenHookAdd(), MatRestrict()
21025dbd56e3SPeter Brune @*/
2103be081cd6SPeter Brune PetscErrorCode DMSubDomainRestrict(DM global,VecScatter oscatter,VecScatter gscatter,DM subdm)
21045dbd56e3SPeter Brune {
21055dbd56e3SPeter Brune   PetscErrorCode      ierr;
2106be081cd6SPeter Brune   DMSubDomainHookLink link;
21075dbd56e3SPeter Brune 
21085dbd56e3SPeter Brune   PetscFunctionBegin;
2109be081cd6SPeter Brune   for (link=global->subdomainhook; link; link=link->next) {
21108865f1eaSKarl Rupp     if (link->restricthook) {
21118865f1eaSKarl Rupp       ierr = (*link->restricthook)(global,oscatter,gscatter,subdm,link->ctx);CHKERRQ(ierr);
21128865f1eaSKarl Rupp     }
21135dbd56e3SPeter Brune   }
21145dbd56e3SPeter Brune   PetscFunctionReturn(0);
21155dbd56e3SPeter Brune }
21165dbd56e3SPeter Brune 
21175dbd56e3SPeter Brune #undef __FUNCT__
21185fe1f584SPeter Brune #define __FUNCT__ "DMGetCoarsenLevel"
21195fe1f584SPeter Brune /*@
21206a7d9d85SPeter Brune     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
21215fe1f584SPeter Brune 
21225fe1f584SPeter Brune     Not Collective
21235fe1f584SPeter Brune 
21245fe1f584SPeter Brune     Input Parameter:
21255fe1f584SPeter Brune .   dm - the DM object
21265fe1f584SPeter Brune 
21275fe1f584SPeter Brune     Output Parameter:
21286a7d9d85SPeter Brune .   level - number of coarsenings
21295fe1f584SPeter Brune 
21305fe1f584SPeter Brune     Level: developer
21315fe1f584SPeter Brune 
21326a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
21335fe1f584SPeter Brune 
21345fe1f584SPeter Brune @*/
21355fe1f584SPeter Brune PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
21365fe1f584SPeter Brune {
21375fe1f584SPeter Brune   PetscFunctionBegin;
21385fe1f584SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
21395fe1f584SPeter Brune   *level = dm->leveldown;
21405fe1f584SPeter Brune   PetscFunctionReturn(0);
21415fe1f584SPeter Brune }
21425fe1f584SPeter Brune 
21435fe1f584SPeter Brune 
21445fe1f584SPeter Brune 
21455fe1f584SPeter Brune #undef __FUNCT__
214647c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy"
214747c6ae99SBarry Smith /*@C
214847c6ae99SBarry Smith     DMRefineHierarchy - Refines a DM object, all levels at once
214947c6ae99SBarry Smith 
215047c6ae99SBarry Smith     Collective on DM
215147c6ae99SBarry Smith 
215247c6ae99SBarry Smith     Input Parameter:
215347c6ae99SBarry Smith +   dm - the DM object
215447c6ae99SBarry Smith -   nlevels - the number of levels of refinement
215547c6ae99SBarry Smith 
215647c6ae99SBarry Smith     Output Parameter:
215747c6ae99SBarry Smith .   dmf - the refined DM hierarchy
215847c6ae99SBarry Smith 
215947c6ae99SBarry Smith     Level: developer
216047c6ae99SBarry Smith 
2161e727c939SJed Brown .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
216247c6ae99SBarry Smith 
216347c6ae99SBarry Smith @*/
21647087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
216547c6ae99SBarry Smith {
216647c6ae99SBarry Smith   PetscErrorCode ierr;
216747c6ae99SBarry Smith 
216847c6ae99SBarry Smith   PetscFunctionBegin;
2169171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2170ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
217147c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
217247c6ae99SBarry Smith   if (dm->ops->refinehierarchy) {
217347c6ae99SBarry Smith     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
217447c6ae99SBarry Smith   } else if (dm->ops->refine) {
217547c6ae99SBarry Smith     PetscInt i;
217647c6ae99SBarry Smith 
2177ce94432eSBarry Smith     ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr);
217847c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
2179ce94432eSBarry Smith       ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr);
218047c6ae99SBarry Smith     }
2181ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
218247c6ae99SBarry Smith   PetscFunctionReturn(0);
218347c6ae99SBarry Smith }
218447c6ae99SBarry Smith 
218547c6ae99SBarry Smith #undef __FUNCT__
218647c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy"
218747c6ae99SBarry Smith /*@C
218847c6ae99SBarry Smith     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
218947c6ae99SBarry Smith 
219047c6ae99SBarry Smith     Collective on DM
219147c6ae99SBarry Smith 
219247c6ae99SBarry Smith     Input Parameter:
219347c6ae99SBarry Smith +   dm - the DM object
219447c6ae99SBarry Smith -   nlevels - the number of levels of coarsening
219547c6ae99SBarry Smith 
219647c6ae99SBarry Smith     Output Parameter:
219747c6ae99SBarry Smith .   dmc - the coarsened DM hierarchy
219847c6ae99SBarry Smith 
219947c6ae99SBarry Smith     Level: developer
220047c6ae99SBarry Smith 
2201e727c939SJed Brown .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
220247c6ae99SBarry Smith 
220347c6ae99SBarry Smith @*/
22047087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
220547c6ae99SBarry Smith {
220647c6ae99SBarry Smith   PetscErrorCode ierr;
220747c6ae99SBarry Smith 
220847c6ae99SBarry Smith   PetscFunctionBegin;
2209171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2210ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
221147c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
221247c6ae99SBarry Smith   PetscValidPointer(dmc,3);
221347c6ae99SBarry Smith   if (dm->ops->coarsenhierarchy) {
221447c6ae99SBarry Smith     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
221547c6ae99SBarry Smith   } else if (dm->ops->coarsen) {
221647c6ae99SBarry Smith     PetscInt i;
221747c6ae99SBarry Smith 
2218ce94432eSBarry Smith     ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr);
221947c6ae99SBarry Smith     for (i=1; i<nlevels; i++) {
2220ce94432eSBarry Smith       ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr);
222147c6ae99SBarry Smith     }
2222ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
222347c6ae99SBarry Smith   PetscFunctionReturn(0);
222447c6ae99SBarry Smith }
222547c6ae99SBarry Smith 
222647c6ae99SBarry Smith #undef __FUNCT__
2227e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates"
222847c6ae99SBarry Smith /*@
2229e727c939SJed Brown    DMCreateAggregates - Gets the aggregates that map between
223047c6ae99SBarry Smith    grids associated with two DMs.
223147c6ae99SBarry Smith 
223247c6ae99SBarry Smith    Collective on DM
223347c6ae99SBarry Smith 
223447c6ae99SBarry Smith    Input Parameters:
223547c6ae99SBarry Smith +  dmc - the coarse grid DM
223647c6ae99SBarry Smith -  dmf - the fine grid DM
223747c6ae99SBarry Smith 
223847c6ae99SBarry Smith    Output Parameters:
223947c6ae99SBarry Smith .  rest - the restriction matrix (transpose of the projection matrix)
224047c6ae99SBarry Smith 
224147c6ae99SBarry Smith    Level: intermediate
224247c6ae99SBarry Smith 
224347c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid
224447c6ae99SBarry Smith 
2245e727c939SJed Brown .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
224647c6ae99SBarry Smith @*/
2247e727c939SJed Brown PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
224847c6ae99SBarry Smith {
224947c6ae99SBarry Smith   PetscErrorCode ierr;
225047c6ae99SBarry Smith 
225147c6ae99SBarry Smith   PetscFunctionBegin;
2252171400e9SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
2253171400e9SBarry Smith   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
225447c6ae99SBarry Smith   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
225547c6ae99SBarry Smith   PetscFunctionReturn(0);
225647c6ae99SBarry Smith }
225747c6ae99SBarry Smith 
225847c6ae99SBarry Smith #undef __FUNCT__
22591a266240SBarry Smith #define __FUNCT__ "DMSetApplicationContextDestroy"
22601a266240SBarry Smith /*@C
22611a266240SBarry Smith     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
22621a266240SBarry Smith 
22631a266240SBarry Smith     Not Collective
22641a266240SBarry Smith 
22651a266240SBarry Smith     Input Parameters:
22661a266240SBarry Smith +   dm - the DM object
22671a266240SBarry Smith -   destroy - the destroy function
22681a266240SBarry Smith 
22691a266240SBarry Smith     Level: intermediate
22701a266240SBarry Smith 
2271e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
22721a266240SBarry Smith 
2273f07f9ceaSJed Brown @*/
22741a266240SBarry Smith PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
22751a266240SBarry Smith {
22761a266240SBarry Smith   PetscFunctionBegin;
2277171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
22781a266240SBarry Smith   dm->ctxdestroy = destroy;
22791a266240SBarry Smith   PetscFunctionReturn(0);
22801a266240SBarry Smith }
22811a266240SBarry Smith 
22821a266240SBarry Smith #undef __FUNCT__
22831b2093e4SBarry Smith #define __FUNCT__ "DMSetApplicationContext"
2284b07ff414SBarry Smith /*@
22851b2093e4SBarry Smith     DMSetApplicationContext - Set a user context into a DM object
228647c6ae99SBarry Smith 
228747c6ae99SBarry Smith     Not Collective
228847c6ae99SBarry Smith 
228947c6ae99SBarry Smith     Input Parameters:
229047c6ae99SBarry Smith +   dm - the DM object
229147c6ae99SBarry Smith -   ctx - the user context
229247c6ae99SBarry Smith 
229347c6ae99SBarry Smith     Level: intermediate
229447c6ae99SBarry Smith 
2295e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
229647c6ae99SBarry Smith 
229747c6ae99SBarry Smith @*/
22981b2093e4SBarry Smith PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
229947c6ae99SBarry Smith {
230047c6ae99SBarry Smith   PetscFunctionBegin;
2301171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
230247c6ae99SBarry Smith   dm->ctx = ctx;
230347c6ae99SBarry Smith   PetscFunctionReturn(0);
230447c6ae99SBarry Smith }
230547c6ae99SBarry Smith 
230647c6ae99SBarry Smith #undef __FUNCT__
23071b2093e4SBarry Smith #define __FUNCT__ "DMGetApplicationContext"
230847c6ae99SBarry Smith /*@
23091b2093e4SBarry Smith     DMGetApplicationContext - Gets a user context from a DM object
231047c6ae99SBarry Smith 
231147c6ae99SBarry Smith     Not Collective
231247c6ae99SBarry Smith 
231347c6ae99SBarry Smith     Input Parameter:
231447c6ae99SBarry Smith .   dm - the DM object
231547c6ae99SBarry Smith 
231647c6ae99SBarry Smith     Output Parameter:
231747c6ae99SBarry Smith .   ctx - the user context
231847c6ae99SBarry Smith 
231947c6ae99SBarry Smith     Level: intermediate
232047c6ae99SBarry Smith 
2321e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
232247c6ae99SBarry Smith 
232347c6ae99SBarry Smith @*/
23241b2093e4SBarry Smith PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
232547c6ae99SBarry Smith {
232647c6ae99SBarry Smith   PetscFunctionBegin;
2327171400e9SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
23281b2093e4SBarry Smith   *(void**)ctx = dm->ctx;
232947c6ae99SBarry Smith   PetscFunctionReturn(0);
233047c6ae99SBarry Smith }
233147c6ae99SBarry Smith 
233247c6ae99SBarry Smith #undef __FUNCT__
233308da532bSDmitry Karpeev #define __FUNCT__ "DMSetVariableBounds"
233408da532bSDmitry Karpeev /*@C
233508da532bSDmitry Karpeev     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
233608da532bSDmitry Karpeev 
233708da532bSDmitry Karpeev     Logically Collective on DM
233808da532bSDmitry Karpeev 
233908da532bSDmitry Karpeev     Input Parameter:
234008da532bSDmitry Karpeev +   dm - the DM object
23410298fd71SBarry Smith -   f - the function that computes variable bounds used by SNESVI (use NULL to cancel a previous function that was set)
234208da532bSDmitry Karpeev 
234308da532bSDmitry Karpeev     Level: intermediate
234408da532bSDmitry Karpeev 
2345835c3ec7SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(),
234608da532bSDmitry Karpeev          DMSetJacobian()
234708da532bSDmitry Karpeev 
234808da532bSDmitry Karpeev @*/
234908da532bSDmitry Karpeev PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
235008da532bSDmitry Karpeev {
235108da532bSDmitry Karpeev   PetscFunctionBegin;
235208da532bSDmitry Karpeev   dm->ops->computevariablebounds = f;
235308da532bSDmitry Karpeev   PetscFunctionReturn(0);
235408da532bSDmitry Karpeev }
235508da532bSDmitry Karpeev 
235608da532bSDmitry Karpeev #undef __FUNCT__
235708da532bSDmitry Karpeev #define __FUNCT__ "DMHasVariableBounds"
235808da532bSDmitry Karpeev /*@
235908da532bSDmitry Karpeev     DMHasVariableBounds - does the DM object have a variable bounds function?
236008da532bSDmitry Karpeev 
236108da532bSDmitry Karpeev     Not Collective
236208da532bSDmitry Karpeev 
236308da532bSDmitry Karpeev     Input Parameter:
236408da532bSDmitry Karpeev .   dm - the DM object to destroy
236508da532bSDmitry Karpeev 
236608da532bSDmitry Karpeev     Output Parameter:
236708da532bSDmitry Karpeev .   flg - PETSC_TRUE if the variable bounds function exists
236808da532bSDmitry Karpeev 
236908da532bSDmitry Karpeev     Level: developer
237008da532bSDmitry Karpeev 
237174e1e8c1SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
237208da532bSDmitry Karpeev 
237308da532bSDmitry Karpeev @*/
237408da532bSDmitry Karpeev PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
237508da532bSDmitry Karpeev {
237608da532bSDmitry Karpeev   PetscFunctionBegin;
237708da532bSDmitry Karpeev   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
237808da532bSDmitry Karpeev   PetscFunctionReturn(0);
237908da532bSDmitry Karpeev }
238008da532bSDmitry Karpeev 
238108da532bSDmitry Karpeev #undef __FUNCT__
238208da532bSDmitry Karpeev #define __FUNCT__ "DMComputeVariableBounds"
238308da532bSDmitry Karpeev /*@C
238408da532bSDmitry Karpeev     DMComputeVariableBounds - compute variable bounds used by SNESVI.
238508da532bSDmitry Karpeev 
238608da532bSDmitry Karpeev     Logically Collective on DM
238708da532bSDmitry Karpeev 
238808da532bSDmitry Karpeev     Input Parameters:
238908da532bSDmitry Karpeev +   dm - the DM object to destroy
239008da532bSDmitry Karpeev -   x  - current solution at which the bounds are computed
239108da532bSDmitry Karpeev 
239208da532bSDmitry Karpeev     Output parameters:
239308da532bSDmitry Karpeev +   xl - lower bound
239408da532bSDmitry Karpeev -   xu - upper bound
239508da532bSDmitry Karpeev 
239608da532bSDmitry Karpeev     Level: intermediate
239708da532bSDmitry Karpeev 
239874e1e8c1SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
239908da532bSDmitry Karpeev 
240008da532bSDmitry Karpeev @*/
240108da532bSDmitry Karpeev PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
240208da532bSDmitry Karpeev {
240308da532bSDmitry Karpeev   PetscErrorCode ierr;
24045fd66863SKarl Rupp 
240508da532bSDmitry Karpeev   PetscFunctionBegin;
240608da532bSDmitry Karpeev   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
240708da532bSDmitry Karpeev   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
240808da532bSDmitry Karpeev   if (dm->ops->computevariablebounds) {
240908da532bSDmitry Karpeev     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu);CHKERRQ(ierr);
24108865f1eaSKarl Rupp   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "This DM is incapable of computing variable bounds.");
241108da532bSDmitry Karpeev   PetscFunctionReturn(0);
241208da532bSDmitry Karpeev }
241308da532bSDmitry Karpeev 
241408da532bSDmitry Karpeev #undef __FUNCT__
2415b0ae01b7SPeter Brune #define __FUNCT__ "DMHasColoring"
2416b0ae01b7SPeter Brune /*@
2417b0ae01b7SPeter Brune     DMHasColoring - does the DM object have a method of providing a coloring?
2418b0ae01b7SPeter Brune 
2419b0ae01b7SPeter Brune     Not Collective
2420b0ae01b7SPeter Brune 
2421b0ae01b7SPeter Brune     Input Parameter:
2422b0ae01b7SPeter Brune .   dm - the DM object
2423b0ae01b7SPeter Brune 
2424b0ae01b7SPeter Brune     Output Parameter:
2425b0ae01b7SPeter Brune .   flg - PETSC_TRUE if the DM has facilities for DMCreateColoring().
2426b0ae01b7SPeter Brune 
2427b0ae01b7SPeter Brune     Level: developer
2428b0ae01b7SPeter Brune 
2429b0ae01b7SPeter Brune .seealso DMHasFunction(), DMCreateColoring()
2430b0ae01b7SPeter Brune 
2431b0ae01b7SPeter Brune @*/
2432b0ae01b7SPeter Brune PetscErrorCode  DMHasColoring(DM dm,PetscBool  *flg)
2433b0ae01b7SPeter Brune {
2434b0ae01b7SPeter Brune   PetscFunctionBegin;
2435b0ae01b7SPeter Brune   *flg =  (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE;
2436b0ae01b7SPeter Brune   PetscFunctionReturn(0);
2437b0ae01b7SPeter Brune }
2438b0ae01b7SPeter Brune 
2439b0ae01b7SPeter Brune #undef  __FUNCT__
244008da532bSDmitry Karpeev #define __FUNCT__ "DMSetVec"
2441748fac09SDmitry Karpeev /*@C
244208da532bSDmitry Karpeev     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
244308da532bSDmitry Karpeev 
244408da532bSDmitry Karpeev     Collective on DM
244508da532bSDmitry Karpeev 
244608da532bSDmitry Karpeev     Input Parameter:
244708da532bSDmitry Karpeev +   dm - the DM object
24480298fd71SBarry Smith -   x - location to compute residual and Jacobian, if NULL is passed to those routines; will be NULL for linear problems.
244908da532bSDmitry Karpeev 
245008da532bSDmitry Karpeev     Level: developer
245108da532bSDmitry Karpeev 
245274e1e8c1SBarry Smith .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
245308da532bSDmitry Karpeev 
245408da532bSDmitry Karpeev @*/
245508da532bSDmitry Karpeev PetscErrorCode  DMSetVec(DM dm,Vec x)
245608da532bSDmitry Karpeev {
245708da532bSDmitry Karpeev   PetscErrorCode ierr;
24585fd66863SKarl Rupp 
245908da532bSDmitry Karpeev   PetscFunctionBegin;
246008da532bSDmitry Karpeev   if (x) {
246108da532bSDmitry Karpeev     if (!dm->x) {
246208da532bSDmitry Karpeev       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
246308da532bSDmitry Karpeev     }
246408da532bSDmitry Karpeev     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
24658865f1eaSKarl Rupp   } else if (dm->x) {
246608da532bSDmitry Karpeev     ierr = VecDestroy(&dm->x);CHKERRQ(ierr);
246708da532bSDmitry Karpeev   }
246808da532bSDmitry Karpeev   PetscFunctionReturn(0);
246908da532bSDmitry Karpeev }
247008da532bSDmitry Karpeev 
24710298fd71SBarry Smith PetscFunctionList DMList              = NULL;
2472264ace61SBarry Smith PetscBool         DMRegisterAllCalled = PETSC_FALSE;
2473264ace61SBarry Smith 
2474264ace61SBarry Smith #undef __FUNCT__
2475264ace61SBarry Smith #define __FUNCT__ "DMSetType"
2476264ace61SBarry Smith /*@C
2477264ace61SBarry Smith   DMSetType - Builds a DM, for a particular DM implementation.
2478264ace61SBarry Smith 
2479264ace61SBarry Smith   Collective on DM
2480264ace61SBarry Smith 
2481264ace61SBarry Smith   Input Parameters:
2482264ace61SBarry Smith + dm     - The DM object
2483264ace61SBarry Smith - method - The name of the DM type
2484264ace61SBarry Smith 
2485264ace61SBarry Smith   Options Database Key:
2486264ace61SBarry Smith . -dm_type <type> - Sets the DM type; use -help for a list of available types
2487264ace61SBarry Smith 
2488264ace61SBarry Smith   Notes:
2489e1589f56SBarry Smith   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
2490264ace61SBarry Smith 
2491264ace61SBarry Smith   Level: intermediate
2492264ace61SBarry Smith 
2493264ace61SBarry Smith .keywords: DM, set, type
2494264ace61SBarry Smith .seealso: DMGetType(), DMCreate()
2495264ace61SBarry Smith @*/
249619fd82e9SBarry Smith PetscErrorCode  DMSetType(DM dm, DMType method)
2497264ace61SBarry Smith {
2498264ace61SBarry Smith   PetscErrorCode (*r)(DM);
2499264ace61SBarry Smith   PetscBool      match;
2500264ace61SBarry Smith   PetscErrorCode ierr;
2501264ace61SBarry Smith 
2502264ace61SBarry Smith   PetscFunctionBegin;
2503264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2504251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
2505264ace61SBarry Smith   if (match) PetscFunctionReturn(0);
2506264ace61SBarry Smith 
2507607a6623SBarry Smith   if (!DMRegisterAllCalled) {ierr = DMRegisterAll();CHKERRQ(ierr);}
25081c9cd337SJed Brown   ierr = PetscFunctionListFind(DMList,method,&r);CHKERRQ(ierr);
2509ce94432eSBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
2510264ace61SBarry Smith 
2511264ace61SBarry Smith   if (dm->ops->destroy) {
2512264ace61SBarry Smith     ierr             = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
25130298fd71SBarry Smith     dm->ops->destroy = NULL;
2514264ace61SBarry Smith   }
2515264ace61SBarry Smith   ierr = (*r)(dm);CHKERRQ(ierr);
2516264ace61SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
2517264ace61SBarry Smith   PetscFunctionReturn(0);
2518264ace61SBarry Smith }
2519264ace61SBarry Smith 
2520264ace61SBarry Smith #undef __FUNCT__
2521264ace61SBarry Smith #define __FUNCT__ "DMGetType"
2522264ace61SBarry Smith /*@C
2523264ace61SBarry Smith   DMGetType - Gets the DM type name (as a string) from the DM.
2524264ace61SBarry Smith 
2525264ace61SBarry Smith   Not Collective
2526264ace61SBarry Smith 
2527264ace61SBarry Smith   Input Parameter:
2528264ace61SBarry Smith . dm  - The DM
2529264ace61SBarry Smith 
2530264ace61SBarry Smith   Output Parameter:
2531264ace61SBarry Smith . type - The DM type name
2532264ace61SBarry Smith 
2533264ace61SBarry Smith   Level: intermediate
2534264ace61SBarry Smith 
2535264ace61SBarry Smith .keywords: DM, get, type, name
2536264ace61SBarry Smith .seealso: DMSetType(), DMCreate()
2537264ace61SBarry Smith @*/
253819fd82e9SBarry Smith PetscErrorCode  DMGetType(DM dm, DMType *type)
2539264ace61SBarry Smith {
2540264ace61SBarry Smith   PetscErrorCode ierr;
2541264ace61SBarry Smith 
2542264ace61SBarry Smith   PetscFunctionBegin;
2543264ace61SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2544264ace61SBarry Smith   PetscValidCharPointer(type,2);
2545264ace61SBarry Smith   if (!DMRegisterAllCalled) {
2546607a6623SBarry Smith     ierr = DMRegisterAll();CHKERRQ(ierr);
2547264ace61SBarry Smith   }
2548264ace61SBarry Smith   *type = ((PetscObject)dm)->type_name;
2549264ace61SBarry Smith   PetscFunctionReturn(0);
2550264ace61SBarry Smith }
2551264ace61SBarry Smith 
255267a56275SMatthew G Knepley #undef __FUNCT__
255367a56275SMatthew G Knepley #define __FUNCT__ "DMConvert"
255467a56275SMatthew G Knepley /*@C
255567a56275SMatthew G Knepley   DMConvert - Converts a DM to another DM, either of the same or different type.
255667a56275SMatthew G Knepley 
255767a56275SMatthew G Knepley   Collective on DM
255867a56275SMatthew G Knepley 
255967a56275SMatthew G Knepley   Input Parameters:
256067a56275SMatthew G Knepley + dm - the DM
256167a56275SMatthew G Knepley - newtype - new DM type (use "same" for the same type)
256267a56275SMatthew G Knepley 
256367a56275SMatthew G Knepley   Output Parameter:
256467a56275SMatthew G Knepley . M - pointer to new DM
256567a56275SMatthew G Knepley 
256667a56275SMatthew G Knepley   Notes:
256767a56275SMatthew G Knepley   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
256867a56275SMatthew G Knepley   the MPI communicator of the generated DM is always the same as the communicator
256967a56275SMatthew G Knepley   of the input DM.
257067a56275SMatthew G Knepley 
257167a56275SMatthew G Knepley   Level: intermediate
257267a56275SMatthew G Knepley 
257367a56275SMatthew G Knepley .seealso: DMCreate()
257467a56275SMatthew G Knepley @*/
257519fd82e9SBarry Smith PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M)
257667a56275SMatthew G Knepley {
257767a56275SMatthew G Knepley   DM             B;
257867a56275SMatthew G Knepley   char           convname[256];
257967a56275SMatthew G Knepley   PetscBool      sametype, issame;
258067a56275SMatthew G Knepley   PetscErrorCode ierr;
258167a56275SMatthew G Knepley 
258267a56275SMatthew G Knepley   PetscFunctionBegin;
258367a56275SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
258467a56275SMatthew G Knepley   PetscValidType(dm,1);
258567a56275SMatthew G Knepley   PetscValidPointer(M,3);
2586251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
258767a56275SMatthew G Knepley   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
258867a56275SMatthew G Knepley   {
25890298fd71SBarry Smith     PetscErrorCode (*conv)(DM, DMType, DM*) = NULL;
259067a56275SMatthew G Knepley 
259167a56275SMatthew G Knepley     /*
259267a56275SMatthew G Knepley        Order of precedence:
259367a56275SMatthew G Knepley        1) See if a specialized converter is known to the current DM.
259467a56275SMatthew G Knepley        2) See if a specialized converter is known to the desired DM class.
259567a56275SMatthew G Knepley        3) See if a good general converter is registered for the desired class
259667a56275SMatthew G Knepley        4) See if a good general converter is known for the current matrix.
259767a56275SMatthew G Knepley        5) Use a really basic converter.
259867a56275SMatthew G Knepley     */
259967a56275SMatthew G Knepley 
260067a56275SMatthew G Knepley     /* 1) See if a specialized converter is known to the current DM and the desired class */
260167a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
260267a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
260367a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
260467a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
260567a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
26060005d66cSJed Brown     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,&conv);CHKERRQ(ierr);
260767a56275SMatthew G Knepley     if (conv) goto foundconv;
260867a56275SMatthew G Knepley 
260967a56275SMatthew G Knepley     /* 2)  See if a specialized converter is known to the desired DM class. */
261082f516ccSBarry Smith     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &B);CHKERRQ(ierr);
261167a56275SMatthew G Knepley     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
261267a56275SMatthew G Knepley     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
261367a56275SMatthew G Knepley     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
261467a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
261567a56275SMatthew G Knepley     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
261667a56275SMatthew G Knepley     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
26170005d66cSJed Brown     ierr = PetscObjectQueryFunction((PetscObject)B,convname,&conv);CHKERRQ(ierr);
261867a56275SMatthew G Knepley     if (conv) {
2619fcfd50ebSBarry Smith       ierr = DMDestroy(&B);CHKERRQ(ierr);
262067a56275SMatthew G Knepley       goto foundconv;
262167a56275SMatthew G Knepley     }
262267a56275SMatthew G Knepley 
262367a56275SMatthew G Knepley #if 0
262467a56275SMatthew G Knepley     /* 3) See if a good general converter is registered for the desired class */
262567a56275SMatthew G Knepley     conv = B->ops->convertfrom;
2626fcfd50ebSBarry Smith     ierr = DMDestroy(&B);CHKERRQ(ierr);
262767a56275SMatthew G Knepley     if (conv) goto foundconv;
262867a56275SMatthew G Knepley 
262967a56275SMatthew G Knepley     /* 4) See if a good general converter is known for the current matrix */
263067a56275SMatthew G Knepley     if (dm->ops->convert) {
263167a56275SMatthew G Knepley       conv = dm->ops->convert;
263267a56275SMatthew G Knepley     }
263367a56275SMatthew G Knepley     if (conv) goto foundconv;
263467a56275SMatthew G Knepley #endif
263567a56275SMatthew G Knepley 
263667a56275SMatthew G Knepley     /* 5) Use a really basic converter. */
263782f516ccSBarry Smith     SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
263867a56275SMatthew G Knepley 
263967a56275SMatthew G Knepley foundconv:
264067a56275SMatthew G Knepley     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
264167a56275SMatthew G Knepley     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
264267a56275SMatthew G Knepley     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
264367a56275SMatthew G Knepley   }
264467a56275SMatthew G Knepley   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
264567a56275SMatthew G Knepley   PetscFunctionReturn(0);
264667a56275SMatthew G Knepley }
2647264ace61SBarry Smith 
2648264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/
2649264ace61SBarry Smith 
2650264ace61SBarry Smith #undef __FUNCT__
2651264ace61SBarry Smith #define __FUNCT__ "DMRegister"
2652264ace61SBarry Smith /*@C
26531c84c290SBarry Smith   DMRegister -  Adds a new DM component implementation
26541c84c290SBarry Smith 
26551c84c290SBarry Smith   Not Collective
26561c84c290SBarry Smith 
26571c84c290SBarry Smith   Input Parameters:
26581c84c290SBarry Smith + name        - The name of a new user-defined creation routine
26591c84c290SBarry Smith - create_func - The creation routine itself
26601c84c290SBarry Smith 
26611c84c290SBarry Smith   Notes:
26621c84c290SBarry Smith   DMRegister() may be called multiple times to add several user-defined DMs
26631c84c290SBarry Smith 
26641c84c290SBarry Smith 
26651c84c290SBarry Smith   Sample usage:
26661c84c290SBarry Smith .vb
2667bdf89e91SBarry Smith     DMRegister("my_da", MyDMCreate);
26681c84c290SBarry Smith .ve
26691c84c290SBarry Smith 
26701c84c290SBarry Smith   Then, your DM type can be chosen with the procedural interface via
26711c84c290SBarry Smith .vb
26721c84c290SBarry Smith     DMCreate(MPI_Comm, DM *);
26731c84c290SBarry Smith     DMSetType(DM,"my_da");
26741c84c290SBarry Smith .ve
26751c84c290SBarry Smith    or at runtime via the option
26761c84c290SBarry Smith .vb
26771c84c290SBarry Smith     -da_type my_da
26781c84c290SBarry Smith .ve
2679264ace61SBarry Smith 
2680264ace61SBarry Smith   Level: advanced
26811c84c290SBarry Smith 
26821c84c290SBarry Smith .keywords: DM, register
2683bdf89e91SBarry Smith .seealso: DMRegisterAll(), DMRegisterDestroy()
26841c84c290SBarry Smith 
2685264ace61SBarry Smith @*/
2686bdf89e91SBarry Smith PetscErrorCode  DMRegister(const char sname[],PetscErrorCode (*function)(DM))
2687264ace61SBarry Smith {
2688264ace61SBarry Smith   PetscErrorCode ierr;
2689264ace61SBarry Smith 
2690264ace61SBarry Smith   PetscFunctionBegin;
2691a240a19fSJed Brown   ierr = PetscFunctionListAdd(&DMList,sname,function);CHKERRQ(ierr);
2692264ace61SBarry Smith   PetscFunctionReturn(0);
2693264ace61SBarry Smith }
2694264ace61SBarry Smith 
2695b859378eSBarry Smith #undef __FUNCT__
2696b859378eSBarry Smith #define __FUNCT__ "DMLoad"
2697b859378eSBarry Smith /*@C
269855849f57SBarry Smith   DMLoad - Loads a DM that has been stored in binary  with DMView().
2699b859378eSBarry Smith 
2700b859378eSBarry Smith   Collective on PetscViewer
2701b859378eSBarry Smith 
2702b859378eSBarry Smith   Input Parameters:
2703b859378eSBarry Smith + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2704b859378eSBarry Smith            some related function before a call to DMLoad().
2705b859378eSBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2706b859378eSBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2707b859378eSBarry Smith 
2708b859378eSBarry Smith    Level: intermediate
2709b859378eSBarry Smith 
2710b859378eSBarry Smith   Notes:
271155849f57SBarry Smith    The type is determined by the data in the file, any type set into the DM before this call is ignored.
2712b859378eSBarry Smith 
2713b859378eSBarry Smith   Notes for advanced users:
2714b859378eSBarry Smith   Most users should not need to know the details of the binary storage
2715b859378eSBarry Smith   format, since DMLoad() and DMView() completely hide these details.
2716b859378eSBarry Smith   But for anyone who's interested, the standard binary matrix storage
2717b859378eSBarry Smith   format is
2718b859378eSBarry Smith .vb
2719b859378eSBarry Smith      has not yet been determined
2720b859378eSBarry Smith .ve
2721b859378eSBarry Smith 
2722b859378eSBarry Smith .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2723b859378eSBarry Smith @*/
2724b859378eSBarry Smith PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2725b859378eSBarry Smith {
27269331c7a4SMatthew G. Knepley   PetscBool      isbinary, ishdf5;
2727b859378eSBarry Smith   PetscErrorCode ierr;
2728b859378eSBarry Smith 
2729b859378eSBarry Smith   PetscFunctionBegin;
2730b859378eSBarry Smith   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2731b859378eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
273232c0f0efSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
27339331c7a4SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
27349331c7a4SMatthew G. Knepley   if (isbinary) {
27359331c7a4SMatthew G. Knepley     PetscInt classid;
27369331c7a4SMatthew G. Knepley     char     type[256];
2737b859378eSBarry Smith 
273832c0f0efSBarry Smith     ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
27399200755eSBarry Smith     if (classid != DM_FILE_CLASSID) SETERRQ1(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not DM next in file, classid found %d",(int)classid);
274032c0f0efSBarry Smith     ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
274132c0f0efSBarry Smith     ierr = DMSetType(newdm, type);CHKERRQ(ierr);
27429331c7a4SMatthew G. Knepley     if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);}
27439331c7a4SMatthew G. Knepley   } else if (ishdf5) {
27449331c7a4SMatthew G. Knepley     if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);}
27459331c7a4SMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen() or PetscViewerHDF5Open()");
2746b859378eSBarry Smith   PetscFunctionReturn(0);
2747b859378eSBarry Smith }
2748b859378eSBarry Smith 
27497da65231SMatthew G Knepley /******************************** FEM Support **********************************/
27507da65231SMatthew G Knepley 
27517da65231SMatthew G Knepley #undef __FUNCT__
27527da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellVector"
2753a6dfd86eSKarl Rupp PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[])
2754a6dfd86eSKarl Rupp {
27551d47ebbbSSatish Balay   PetscInt       f;
27561b30c384SMatthew G Knepley   PetscErrorCode ierr;
27571b30c384SMatthew G Knepley 
27587da65231SMatthew G Knepley   PetscFunctionBegin;
275974778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
27601d47ebbbSSatish Balay   for (f = 0; f < len; ++f) {
276157622a8eSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", (double)PetscRealPart(x[f]));CHKERRQ(ierr);
27627da65231SMatthew G Knepley   }
27637da65231SMatthew G Knepley   PetscFunctionReturn(0);
27647da65231SMatthew G Knepley }
27657da65231SMatthew G Knepley 
27667da65231SMatthew G Knepley #undef __FUNCT__
27677da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellMatrix"
2768a6dfd86eSKarl Rupp PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
2769a6dfd86eSKarl Rupp {
27701b30c384SMatthew G Knepley   PetscInt       f, g;
27717da65231SMatthew G Knepley   PetscErrorCode ierr;
27727da65231SMatthew G Knepley 
27737da65231SMatthew G Knepley   PetscFunctionBegin;
277474778d6cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
27751d47ebbbSSatish Balay   for (f = 0; f < rows; ++f) {
277674778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
27771d47ebbbSSatish Balay     for (g = 0; g < cols; ++g) {
2778e3556bceSMatthew G. Knepley       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5g", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
27797da65231SMatthew G Knepley     }
278074778d6cSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
27817da65231SMatthew G Knepley   }
27827da65231SMatthew G Knepley   PetscFunctionReturn(0);
27837da65231SMatthew G Knepley }
2784e7c4fc90SDmitry Karpeev 
2785970e74d5SMatthew G Knepley #undef __FUNCT__
2786e759306cSMatthew G. Knepley #define __FUNCT__ "DMPrintLocalVec"
27876113b454SMatthew G. Knepley PetscErrorCode DMPrintLocalVec(DM dm, const char name[], PetscReal tol, Vec X)
2788e759306cSMatthew G. Knepley {
2789e759306cSMatthew G. Knepley   PetscMPIInt    rank, numProcs;
2790e759306cSMatthew G. Knepley   PetscInt       p;
2791e759306cSMatthew G. Knepley   PetscErrorCode ierr;
2792e759306cSMatthew G. Knepley 
2793e759306cSMatthew G. Knepley   PetscFunctionBegin;
2794e759306cSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
2795e759306cSMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
2796e759306cSMatthew G. Knepley   ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s:\n", name);CHKERRQ(ierr);
2797e759306cSMatthew G. Knepley   for (p = 0; p < numProcs; ++p) {
2798e759306cSMatthew G. Knepley     if (p == rank) {
2799e759306cSMatthew G. Knepley       Vec x;
2800e759306cSMatthew G. Knepley 
2801e759306cSMatthew G. Knepley       ierr = VecDuplicate(X, &x);CHKERRQ(ierr);
2802e759306cSMatthew G. Knepley       ierr = VecCopy(X, x);CHKERRQ(ierr);
28036113b454SMatthew G. Knepley       ierr = VecChop(x, tol);CHKERRQ(ierr);
2804e759306cSMatthew G. Knepley       ierr = VecView(x, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2805e759306cSMatthew G. Knepley       ierr = VecDestroy(&x);CHKERRQ(ierr);
2806e759306cSMatthew G. Knepley       ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2807e759306cSMatthew G. Knepley     }
2808e759306cSMatthew G. Knepley     ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
2809e759306cSMatthew G. Knepley   }
2810e759306cSMatthew G. Knepley   PetscFunctionReturn(0);
2811e759306cSMatthew G. Knepley }
2812e759306cSMatthew G. Knepley 
2813e759306cSMatthew G. Knepley #undef __FUNCT__
281488ed4aceSMatthew G Knepley #define __FUNCT__ "DMGetDefaultSection"
281588ed4aceSMatthew G Knepley /*@
281688ed4aceSMatthew G Knepley   DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM.
281788ed4aceSMatthew G Knepley 
281888ed4aceSMatthew G Knepley   Input Parameter:
281988ed4aceSMatthew G Knepley . dm - The DM
282088ed4aceSMatthew G Knepley 
282188ed4aceSMatthew G Knepley   Output Parameter:
282288ed4aceSMatthew G Knepley . section - The PetscSection
282388ed4aceSMatthew G Knepley 
282488ed4aceSMatthew G Knepley   Level: intermediate
282588ed4aceSMatthew G Knepley 
282688ed4aceSMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
282788ed4aceSMatthew G Knepley 
282888ed4aceSMatthew G Knepley .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
282988ed4aceSMatthew G Knepley @*/
28300adebc6cSBarry Smith PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section)
28310adebc6cSBarry Smith {
2832fd59a867SMatthew G. Knepley   PetscErrorCode ierr;
2833fd59a867SMatthew G. Knepley 
283488ed4aceSMatthew G Knepley   PetscFunctionBegin;
283588ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
283688ed4aceSMatthew G Knepley   PetscValidPointer(section, 2);
2837fd59a867SMatthew G. Knepley   if (!dm->defaultSection && dm->ops->createdefaultsection) {ierr = (*dm->ops->createdefaultsection)(dm);CHKERRQ(ierr);}
283888ed4aceSMatthew G Knepley   *section = dm->defaultSection;
283988ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
284088ed4aceSMatthew G Knepley }
284188ed4aceSMatthew G Knepley 
284288ed4aceSMatthew G Knepley #undef __FUNCT__
284388ed4aceSMatthew G Knepley #define __FUNCT__ "DMSetDefaultSection"
284488ed4aceSMatthew G Knepley /*@
284588ed4aceSMatthew G Knepley   DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM.
284688ed4aceSMatthew G Knepley 
284788ed4aceSMatthew G Knepley   Input Parameters:
284888ed4aceSMatthew G Knepley + dm - The DM
284988ed4aceSMatthew G Knepley - section - The PetscSection
285088ed4aceSMatthew G Knepley 
285188ed4aceSMatthew G Knepley   Level: intermediate
285288ed4aceSMatthew G Knepley 
285388ed4aceSMatthew G Knepley   Note: Any existing Section will be destroyed
285488ed4aceSMatthew G Knepley 
285588ed4aceSMatthew G Knepley .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
285688ed4aceSMatthew G Knepley @*/
28570adebc6cSBarry Smith PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section)
28580adebc6cSBarry Smith {
2859af122d2aSMatthew G Knepley   PetscInt       numFields;
2860af122d2aSMatthew G Knepley   PetscInt       f;
286188ed4aceSMatthew G Knepley   PetscErrorCode ierr;
286288ed4aceSMatthew G Knepley 
286388ed4aceSMatthew G Knepley   PetscFunctionBegin;
286488ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28651d799100SJed Brown   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
28661d799100SJed Brown   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
286788ed4aceSMatthew G Knepley   ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr);
286888ed4aceSMatthew G Knepley   dm->defaultSection = section;
2869af122d2aSMatthew G Knepley   ierr = PetscSectionGetNumFields(dm->defaultSection, &numFields);CHKERRQ(ierr);
2870af122d2aSMatthew G Knepley   if (numFields) {
2871af122d2aSMatthew G Knepley     ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr);
2872af122d2aSMatthew G Knepley     for (f = 0; f < numFields; ++f) {
2873af122d2aSMatthew G Knepley       const char *name;
2874af122d2aSMatthew G Knepley 
2875af122d2aSMatthew G Knepley       ierr = PetscSectionGetFieldName(dm->defaultSection, f, &name);CHKERRQ(ierr);
2876decb47aaSMatthew G. Knepley       ierr = PetscObjectSetName((PetscObject) dm->fields[f], name);CHKERRQ(ierr);
2877af122d2aSMatthew G Knepley     }
2878af122d2aSMatthew G Knepley   }
28791d799100SJed Brown   /* The global section will be rebuilt in the next call to DMGetDefaultGlobalSection(). */
28801d799100SJed Brown   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
288188ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
288288ed4aceSMatthew G Knepley }
288388ed4aceSMatthew G Knepley 
288488ed4aceSMatthew G Knepley #undef __FUNCT__
288588ed4aceSMatthew G Knepley #define __FUNCT__ "DMGetDefaultGlobalSection"
288688ed4aceSMatthew G Knepley /*@
288788ed4aceSMatthew G Knepley   DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM.
288888ed4aceSMatthew G Knepley 
28898b1ab98fSJed Brown   Collective on DM
28908b1ab98fSJed Brown 
289188ed4aceSMatthew G Knepley   Input Parameter:
289288ed4aceSMatthew G Knepley . dm - The DM
289388ed4aceSMatthew G Knepley 
289488ed4aceSMatthew G Knepley   Output Parameter:
289588ed4aceSMatthew G Knepley . section - The PetscSection
289688ed4aceSMatthew G Knepley 
289788ed4aceSMatthew G Knepley   Level: intermediate
289888ed4aceSMatthew G Knepley 
289988ed4aceSMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
290088ed4aceSMatthew G Knepley 
290188ed4aceSMatthew G Knepley .seealso: DMSetDefaultSection(), DMGetDefaultSection()
290288ed4aceSMatthew G Knepley @*/
29030adebc6cSBarry Smith PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section)
29040adebc6cSBarry Smith {
290588ed4aceSMatthew G Knepley   PetscErrorCode ierr;
290688ed4aceSMatthew G Knepley 
290788ed4aceSMatthew G Knepley   PetscFunctionBegin;
290888ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
290988ed4aceSMatthew G Knepley   PetscValidPointer(section, 2);
291088ed4aceSMatthew G Knepley   if (!dm->defaultGlobalSection) {
2911fd59a867SMatthew G. Knepley     PetscSection s;
2912fd59a867SMatthew G. Knepley 
2913fd59a867SMatthew G. Knepley     ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
2914fd59a867SMatthew 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");
2915fd59a867SMatthew 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");
2916fd59a867SMatthew G. Knepley     ierr = PetscSectionCreateGlobalSection(s, dm->sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr);
2917cf06b437SMatthew G. Knepley     ierr = PetscLayoutDestroy(&dm->map);CHKERRQ(ierr);
2918ce94432eSBarry Smith     ierr = PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm), dm->defaultGlobalSection, &dm->map);CHKERRQ(ierr);
291988ed4aceSMatthew G Knepley   }
292088ed4aceSMatthew G Knepley   *section = dm->defaultGlobalSection;
292188ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
292288ed4aceSMatthew G Knepley }
292388ed4aceSMatthew G Knepley 
292488ed4aceSMatthew G Knepley #undef __FUNCT__
2925b21d0597SMatthew G Knepley #define __FUNCT__ "DMSetDefaultGlobalSection"
2926b21d0597SMatthew G Knepley /*@
2927b21d0597SMatthew G Knepley   DMSetDefaultGlobalSection - Set the PetscSection encoding the global data layout for the DM.
2928b21d0597SMatthew G Knepley 
2929b21d0597SMatthew G Knepley   Input Parameters:
2930b21d0597SMatthew G Knepley + dm - The DM
29315080bbdbSMatthew G Knepley - section - The PetscSection, or NULL
2932b21d0597SMatthew G Knepley 
2933b21d0597SMatthew G Knepley   Level: intermediate
2934b21d0597SMatthew G Knepley 
2935b21d0597SMatthew G Knepley   Note: Any existing Section will be destroyed
2936b21d0597SMatthew G Knepley 
2937b21d0597SMatthew G Knepley .seealso: DMGetDefaultGlobalSection(), DMSetDefaultSection()
2938b21d0597SMatthew G Knepley @*/
29390adebc6cSBarry Smith PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection section)
29400adebc6cSBarry Smith {
2941b21d0597SMatthew G Knepley   PetscErrorCode ierr;
2942b21d0597SMatthew G Knepley 
2943b21d0597SMatthew G Knepley   PetscFunctionBegin;
2944b21d0597SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
29455080bbdbSMatthew G Knepley   if (section) PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
29461d799100SJed Brown   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
2947b21d0597SMatthew G Knepley   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
2948b21d0597SMatthew G Knepley   dm->defaultGlobalSection = section;
2949b21d0597SMatthew G Knepley   PetscFunctionReturn(0);
2950b21d0597SMatthew G Knepley }
2951b21d0597SMatthew G Knepley 
2952b21d0597SMatthew G Knepley #undef __FUNCT__
295388ed4aceSMatthew G Knepley #define __FUNCT__ "DMGetDefaultSF"
295488ed4aceSMatthew G Knepley /*@
295588ed4aceSMatthew G Knepley   DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set,
295688ed4aceSMatthew G Knepley   it is created from the default PetscSection layouts in the DM.
295788ed4aceSMatthew G Knepley 
295888ed4aceSMatthew G Knepley   Input Parameter:
295988ed4aceSMatthew G Knepley . dm - The DM
296088ed4aceSMatthew G Knepley 
296188ed4aceSMatthew G Knepley   Output Parameter:
296288ed4aceSMatthew G Knepley . sf - The PetscSF
296388ed4aceSMatthew G Knepley 
296488ed4aceSMatthew G Knepley   Level: intermediate
296588ed4aceSMatthew G Knepley 
296688ed4aceSMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
296788ed4aceSMatthew G Knepley 
296888ed4aceSMatthew G Knepley .seealso: DMSetDefaultSF(), DMCreateDefaultSF()
296988ed4aceSMatthew G Knepley @*/
29700adebc6cSBarry Smith PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf)
29710adebc6cSBarry Smith {
297288ed4aceSMatthew G Knepley   PetscInt       nroots;
297388ed4aceSMatthew G Knepley   PetscErrorCode ierr;
297488ed4aceSMatthew G Knepley 
297588ed4aceSMatthew G Knepley   PetscFunctionBegin;
297688ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
297788ed4aceSMatthew G Knepley   PetscValidPointer(sf, 2);
29780298fd71SBarry Smith   ierr = PetscSFGetGraph(dm->defaultSF, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
297988ed4aceSMatthew G Knepley   if (nroots < 0) {
298088ed4aceSMatthew G Knepley     PetscSection section, gSection;
298188ed4aceSMatthew G Knepley 
298288ed4aceSMatthew G Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
298331ea6d37SMatthew G Knepley     if (section) {
298488ed4aceSMatthew G Knepley       ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
298588ed4aceSMatthew G Knepley       ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr);
298631ea6d37SMatthew G Knepley     } else {
29870298fd71SBarry Smith       *sf = NULL;
298831ea6d37SMatthew G Knepley       PetscFunctionReturn(0);
298931ea6d37SMatthew G Knepley     }
299088ed4aceSMatthew G Knepley   }
299188ed4aceSMatthew G Knepley   *sf = dm->defaultSF;
299288ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
299388ed4aceSMatthew G Knepley }
299488ed4aceSMatthew G Knepley 
299588ed4aceSMatthew G Knepley #undef __FUNCT__
299688ed4aceSMatthew G Knepley #define __FUNCT__ "DMSetDefaultSF"
299788ed4aceSMatthew G Knepley /*@
299888ed4aceSMatthew G Knepley   DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM
299988ed4aceSMatthew G Knepley 
300088ed4aceSMatthew G Knepley   Input Parameters:
300188ed4aceSMatthew G Knepley + dm - The DM
300288ed4aceSMatthew G Knepley - sf - The PetscSF
300388ed4aceSMatthew G Knepley 
300488ed4aceSMatthew G Knepley   Level: intermediate
300588ed4aceSMatthew G Knepley 
300688ed4aceSMatthew G Knepley   Note: Any previous SF is destroyed
300788ed4aceSMatthew G Knepley 
300888ed4aceSMatthew G Knepley .seealso: DMGetDefaultSF(), DMCreateDefaultSF()
300988ed4aceSMatthew G Knepley @*/
30100adebc6cSBarry Smith PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf)
30110adebc6cSBarry Smith {
301288ed4aceSMatthew G Knepley   PetscErrorCode ierr;
301388ed4aceSMatthew G Knepley 
301488ed4aceSMatthew G Knepley   PetscFunctionBegin;
301588ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
301688ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
301788ed4aceSMatthew G Knepley   ierr          = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr);
301888ed4aceSMatthew G Knepley   dm->defaultSF = sf;
301988ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
302088ed4aceSMatthew G Knepley }
302188ed4aceSMatthew G Knepley 
302288ed4aceSMatthew G Knepley #undef __FUNCT__
302388ed4aceSMatthew G Knepley #define __FUNCT__ "DMCreateDefaultSF"
302488ed4aceSMatthew G Knepley /*@C
302588ed4aceSMatthew G Knepley   DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections
302688ed4aceSMatthew G Knepley   describing the data layout.
302788ed4aceSMatthew G Knepley 
302888ed4aceSMatthew G Knepley   Input Parameters:
302988ed4aceSMatthew G Knepley + dm - The DM
303088ed4aceSMatthew G Knepley . localSection - PetscSection describing the local data layout
303188ed4aceSMatthew G Knepley - globalSection - PetscSection describing the global data layout
303288ed4aceSMatthew G Knepley 
303388ed4aceSMatthew G Knepley   Level: intermediate
303488ed4aceSMatthew G Knepley 
303588ed4aceSMatthew G Knepley .seealso: DMGetDefaultSF(), DMSetDefaultSF()
303688ed4aceSMatthew G Knepley @*/
303788ed4aceSMatthew G Knepley PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection)
303888ed4aceSMatthew G Knepley {
303982f516ccSBarry Smith   MPI_Comm       comm;
304088ed4aceSMatthew G Knepley   PetscLayout    layout;
304188ed4aceSMatthew G Knepley   const PetscInt *ranges;
304288ed4aceSMatthew G Knepley   PetscInt       *local;
304388ed4aceSMatthew G Knepley   PetscSFNode    *remote;
3044ecd73843SMatthew G. Knepley   PetscInt       pStart, pEnd, p, nroots, nleaves = 0, l;
304588ed4aceSMatthew G Knepley   PetscMPIInt    size, rank;
304688ed4aceSMatthew G Knepley   PetscErrorCode ierr;
304788ed4aceSMatthew G Knepley 
304888ed4aceSMatthew G Knepley   PetscFunctionBegin;
304982f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
305088ed4aceSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
305188ed4aceSMatthew G Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
305288ed4aceSMatthew G Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
305388ed4aceSMatthew G Knepley   ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr);
305488ed4aceSMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr);
305588ed4aceSMatthew G Knepley   ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
305688ed4aceSMatthew G Knepley   ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
305788ed4aceSMatthew G Knepley   ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr);
305888ed4aceSMatthew G Knepley   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
305988ed4aceSMatthew G Knepley   ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
3060ecd73843SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
30616636e97aSMatthew G Knepley     PetscInt gdof, gcdof;
306288ed4aceSMatthew G Knepley 
30636636e97aSMatthew G Knepley     ierr     = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
30646636e97aSMatthew G Knepley     ierr     = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
30656636e97aSMatthew G Knepley     nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
306688ed4aceSMatthew G Knepley   }
3067785e854fSJed Brown   ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr);
3068785e854fSJed Brown   ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr);
306988ed4aceSMatthew G Knepley   for (p = pStart, l = 0; p < pEnd; ++p) {
30701f588964SMatthew G Knepley     const PetscInt *cind;
30716636e97aSMatthew G Knepley     PetscInt       dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
307288ed4aceSMatthew G Knepley 
307388ed4aceSMatthew G Knepley     ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr);
307488ed4aceSMatthew G Knepley     ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr);
307588ed4aceSMatthew G Knepley     ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr);
307688ed4aceSMatthew G Knepley     ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr);
307788ed4aceSMatthew G Knepley     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
30786636e97aSMatthew G Knepley     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
307988ed4aceSMatthew G Knepley     ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
30806636e97aSMatthew G Knepley     if (!gdof) continue; /* Censored point */
30816636e97aSMatthew G Knepley     gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
30826636e97aSMatthew G Knepley     if (gsize != dof-cdof) {
3083057b4bcdSMatthew 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);
30846636e97aSMatthew G Knepley       cdof = 0; /* Ignore constraints */
30856636e97aSMatthew G Knepley     }
308688ed4aceSMatthew G Knepley     for (d = 0, c = 0; d < dof; ++d) {
308788ed4aceSMatthew G Knepley       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
308888ed4aceSMatthew G Knepley       local[l+d-c] = off+d;
308988ed4aceSMatthew G Knepley     }
309088ed4aceSMatthew G Knepley     if (gdof < 0) {
30916636e97aSMatthew G Knepley       for (d = 0; d < gsize; ++d, ++l) {
309288ed4aceSMatthew G Knepley         PetscInt offset = -(goff+1) + d, r;
309388ed4aceSMatthew G Knepley 
309405376888SMatthew G. Knepley         ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr);
309531d3f06eSJed Brown         if (r < 0) r = -(r+2);
309605376888SMatthew 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);
309788ed4aceSMatthew G Knepley         remote[l].rank  = r;
309888ed4aceSMatthew G Knepley         remote[l].index = offset - ranges[r];
309988ed4aceSMatthew G Knepley       }
310088ed4aceSMatthew G Knepley     } else {
31016636e97aSMatthew G Knepley       for (d = 0; d < gsize; ++d, ++l) {
310288ed4aceSMatthew G Knepley         remote[l].rank  = rank;
310388ed4aceSMatthew G Knepley         remote[l].index = goff+d - ranges[rank];
310488ed4aceSMatthew G Knepley       }
310588ed4aceSMatthew G Knepley     }
310688ed4aceSMatthew G Knepley   }
31076636e97aSMatthew G Knepley   if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %d != nleaves %d", l, nleaves);
310888ed4aceSMatthew G Knepley   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
310988ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
311088ed4aceSMatthew G Knepley   PetscFunctionReturn(0);
311188ed4aceSMatthew G Knepley }
3112af122d2aSMatthew G Knepley 
3113af122d2aSMatthew G Knepley #undef __FUNCT__
3114b21d0597SMatthew G Knepley #define __FUNCT__ "DMGetPointSF"
3115b21d0597SMatthew G Knepley /*@
3116b21d0597SMatthew G Knepley   DMGetPointSF - Get the PetscSF encoding the parallel section point overlap for the DM.
3117b21d0597SMatthew G Knepley 
3118b21d0597SMatthew G Knepley   Input Parameter:
3119b21d0597SMatthew G Knepley . dm - The DM
3120b21d0597SMatthew G Knepley 
3121b21d0597SMatthew G Knepley   Output Parameter:
3122b21d0597SMatthew G Knepley . sf - The PetscSF
3123b21d0597SMatthew G Knepley 
3124b21d0597SMatthew G Knepley   Level: intermediate
3125b21d0597SMatthew G Knepley 
3126b21d0597SMatthew G Knepley   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
3127b21d0597SMatthew G Knepley 
3128057b4bcdSMatthew G Knepley .seealso: DMSetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3129b21d0597SMatthew G Knepley @*/
31300adebc6cSBarry Smith PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf)
31310adebc6cSBarry Smith {
3132b21d0597SMatthew G Knepley   PetscFunctionBegin;
3133b21d0597SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3134b21d0597SMatthew G Knepley   PetscValidPointer(sf, 2);
3135b21d0597SMatthew G Knepley   *sf = dm->sf;
3136b21d0597SMatthew G Knepley   PetscFunctionReturn(0);
3137b21d0597SMatthew G Knepley }
3138b21d0597SMatthew G Knepley 
3139b21d0597SMatthew G Knepley #undef __FUNCT__
3140057b4bcdSMatthew G Knepley #define __FUNCT__ "DMSetPointSF"
3141057b4bcdSMatthew G Knepley /*@
3142057b4bcdSMatthew G Knepley   DMSetPointSF - Set the PetscSF encoding the parallel section point overlap for the DM.
3143057b4bcdSMatthew G Knepley 
3144057b4bcdSMatthew G Knepley   Input Parameters:
3145057b4bcdSMatthew G Knepley + dm - The DM
3146057b4bcdSMatthew G Knepley - sf - The PetscSF
3147057b4bcdSMatthew G Knepley 
3148057b4bcdSMatthew G Knepley   Level: intermediate
3149057b4bcdSMatthew G Knepley 
3150057b4bcdSMatthew G Knepley .seealso: DMGetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3151057b4bcdSMatthew G Knepley @*/
31520adebc6cSBarry Smith PetscErrorCode DMSetPointSF(DM dm, PetscSF sf)
31530adebc6cSBarry Smith {
3154057b4bcdSMatthew G Knepley   PetscErrorCode ierr;
3155057b4bcdSMatthew G Knepley 
3156057b4bcdSMatthew G Knepley   PetscFunctionBegin;
3157057b4bcdSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3158057b4bcdSMatthew G Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
3159057b4bcdSMatthew G Knepley   ierr   = PetscSFDestroy(&dm->sf);CHKERRQ(ierr);
3160057b4bcdSMatthew G Knepley   ierr   = PetscObjectReference((PetscObject) sf);CHKERRQ(ierr);
3161057b4bcdSMatthew G Knepley   dm->sf = sf;
3162057b4bcdSMatthew G Knepley   PetscFunctionReturn(0);
3163057b4bcdSMatthew G Knepley }
3164057b4bcdSMatthew G Knepley 
3165057b4bcdSMatthew G Knepley #undef __FUNCT__
3166af122d2aSMatthew G Knepley #define __FUNCT__ "DMGetNumFields"
3167af122d2aSMatthew G Knepley PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields)
3168af122d2aSMatthew G Knepley {
3169af122d2aSMatthew G Knepley   PetscFunctionBegin;
3170af122d2aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3171af122d2aSMatthew G Knepley   PetscValidPointer(numFields, 2);
3172af122d2aSMatthew G Knepley   *numFields = dm->numFields;
3173af122d2aSMatthew G Knepley   PetscFunctionReturn(0);
3174af122d2aSMatthew G Knepley }
3175af122d2aSMatthew G Knepley 
3176af122d2aSMatthew G Knepley #undef __FUNCT__
3177af122d2aSMatthew G Knepley #define __FUNCT__ "DMSetNumFields"
3178af122d2aSMatthew G Knepley PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields)
3179af122d2aSMatthew G Knepley {
3180af122d2aSMatthew G Knepley   PetscInt       f;
3181af122d2aSMatthew G Knepley   PetscErrorCode ierr;
3182af122d2aSMatthew G Knepley 
3183af122d2aSMatthew G Knepley   PetscFunctionBegin;
3184af122d2aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3185decb47aaSMatthew G. Knepley   if (dm->numFields == numFields) PetscFunctionReturn(0);
3186af122d2aSMatthew G Knepley   for (f = 0; f < dm->numFields; ++f) {
3187decb47aaSMatthew G. Knepley     ierr = PetscObjectDestroy((PetscObject *) &dm->fields[f]);CHKERRQ(ierr);
3188af122d2aSMatthew G Knepley   }
3189af122d2aSMatthew G Knepley   ierr          = PetscFree(dm->fields);CHKERRQ(ierr);
3190af122d2aSMatthew G Knepley   dm->numFields = numFields;
3191785e854fSJed Brown   ierr          = PetscMalloc1(dm->numFields, &dm->fields);CHKERRQ(ierr);
3192af122d2aSMatthew G Knepley   for (f = 0; f < dm->numFields; ++f) {
319382f516ccSBarry Smith     ierr = PetscContainerCreate(PetscObjectComm((PetscObject) dm), (PetscContainer *) &dm->fields[f]);CHKERRQ(ierr);
3194af122d2aSMatthew G Knepley   }
3195af122d2aSMatthew G Knepley   PetscFunctionReturn(0);
3196af122d2aSMatthew G Knepley }
3197af122d2aSMatthew G Knepley 
3198af122d2aSMatthew G Knepley #undef __FUNCT__
3199af122d2aSMatthew G Knepley #define __FUNCT__ "DMGetField"
3200af122d2aSMatthew G Knepley PetscErrorCode DMGetField(DM dm, PetscInt f, PetscObject *field)
3201af122d2aSMatthew G Knepley {
3202af122d2aSMatthew G Knepley   PetscFunctionBegin;
3203af122d2aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3204decb47aaSMatthew G. Knepley   PetscValidPointer(field, 3);
320582f516ccSBarry Smith   if (!dm->fields) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Fields have not been setup in this DM. Call DMSetNumFields()");
320682f516ccSBarry Smith   if ((f < 0) || (f >= dm->numFields)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Field %d should be in [%d,%d)", f, 0, dm->numFields);
3207decb47aaSMatthew G. Knepley   *field = (PetscObject) dm->fields[f];
3208decb47aaSMatthew G. Knepley   PetscFunctionReturn(0);
3209decb47aaSMatthew G. Knepley }
3210decb47aaSMatthew G. Knepley 
3211decb47aaSMatthew G. Knepley #undef __FUNCT__
3212decb47aaSMatthew G. Knepley #define __FUNCT__ "DMSetField"
3213decb47aaSMatthew G. Knepley PetscErrorCode DMSetField(DM dm, PetscInt f, PetscObject field)
3214decb47aaSMatthew G. Knepley {
3215decb47aaSMatthew G. Knepley   PetscErrorCode ierr;
3216decb47aaSMatthew G. Knepley 
3217decb47aaSMatthew G. Knepley   PetscFunctionBegin;
3218decb47aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3219decb47aaSMatthew G. Knepley   if (!dm->fields) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "Fields have not been setup in this DM. Call DMSetNumFields()");
3220decb47aaSMatthew G. Knepley   if ((f < 0) || (f >= dm->numFields)) SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Field %d should be in [%d,%d)", f, 0, dm->numFields);
32219e8c3d3aSMatthew G. Knepley   if (((PetscObject) dm->fields[f]) == field) PetscFunctionReturn(0);
3222decb47aaSMatthew G. Knepley   ierr = PetscObjectDestroy((PetscObject *) &dm->fields[f]);CHKERRQ(ierr);
3223decb47aaSMatthew G. Knepley   dm->fields[f] = (PetscFE) field;
3224decb47aaSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) dm->fields[f]);CHKERRQ(ierr);
3225af122d2aSMatthew G Knepley   PetscFunctionReturn(0);
3226af122d2aSMatthew G Knepley }
32276636e97aSMatthew G Knepley 
32286636e97aSMatthew G Knepley #undef __FUNCT__
3229b64e0483SPeter Brune #define __FUNCT__ "DMRestrictHook_Coordinates"
3230b64e0483SPeter Brune PetscErrorCode DMRestrictHook_Coordinates(DM dm,DM dmc,void *ctx)
3231b64e0483SPeter Brune {
3232b64e0483SPeter Brune   DM dm_coord,dmc_coord;
3233b64e0483SPeter Brune   PetscErrorCode ierr;
3234b64e0483SPeter Brune   Vec coords,ccoords;
3235b64e0483SPeter Brune   VecScatter scat;
3236b64e0483SPeter Brune   PetscFunctionBegin;
3237b64e0483SPeter Brune   ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr);
3238b64e0483SPeter Brune   ierr = DMGetCoordinateDM(dmc,&dmc_coord);CHKERRQ(ierr);
3239b64e0483SPeter Brune   ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr);
3240b64e0483SPeter Brune   ierr = DMGetCoordinates(dmc,&ccoords);CHKERRQ(ierr);
3241b64e0483SPeter Brune   if (coords && !ccoords) {
3242b64e0483SPeter Brune     ierr = DMCreateGlobalVector(dmc_coord,&ccoords);CHKERRQ(ierr);
3243b64e0483SPeter Brune     ierr = DMCreateInjection(dmc_coord,dm_coord,&scat);CHKERRQ(ierr);
3244b64e0483SPeter Brune     ierr = VecScatterBegin(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3245b64e0483SPeter Brune     ierr = VecScatterEnd(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3246b64e0483SPeter Brune     ierr = DMSetCoordinates(dmc,ccoords);CHKERRQ(ierr);
3247b64e0483SPeter Brune     ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
3248b64e0483SPeter Brune     ierr = VecDestroy(&ccoords);CHKERRQ(ierr);
3249b64e0483SPeter Brune   }
3250b64e0483SPeter Brune   PetscFunctionReturn(0);
3251b64e0483SPeter Brune }
3252b64e0483SPeter Brune 
3253b64e0483SPeter Brune #undef __FUNCT__
325403dadc2fSPeter Brune #define __FUNCT__ "DMSubDomainHook_Coordinates"
325503dadc2fSPeter Brune static PetscErrorCode DMSubDomainHook_Coordinates(DM dm,DM subdm,void *ctx)
325603dadc2fSPeter Brune {
325703dadc2fSPeter Brune   DM dm_coord,subdm_coord;
325803dadc2fSPeter Brune   PetscErrorCode ierr;
325903dadc2fSPeter Brune   Vec coords,ccoords,clcoords;
326003dadc2fSPeter Brune   VecScatter *scat_i,*scat_g;
326103dadc2fSPeter Brune   PetscFunctionBegin;
326203dadc2fSPeter Brune   ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr);
326303dadc2fSPeter Brune   ierr = DMGetCoordinateDM(subdm,&subdm_coord);CHKERRQ(ierr);
326403dadc2fSPeter Brune   ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr);
326503dadc2fSPeter Brune   ierr = DMGetCoordinates(subdm,&ccoords);CHKERRQ(ierr);
326603dadc2fSPeter Brune   if (coords && !ccoords) {
326703dadc2fSPeter Brune     ierr = DMCreateGlobalVector(subdm_coord,&ccoords);CHKERRQ(ierr);
326803dadc2fSPeter Brune     ierr = DMCreateLocalVector(subdm_coord,&clcoords);CHKERRQ(ierr);
326903dadc2fSPeter Brune     ierr = DMCreateDomainDecompositionScatters(dm_coord,1,&subdm_coord,NULL,&scat_i,&scat_g);CHKERRQ(ierr);
327003dadc2fSPeter Brune     ierr = VecScatterBegin(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
327103dadc2fSPeter Brune     ierr = VecScatterBegin(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
327203dadc2fSPeter Brune     ierr = VecScatterEnd(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
327303dadc2fSPeter Brune     ierr = VecScatterEnd(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
327403dadc2fSPeter Brune     ierr = DMSetCoordinates(subdm,ccoords);CHKERRQ(ierr);
327503dadc2fSPeter Brune     ierr = DMSetCoordinatesLocal(subdm,clcoords);CHKERRQ(ierr);
327603dadc2fSPeter Brune     ierr = VecScatterDestroy(&scat_i[0]);CHKERRQ(ierr);
327703dadc2fSPeter Brune     ierr = VecScatterDestroy(&scat_g[0]);CHKERRQ(ierr);
327803dadc2fSPeter Brune     ierr = VecDestroy(&ccoords);CHKERRQ(ierr);
327903dadc2fSPeter Brune     ierr = VecDestroy(&clcoords);CHKERRQ(ierr);
328003dadc2fSPeter Brune     ierr = PetscFree(scat_i);CHKERRQ(ierr);
328103dadc2fSPeter Brune     ierr = PetscFree(scat_g);CHKERRQ(ierr);
328203dadc2fSPeter Brune   }
328303dadc2fSPeter Brune   PetscFunctionReturn(0);
328403dadc2fSPeter Brune }
328503dadc2fSPeter Brune 
328603dadc2fSPeter Brune #undef __FUNCT__
32876636e97aSMatthew G Knepley #define __FUNCT__ "DMSetCoordinates"
32886636e97aSMatthew G Knepley /*@
32896636e97aSMatthew G Knepley   DMSetCoordinates - Sets into the DM a global vector that holds the coordinates
32906636e97aSMatthew G Knepley 
32916636e97aSMatthew G Knepley   Collective on DM
32926636e97aSMatthew G Knepley 
32936636e97aSMatthew G Knepley   Input Parameters:
32946636e97aSMatthew G Knepley + dm - the DM
32956636e97aSMatthew G Knepley - c - coordinate vector
32966636e97aSMatthew G Knepley 
32976636e97aSMatthew G Knepley   Note:
32986636e97aSMatthew G Knepley   The coordinates do include those for ghost points, which are in the local vector
32996636e97aSMatthew G Knepley 
33006636e97aSMatthew G Knepley   Level: intermediate
33016636e97aSMatthew G Knepley 
33026636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
33036636e97aSMatthew G Knepley .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLoca(), DMGetCoordinateDM()
33046636e97aSMatthew G Knepley @*/
33056636e97aSMatthew G Knepley PetscErrorCode DMSetCoordinates(DM dm, Vec c)
33066636e97aSMatthew G Knepley {
33076636e97aSMatthew G Knepley   PetscErrorCode ierr;
33086636e97aSMatthew G Knepley 
33096636e97aSMatthew G Knepley   PetscFunctionBegin;
33106636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
33116636e97aSMatthew G Knepley   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
33126636e97aSMatthew G Knepley   ierr            = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
33136636e97aSMatthew G Knepley   ierr            = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
33146636e97aSMatthew G Knepley   dm->coordinates = c;
33156636e97aSMatthew G Knepley   ierr            = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
3316b64e0483SPeter Brune   ierr            = DMCoarsenHookAdd(dm,DMRestrictHook_Coordinates,NULL,NULL);CHKERRQ(ierr);
331703dadc2fSPeter Brune   ierr            = DMSubDomainHookAdd(dm,DMSubDomainHook_Coordinates,NULL,NULL);CHKERRQ(ierr);
33186636e97aSMatthew G Knepley   PetscFunctionReturn(0);
33196636e97aSMatthew G Knepley }
33206636e97aSMatthew G Knepley 
33216636e97aSMatthew G Knepley #undef __FUNCT__
33226636e97aSMatthew G Knepley #define __FUNCT__ "DMSetCoordinatesLocal"
33236636e97aSMatthew G Knepley /*@
33246636e97aSMatthew G Knepley   DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates
33256636e97aSMatthew G Knepley 
33266636e97aSMatthew G Knepley   Collective on DM
33276636e97aSMatthew G Knepley 
33286636e97aSMatthew G Knepley    Input Parameters:
33296636e97aSMatthew G Knepley +  dm - the DM
33306636e97aSMatthew G Knepley -  c - coordinate vector
33316636e97aSMatthew G Knepley 
33326636e97aSMatthew G Knepley   Note:
33336636e97aSMatthew G Knepley   The coordinates of ghost points can be set using DMSetCoordinates()
33346636e97aSMatthew G Knepley   followed by DMGetCoordinatesLocal(). This is intended to enable the
33356636e97aSMatthew G Knepley   setting of ghost coordinates outside of the domain.
33366636e97aSMatthew G Knepley 
33376636e97aSMatthew G Knepley   Level: intermediate
33386636e97aSMatthew G Knepley 
33396636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
33406636e97aSMatthew G Knepley .seealso: DMGetCoordinatesLocal(), DMSetCoordinates(), DMGetCoordinates(), DMGetCoordinateDM()
33416636e97aSMatthew G Knepley @*/
33426636e97aSMatthew G Knepley PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
33436636e97aSMatthew G Knepley {
33446636e97aSMatthew G Knepley   PetscErrorCode ierr;
33456636e97aSMatthew G Knepley 
33466636e97aSMatthew G Knepley   PetscFunctionBegin;
33476636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
33486636e97aSMatthew G Knepley   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
33496636e97aSMatthew G Knepley   ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
33506636e97aSMatthew G Knepley   ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
33518865f1eaSKarl Rupp 
33526636e97aSMatthew G Knepley   dm->coordinatesLocal = c;
33538865f1eaSKarl Rupp 
33546636e97aSMatthew G Knepley   ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
33556636e97aSMatthew G Knepley   PetscFunctionReturn(0);
33566636e97aSMatthew G Knepley }
33576636e97aSMatthew G Knepley 
33586636e97aSMatthew G Knepley #undef __FUNCT__
33596636e97aSMatthew G Knepley #define __FUNCT__ "DMGetCoordinates"
33606636e97aSMatthew G Knepley /*@
33616636e97aSMatthew G Knepley   DMGetCoordinates - Gets a global vector with the coordinates associated with the DM.
33626636e97aSMatthew G Knepley 
33636636e97aSMatthew G Knepley   Not Collective
33646636e97aSMatthew G Knepley 
33656636e97aSMatthew G Knepley   Input Parameter:
33666636e97aSMatthew G Knepley . dm - the DM
33676636e97aSMatthew G Knepley 
33686636e97aSMatthew G Knepley   Output Parameter:
33696636e97aSMatthew G Knepley . c - global coordinate vector
33706636e97aSMatthew G Knepley 
33716636e97aSMatthew G Knepley   Note:
33726636e97aSMatthew G Knepley   This is a borrowed reference, so the user should NOT destroy this vector
33736636e97aSMatthew G Knepley 
33746636e97aSMatthew G Knepley   Each process has only the local coordinates (does NOT have the ghost coordinates).
33756636e97aSMatthew G Knepley 
33766636e97aSMatthew G Knepley   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
33776636e97aSMatthew G Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
33786636e97aSMatthew G Knepley 
33796636e97aSMatthew G Knepley   Level: intermediate
33806636e97aSMatthew G Knepley 
33816636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
33826636e97aSMatthew G Knepley .seealso: DMSetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM()
33836636e97aSMatthew G Knepley @*/
33846636e97aSMatthew G Knepley PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
33856636e97aSMatthew G Knepley {
33866636e97aSMatthew G Knepley   PetscErrorCode ierr;
33876636e97aSMatthew G Knepley 
33886636e97aSMatthew G Knepley   PetscFunctionBegin;
33896636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
33906636e97aSMatthew G Knepley   PetscValidPointer(c,2);
33911f588964SMatthew G Knepley   if (!dm->coordinates && dm->coordinatesLocal) {
33920298fd71SBarry Smith     DM cdm = NULL;
33936636e97aSMatthew G Knepley 
33946636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
33956636e97aSMatthew G Knepley     ierr = DMCreateGlobalVector(cdm, &dm->coordinates);CHKERRQ(ierr);
33966636e97aSMatthew G Knepley     ierr = PetscObjectSetName((PetscObject) dm->coordinates, "coordinates");CHKERRQ(ierr);
33976636e97aSMatthew G Knepley     ierr = DMLocalToGlobalBegin(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
33986636e97aSMatthew G Knepley     ierr = DMLocalToGlobalEnd(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
33996636e97aSMatthew G Knepley   }
34006636e97aSMatthew G Knepley   *c = dm->coordinates;
34016636e97aSMatthew G Knepley   PetscFunctionReturn(0);
34026636e97aSMatthew G Knepley }
34036636e97aSMatthew G Knepley 
34046636e97aSMatthew G Knepley #undef __FUNCT__
34056636e97aSMatthew G Knepley #define __FUNCT__ "DMGetCoordinatesLocal"
34066636e97aSMatthew G Knepley /*@
34076636e97aSMatthew G Knepley   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM.
34086636e97aSMatthew G Knepley 
34096636e97aSMatthew G Knepley   Collective on DM
34106636e97aSMatthew G Knepley 
34116636e97aSMatthew G Knepley   Input Parameter:
34126636e97aSMatthew G Knepley . dm - the DM
34136636e97aSMatthew G Knepley 
34146636e97aSMatthew G Knepley   Output Parameter:
34156636e97aSMatthew G Knepley . c - coordinate vector
34166636e97aSMatthew G Knepley 
34176636e97aSMatthew G Knepley   Note:
34186636e97aSMatthew G Knepley   This is a borrowed reference, so the user should NOT destroy this vector
34196636e97aSMatthew G Knepley 
34206636e97aSMatthew G Knepley   Each process has the local and ghost coordinates
34216636e97aSMatthew G Knepley 
34226636e97aSMatthew G Knepley   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
34236636e97aSMatthew G Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
34246636e97aSMatthew G Knepley 
34256636e97aSMatthew G Knepley   Level: intermediate
34266636e97aSMatthew G Knepley 
34276636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
34286636e97aSMatthew G Knepley .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMSetCoordinates(), DMGetCoordinateDM()
34296636e97aSMatthew G Knepley @*/
34306636e97aSMatthew G Knepley PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
34316636e97aSMatthew G Knepley {
34326636e97aSMatthew G Knepley   PetscErrorCode ierr;
34336636e97aSMatthew G Knepley 
34346636e97aSMatthew G Knepley   PetscFunctionBegin;
34356636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
34366636e97aSMatthew G Knepley   PetscValidPointer(c,2);
34371f588964SMatthew G Knepley   if (!dm->coordinatesLocal && dm->coordinates) {
34380298fd71SBarry Smith     DM cdm = NULL;
34396636e97aSMatthew G Knepley 
34406636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
34416636e97aSMatthew G Knepley     ierr = DMCreateLocalVector(cdm, &dm->coordinatesLocal);CHKERRQ(ierr);
34426636e97aSMatthew G Knepley     ierr = PetscObjectSetName((PetscObject) dm->coordinatesLocal, "coordinates");CHKERRQ(ierr);
34436636e97aSMatthew G Knepley     ierr = DMGlobalToLocalBegin(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
34446636e97aSMatthew G Knepley     ierr = DMGlobalToLocalEnd(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
34456636e97aSMatthew G Knepley   }
34466636e97aSMatthew G Knepley   *c = dm->coordinatesLocal;
34476636e97aSMatthew G Knepley   PetscFunctionReturn(0);
34486636e97aSMatthew G Knepley }
34496636e97aSMatthew G Knepley 
34506636e97aSMatthew G Knepley #undef __FUNCT__
34516636e97aSMatthew G Knepley #define __FUNCT__ "DMGetCoordinateDM"
34526636e97aSMatthew G Knepley /*@
34531cfe2091SMatthew G. Knepley   DMGetCoordinateDM - Gets the DM that prescribes coordinate layout and scatters between global and local coordinates
34546636e97aSMatthew G Knepley 
34556636e97aSMatthew G Knepley   Collective on DM
34566636e97aSMatthew G Knepley 
34576636e97aSMatthew G Knepley   Input Parameter:
34586636e97aSMatthew G Knepley . dm - the DM
34596636e97aSMatthew G Knepley 
34606636e97aSMatthew G Knepley   Output Parameter:
34616636e97aSMatthew G Knepley . cdm - coordinate DM
34626636e97aSMatthew G Knepley 
34636636e97aSMatthew G Knepley   Level: intermediate
34646636e97aSMatthew G Knepley 
34656636e97aSMatthew G Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
34661cfe2091SMatthew G. Knepley .seealso: DMSetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
34676636e97aSMatthew G Knepley @*/
34686636e97aSMatthew G Knepley PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
34696636e97aSMatthew G Knepley {
34706636e97aSMatthew G Knepley   PetscErrorCode ierr;
34716636e97aSMatthew G Knepley 
34726636e97aSMatthew G Knepley   PetscFunctionBegin;
34736636e97aSMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
34746636e97aSMatthew G Knepley   PetscValidPointer(cdm,2);
34756636e97aSMatthew G Knepley   if (!dm->coordinateDM) {
347682f516ccSBarry Smith     if (!dm->ops->createcoordinatedm) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unable to create coordinates for this DM");
34776636e97aSMatthew G Knepley     ierr = (*dm->ops->createcoordinatedm)(dm, &dm->coordinateDM);CHKERRQ(ierr);
34786636e97aSMatthew G Knepley   }
34796636e97aSMatthew G Knepley   *cdm = dm->coordinateDM;
34806636e97aSMatthew G Knepley   PetscFunctionReturn(0);
34816636e97aSMatthew G Knepley }
3482e87bb0d3SMatthew G Knepley 
3483e87bb0d3SMatthew G Knepley #undef __FUNCT__
34841cfe2091SMatthew G. Knepley #define __FUNCT__ "DMSetCoordinateDM"
34851cfe2091SMatthew G. Knepley /*@
34861cfe2091SMatthew G. Knepley   DMSetCoordinateDM - Sets the DM that prescribes coordinate layout and scatters between global and local coordinates
34871cfe2091SMatthew G. Knepley 
34881cfe2091SMatthew G. Knepley   Logically Collective on DM
34891cfe2091SMatthew G. Knepley 
34901cfe2091SMatthew G. Knepley   Input Parameters:
34911cfe2091SMatthew G. Knepley + dm - the DM
34921cfe2091SMatthew G. Knepley - cdm - coordinate DM
34931cfe2091SMatthew G. Knepley 
34941cfe2091SMatthew G. Knepley   Level: intermediate
34951cfe2091SMatthew G. Knepley 
34961cfe2091SMatthew G. Knepley .keywords: distributed array, get, corners, nodes, local indices, coordinates
34971cfe2091SMatthew G. Knepley .seealso: DMGetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
34981cfe2091SMatthew G. Knepley @*/
34991cfe2091SMatthew G. Knepley PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm)
35001cfe2091SMatthew G. Knepley {
35011cfe2091SMatthew G. Knepley   PetscErrorCode ierr;
35021cfe2091SMatthew G. Knepley 
35031cfe2091SMatthew G. Knepley   PetscFunctionBegin;
35041cfe2091SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
35051cfe2091SMatthew G. Knepley   PetscValidHeaderSpecific(cdm,DM_CLASSID,2);
35061cfe2091SMatthew G. Knepley   ierr = DMDestroy(&dm->coordinateDM);CHKERRQ(ierr);
35071cfe2091SMatthew G. Knepley   dm->coordinateDM = cdm;
35081cfe2091SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) dm->coordinateDM);CHKERRQ(ierr);
35091cfe2091SMatthew G. Knepley   PetscFunctionReturn(0);
35101cfe2091SMatthew G. Knepley }
35111cfe2091SMatthew G. Knepley 
35121cfe2091SMatthew G. Knepley #undef __FUNCT__
3513e8abe2deSMatthew G. Knepley #define __FUNCT__ "DMGetCoordinateSection"
3514e8abe2deSMatthew G. Knepley /*@
3515e8abe2deSMatthew G. Knepley   DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh.
3516e8abe2deSMatthew G. Knepley 
3517e8abe2deSMatthew G. Knepley   Not Collective
3518e8abe2deSMatthew G. Knepley 
3519e8abe2deSMatthew G. Knepley   Input Parameter:
3520e8abe2deSMatthew G. Knepley . dm - The DM object
3521e8abe2deSMatthew G. Knepley 
3522e8abe2deSMatthew G. Knepley   Output Parameter:
3523e8abe2deSMatthew G. Knepley . section - The PetscSection object
3524e8abe2deSMatthew G. Knepley 
3525e8abe2deSMatthew G. Knepley   Level: intermediate
3526e8abe2deSMatthew G. Knepley 
3527e8abe2deSMatthew G. Knepley .keywords: mesh, coordinates
3528e8abe2deSMatthew G. Knepley .seealso: DMGetCoordinateDM(), DMGetDefaultSection(), DMSetDefaultSection()
3529e8abe2deSMatthew G. Knepley @*/
3530e8abe2deSMatthew G. Knepley PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section)
3531e8abe2deSMatthew G. Knepley {
3532e8abe2deSMatthew G. Knepley   DM             cdm;
3533e8abe2deSMatthew G. Knepley   PetscErrorCode ierr;
3534e8abe2deSMatthew G. Knepley 
3535e8abe2deSMatthew G. Knepley   PetscFunctionBegin;
3536e8abe2deSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3537e8abe2deSMatthew G. Knepley   PetscValidPointer(section, 2);
3538e8abe2deSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3539e8abe2deSMatthew G. Knepley   ierr = DMGetDefaultSection(cdm, section);CHKERRQ(ierr);
3540e8abe2deSMatthew G. Knepley   PetscFunctionReturn(0);
3541e8abe2deSMatthew G. Knepley }
3542e8abe2deSMatthew G. Knepley 
3543e8abe2deSMatthew G. Knepley #undef __FUNCT__
3544e8abe2deSMatthew G. Knepley #define __FUNCT__ "DMSetCoordinateSection"
3545e8abe2deSMatthew G. Knepley /*@
3546e8abe2deSMatthew G. Knepley   DMSetCoordinateSection - Set the layout of coordinate values over the mesh.
3547e8abe2deSMatthew G. Knepley 
3548e8abe2deSMatthew G. Knepley   Not Collective
3549e8abe2deSMatthew G. Knepley 
3550e8abe2deSMatthew G. Knepley   Input Parameters:
3551e8abe2deSMatthew G. Knepley + dm      - The DM object
3552e8abe2deSMatthew G. Knepley - section - The PetscSection object
3553e8abe2deSMatthew G. Knepley 
3554e8abe2deSMatthew G. Knepley   Level: intermediate
3555e8abe2deSMatthew G. Knepley 
3556e8abe2deSMatthew G. Knepley .keywords: mesh, coordinates
3557e8abe2deSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetDefaultSection(), DMSetDefaultSection()
3558e8abe2deSMatthew G. Knepley @*/
3559e8abe2deSMatthew G. Knepley PetscErrorCode DMSetCoordinateSection(DM dm, PetscSection section)
3560e8abe2deSMatthew G. Knepley {
3561e8abe2deSMatthew G. Knepley   DM             cdm;
3562e8abe2deSMatthew G. Knepley   PetscErrorCode ierr;
3563e8abe2deSMatthew G. Knepley 
3564e8abe2deSMatthew G. Knepley   PetscFunctionBegin;
3565e8abe2deSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3566e8abe2deSMatthew G. Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
3567e8abe2deSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3568e8abe2deSMatthew G. Knepley   ierr = DMSetDefaultSection(cdm, section);CHKERRQ(ierr);
3569e8abe2deSMatthew G. Knepley   PetscFunctionReturn(0);
3570e8abe2deSMatthew G. Knepley }
3571e8abe2deSMatthew G. Knepley 
3572e8abe2deSMatthew G. Knepley #undef __FUNCT__
3573c6b900c6SMatthew G. Knepley #define __FUNCT__ "DMGetPeriodicity"
3574c6b900c6SMatthew G. Knepley PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal **maxCell, const PetscReal **L)
3575c6b900c6SMatthew G. Knepley {
3576c6b900c6SMatthew G. Knepley   PetscFunctionBegin;
3577c6b900c6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3578c6b900c6SMatthew G. Knepley   if (L)       *L       = dm->L;
3579c6b900c6SMatthew G. Knepley   if (maxCell) *maxCell = dm->maxCell;
3580c6b900c6SMatthew G. Knepley   PetscFunctionReturn(0);
3581c6b900c6SMatthew G. Knepley }
3582c6b900c6SMatthew G. Knepley 
3583c6b900c6SMatthew G. Knepley #undef __FUNCT__
3584c6b900c6SMatthew G. Knepley #define __FUNCT__ "DMSetPeriodicity"
3585c6b900c6SMatthew G. Knepley PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal L[])
3586c6b900c6SMatthew G. Knepley {
3587c6b900c6SMatthew G. Knepley   Vec            coordinates;
3588c6b900c6SMatthew G. Knepley   PetscInt       dim, d;
3589c6b900c6SMatthew G. Knepley   PetscErrorCode ierr;
3590c6b900c6SMatthew G. Knepley 
3591c6b900c6SMatthew G. Knepley   PetscFunctionBegin;
3592c6b900c6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3593c6b900c6SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
3594c6b900c6SMatthew G. Knepley   ierr = VecGetBlockSize(coordinates, &dim);CHKERRQ(ierr);
3595c6b900c6SMatthew G. Knepley   ierr = PetscMalloc1(dim,&dm->L);CHKERRQ(ierr);
3596c6b900c6SMatthew G. Knepley   ierr = PetscMalloc1(dim,&dm->maxCell);CHKERRQ(ierr);
3597c6b900c6SMatthew G. Knepley   for (d = 0; d < dim; ++d) {dm->L[d] = L[d]; dm->maxCell[d] = maxCell[d];}
3598c6b900c6SMatthew G. Knepley   PetscFunctionReturn(0);
3599c6b900c6SMatthew G. Knepley }
3600c6b900c6SMatthew G. Knepley 
3601c6b900c6SMatthew G. Knepley #undef __FUNCT__
3602e87bb0d3SMatthew G Knepley #define __FUNCT__ "DMLocatePoints"
3603e87bb0d3SMatthew G Knepley /*@
3604e87bb0d3SMatthew G Knepley   DMLocatePoints - Locate the points in v in the mesh and return an IS of the containing cells
3605e87bb0d3SMatthew G Knepley 
3606e87bb0d3SMatthew G Knepley   Not collective
3607e87bb0d3SMatthew G Knepley 
3608e87bb0d3SMatthew G Knepley   Input Parameters:
3609e87bb0d3SMatthew G Knepley + dm - The DM
3610e87bb0d3SMatthew G Knepley - v - The Vec of points
3611e87bb0d3SMatthew G Knepley 
361261e3bb9bSMatthew G Knepley   Output Parameter:
3613e87bb0d3SMatthew G Knepley . cells - The local cell numbers for cells which contain the points
3614e87bb0d3SMatthew G Knepley 
3615e87bb0d3SMatthew G Knepley   Level: developer
361661e3bb9bSMatthew G Knepley 
361761e3bb9bSMatthew G Knepley .keywords: point location, mesh
361861e3bb9bSMatthew G Knepley .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
361961e3bb9bSMatthew G Knepley @*/
3620e87bb0d3SMatthew G Knepley PetscErrorCode DMLocatePoints(DM dm, Vec v, IS *cells)
3621e87bb0d3SMatthew G Knepley {
3622735aa83eSMatthew G Knepley   PetscErrorCode ierr;
3623735aa83eSMatthew G Knepley 
3624e87bb0d3SMatthew G Knepley   PetscFunctionBegin;
3625e87bb0d3SMatthew G Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3626e87bb0d3SMatthew G Knepley   PetscValidHeaderSpecific(v,VEC_CLASSID,2);
3627e87bb0d3SMatthew G Knepley   PetscValidPointer(cells,3);
3628735aa83eSMatthew G Knepley   if (dm->ops->locatepoints) {
3629735aa83eSMatthew G Knepley     ierr = (*dm->ops->locatepoints)(dm,v,cells);CHKERRQ(ierr);
363082f516ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Point location not available for this DM");
3631e87bb0d3SMatthew G Knepley   PetscFunctionReturn(0);
3632e87bb0d3SMatthew G Knepley }
363314f150ffSMatthew G. Knepley 
363414f150ffSMatthew G. Knepley #undef __FUNCT__
363514f150ffSMatthew G. Knepley #define __FUNCT__ "DMGetOutputDM"
3636f4d763aaSMatthew G. Knepley /*@
3637f4d763aaSMatthew G. Knepley   DMGetOutputDM - Retrieve the DM associated with the layout for output
3638f4d763aaSMatthew G. Knepley 
3639f4d763aaSMatthew G. Knepley   Input Parameter:
3640f4d763aaSMatthew G. Knepley . dm - The original DM
3641f4d763aaSMatthew G. Knepley 
3642f4d763aaSMatthew G. Knepley   Output Parameter:
3643f4d763aaSMatthew G. Knepley . odm - The DM which provides the layout for output
3644f4d763aaSMatthew G. Knepley 
3645f4d763aaSMatthew G. Knepley   Level: intermediate
3646f4d763aaSMatthew G. Knepley 
3647f4d763aaSMatthew G. Knepley .seealso: VecView(), DMGetDefaultGlobalSection()
3648f4d763aaSMatthew G. Knepley @*/
364914f150ffSMatthew G. Knepley PetscErrorCode DMGetOutputDM(DM dm, DM *odm)
365014f150ffSMatthew G. Knepley {
3651c26acbdeSMatthew G. Knepley   PetscSection   section;
3652c26acbdeSMatthew G. Knepley   PetscBool      hasConstraints;
365314f150ffSMatthew G. Knepley   PetscErrorCode ierr;
365414f150ffSMatthew G. Knepley 
365514f150ffSMatthew G. Knepley   PetscFunctionBegin;
365614f150ffSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
365714f150ffSMatthew G. Knepley   PetscValidPointer(odm,2);
3658c26acbdeSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
3659c26acbdeSMatthew G. Knepley   ierr = PetscSectionHasConstraints(section, &hasConstraints);CHKERRQ(ierr);
3660c26acbdeSMatthew G. Knepley   if (!hasConstraints) {
3661c26acbdeSMatthew G. Knepley     *odm = dm;
3662c26acbdeSMatthew G. Knepley     PetscFunctionReturn(0);
3663c26acbdeSMatthew G. Knepley   }
366414f150ffSMatthew G. Knepley   if (!dm->dmBC) {
3665c26acbdeSMatthew G. Knepley     PetscSection newSection, gsection;
366614f150ffSMatthew G. Knepley     PetscSF      sf;
366714f150ffSMatthew G. Knepley 
366814f150ffSMatthew G. Knepley     ierr = DMClone(dm, &dm->dmBC);CHKERRQ(ierr);
366914f150ffSMatthew G. Knepley     ierr = PetscSectionClone(section, &newSection);CHKERRQ(ierr);
367014f150ffSMatthew G. Knepley     ierr = DMSetDefaultSection(dm->dmBC, newSection);CHKERRQ(ierr);
367114f150ffSMatthew G. Knepley     ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);
367214f150ffSMatthew G. Knepley     ierr = DMGetPointSF(dm->dmBC, &sf);CHKERRQ(ierr);
367314f150ffSMatthew G. Knepley     ierr = PetscSectionCreateGlobalSection(section, sf, PETSC_TRUE, &gsection);CHKERRQ(ierr);
367414f150ffSMatthew G. Knepley     ierr = DMSetDefaultGlobalSection(dm->dmBC, gsection);CHKERRQ(ierr);
367514f150ffSMatthew G. Knepley     ierr = PetscSectionDestroy(&gsection);CHKERRQ(ierr);
367614f150ffSMatthew G. Knepley   }
367714f150ffSMatthew G. Knepley   *odm = dm->dmBC;
367814f150ffSMatthew G. Knepley   PetscFunctionReturn(0);
367914f150ffSMatthew G. Knepley }
3680f4d763aaSMatthew G. Knepley 
3681f4d763aaSMatthew G. Knepley #undef __FUNCT__
3682f4d763aaSMatthew G. Knepley #define __FUNCT__ "DMGetOutputSequenceNumber"
3683f4d763aaSMatthew G. Knepley /*@
3684f4d763aaSMatthew G. Knepley   DMGetOutputSequenceNumber - Retrieve the sequence number for output
3685f4d763aaSMatthew G. Knepley 
3686f4d763aaSMatthew G. Knepley   Input Parameter:
3687f4d763aaSMatthew G. Knepley . dm - The original DM
3688f4d763aaSMatthew G. Knepley 
3689f4d763aaSMatthew G. Knepley   Output Parameter:
3690f4d763aaSMatthew G. Knepley . num - The output sequence number
3691f4d763aaSMatthew G. Knepley 
3692f4d763aaSMatthew G. Knepley   Level: intermediate
3693f4d763aaSMatthew G. Knepley 
3694f4d763aaSMatthew G. Knepley   Note: This is intended for output that should appear in sequence, for instance
3695f4d763aaSMatthew G. Knepley   a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system.
3696f4d763aaSMatthew G. Knepley 
3697f4d763aaSMatthew G. Knepley .seealso: VecView()
3698f4d763aaSMatthew G. Knepley @*/
3699f4d763aaSMatthew G. Knepley PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num)
3700f4d763aaSMatthew G. Knepley {
3701f4d763aaSMatthew G. Knepley   PetscFunctionBegin;
3702f4d763aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3703f4d763aaSMatthew G. Knepley   PetscValidPointer(num,2);
3704f4d763aaSMatthew G. Knepley   *num = dm->outputSequenceNum;
3705f4d763aaSMatthew G. Knepley   PetscFunctionReturn(0);
3706f4d763aaSMatthew G. Knepley }
3707f4d763aaSMatthew G. Knepley 
3708f4d763aaSMatthew G. Knepley #undef __FUNCT__
3709f4d763aaSMatthew G. Knepley #define __FUNCT__ "DMSetOutputSequenceNumber"
3710f4d763aaSMatthew G. Knepley /*@
3711f4d763aaSMatthew G. Knepley   DMSetOutputSequenceNumber - Set the sequence number for output
3712f4d763aaSMatthew G. Knepley 
3713f4d763aaSMatthew G. Knepley   Input Parameters:
3714f4d763aaSMatthew G. Knepley + dm - The original DM
3715f4d763aaSMatthew G. Knepley - num - The output sequence number
3716f4d763aaSMatthew G. Knepley 
3717f4d763aaSMatthew G. Knepley   Level: intermediate
3718f4d763aaSMatthew G. Knepley 
3719f4d763aaSMatthew G. Knepley   Note: This is intended for output that should appear in sequence, for instance
3720f4d763aaSMatthew G. Knepley   a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system.
3721f4d763aaSMatthew G. Knepley 
3722f4d763aaSMatthew G. Knepley .seealso: VecView()
3723f4d763aaSMatthew G. Knepley @*/
3724f4d763aaSMatthew G. Knepley PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num)
3725f4d763aaSMatthew G. Knepley {
3726f4d763aaSMatthew G. Knepley   PetscFunctionBegin;
3727f4d763aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3728f4d763aaSMatthew G. Knepley   dm->outputSequenceNum = num;
3729f4d763aaSMatthew G. Knepley   PetscFunctionReturn(0);
3730f4d763aaSMatthew G. Knepley }
3731