xref: /petsc/src/dm/interface/dm.c (revision f089877a000c0949609bd4a1bfd01a0aa8158bb4)
1 #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/
2 #include <petscsf.h>
3 
4 PetscClassId  DM_CLASSID;
5 PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocalToLocal;
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "DMViewFromOptions"
9 /*
10   DMViewFromOptions - Processes command line options to determine if/how a DM is to be viewed.
11 
12   Collective on Vec
13 
14   Input Parameters:
15 + dm   - the DM
16 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
17 - optionname - option to activate viewing
18 
19   Level: intermediate
20 
21 .keywords: DM, view, options, database
22 .seealso: VecViewFromOptions(), MatViewFromOptions()
23 */
24 PetscErrorCode  DMViewFromOptions(DM dm,const char prefix[],const char optionname[])
25 {
26   PetscErrorCode    ierr;
27   PetscBool         flg;
28   PetscViewer       viewer;
29   PetscViewerFormat format;
30 
31   PetscFunctionBegin;
32   if (prefix) {
33     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)dm),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
34   } else {
35     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
36   }
37   if (flg) {
38     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
39     ierr = DMView(dm,viewer);CHKERRQ(ierr);
40     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
41     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "DMCreate"
49 /*@
50   DMCreate - Creates an empty DM object. The type can then be set with DMSetType().
51 
52    If you never  call DMSetType()  it will generate an
53    error when you try to use the vector.
54 
55   Collective on MPI_Comm
56 
57   Input Parameter:
58 . comm - The communicator for the DM object
59 
60   Output Parameter:
61 . dm - The DM object
62 
63   Level: beginner
64 
65 .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE
66 @*/
67 PetscErrorCode  DMCreate(MPI_Comm comm,DM *dm)
68 {
69   DM             v;
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   PetscValidPointer(dm,2);
74   *dm = NULL;
75 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
76   ierr = VecInitializePackage();CHKERRQ(ierr);
77   ierr = MatInitializePackage();CHKERRQ(ierr);
78   ierr = DMInitializePackage();CHKERRQ(ierr);
79 #endif
80 
81   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
82   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
83 
84 
85   v->ltogmap              = NULL;
86   v->ltogmapb             = NULL;
87   v->bs                   = 1;
88   v->coloringtype         = IS_COLORING_GLOBAL;
89   ierr                    = PetscSFCreate(comm, &v->sf);CHKERRQ(ierr);
90   ierr                    = PetscSFCreate(comm, &v->defaultSF);CHKERRQ(ierr);
91   v->defaultSection       = NULL;
92   v->defaultGlobalSection = NULL;
93   {
94     PetscInt i;
95     for (i = 0; i < 10; ++i) {
96       v->nullspaceConstructors[i] = NULL;
97     }
98   }
99   v->numFields = 0;
100   v->fields    = NULL;
101 
102   *dm = v;
103   PetscFunctionReturn(0);
104 }
105 
106 #undef __FUNCT__
107 #define __FUNCT__ "DMSetVecType"
108 /*@C
109        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
110 
111    Logically Collective on DMDA
112 
113    Input Parameter:
114 +  da - initial distributed array
115 .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
116 
117    Options Database:
118 .   -dm_vec_type ctype
119 
120    Level: intermediate
121 
122 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
123 @*/
124 PetscErrorCode  DMSetVecType(DM da,VecType ctype)
125 {
126   PetscErrorCode ierr;
127 
128   PetscFunctionBegin;
129   PetscValidHeaderSpecific(da,DM_CLASSID,1);
130   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
131   ierr = PetscStrallocpy(ctype,(char**)&da->vectype);CHKERRQ(ierr);
132   PetscFunctionReturn(0);
133 }
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "VecGetDM"
137 /*@
138   VecGetDM - Gets the DM defining the data layout of the vector
139 
140   Not collective
141 
142   Input Parameter:
143 . v - The Vec
144 
145   Output Parameter:
146 . dm - The DM
147 
148   Level: intermediate
149 
150 .seealso: VecSetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType()
151 @*/
152 PetscErrorCode VecGetDM(Vec v, DM *dm)
153 {
154   PetscErrorCode ierr;
155 
156   PetscFunctionBegin;
157   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
158   PetscValidPointer(dm,2);
159   ierr = PetscObjectQuery((PetscObject) v, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "VecSetDM"
165 /*@
166   VecSetDM - Sets the DM defining the data layout of the vector
167 
168   Not collective
169 
170   Input Parameters:
171 + v - The Vec
172 - dm - The DM
173 
174   Level: intermediate
175 
176 .seealso: VecGetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType()
177 @*/
178 PetscErrorCode VecSetDM(Vec v, DM dm)
179 {
180   PetscErrorCode ierr;
181 
182   PetscFunctionBegin;
183   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
184   if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2);
185   ierr = PetscObjectCompose((PetscObject) v, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr);
186   PetscFunctionReturn(0);
187 }
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "DMSetMatType"
191 /*@C
192        DMSetMatType - Sets the type of matrix created with DMCreateMatrix()
193 
194    Logically Collective on DM
195 
196    Input Parameter:
197 +  dm - the DM context
198 .  ctype - the matrix type
199 
200    Options Database:
201 .   -dm_mat_type ctype
202 
203    Level: intermediate
204 
205 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType
206 @*/
207 PetscErrorCode  DMSetMatType(DM dm,MatType ctype)
208 {
209   PetscErrorCode ierr;
210 
211   PetscFunctionBegin;
212   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
213   ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
214   ierr = PetscStrallocpy(ctype,(char**)&dm->mattype);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "MatGetDM"
220 /*@
221   MatGetDM - Gets the DM defining the data layout of the matrix
222 
223   Not collective
224 
225   Input Parameter:
226 . A - The Mat
227 
228   Output Parameter:
229 . dm - The DM
230 
231   Level: intermediate
232 
233 .seealso: MatSetDM(), DMCreateMatrix(), DMSetMatType()
234 @*/
235 PetscErrorCode MatGetDM(Mat A, DM *dm)
236 {
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
241   PetscValidPointer(dm,2);
242   ierr = PetscObjectQuery((PetscObject) A, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "MatSetDM"
248 /*@
249   MatSetDM - Sets the DM defining the data layout of the matrix
250 
251   Not collective
252 
253   Input Parameters:
254 + A - The Mat
255 - dm - The DM
256 
257   Level: intermediate
258 
259 .seealso: MatGetDM(), DMCreateMatrix(), DMSetMatType()
260 @*/
261 PetscErrorCode MatSetDM(Mat A, DM dm)
262 {
263   PetscErrorCode ierr;
264 
265   PetscFunctionBegin;
266   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
267   if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2);
268   ierr = PetscObjectCompose((PetscObject) A, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr);
269   PetscFunctionReturn(0);
270 }
271 
272 #undef __FUNCT__
273 #define __FUNCT__ "DMSetOptionsPrefix"
274 /*@C
275    DMSetOptionsPrefix - Sets the prefix used for searching for all
276    DMDA options in the database.
277 
278    Logically Collective on DMDA
279 
280    Input Parameter:
281 +  da - the DMDA context
282 -  prefix - the prefix to prepend to all option names
283 
284    Notes:
285    A hyphen (-) must NOT be given at the beginning of the prefix name.
286    The first character of all runtime options is AUTOMATICALLY the hyphen.
287 
288    Level: advanced
289 
290 .keywords: DMDA, set, options, prefix, database
291 
292 .seealso: DMSetFromOptions()
293 @*/
294 PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
295 {
296   PetscErrorCode ierr;
297 
298   PetscFunctionBegin;
299   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
300   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "DMDestroy"
306 /*@
307     DMDestroy - Destroys a vector packer or DMDA.
308 
309     Collective on DM
310 
311     Input Parameter:
312 .   dm - the DM object to destroy
313 
314     Level: developer
315 
316 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
317 
318 @*/
319 PetscErrorCode  DMDestroy(DM *dm)
320 {
321   PetscInt       i, cnt = 0, f;
322   DMNamedVecLink nlink,nnext;
323   PetscErrorCode ierr;
324 
325   PetscFunctionBegin;
326   if (!*dm) PetscFunctionReturn(0);
327   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
328 
329   /* count all the circular references of DM and its contained Vecs */
330   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
331     if ((*dm)->localin[i])  cnt++;
332     if ((*dm)->globalin[i]) cnt++;
333   }
334   for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++;
335   if ((*dm)->x) {
336     DM obj;
337     ierr = VecGetDM((*dm)->x, &obj);CHKERRQ(ierr);
338     if (obj == *dm) cnt++;
339   }
340   for (nlink=(*dm)->namedlocal; nlink; nlink=nlink->next) cnt++;
341   if ((*dm)->x) {
342     DM obj;
343     ierr = VecGetDM((*dm)->x, &obj);CHKERRQ(ierr);
344     if (obj == *dm) cnt++;
345   }
346 
347   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
348   /*
349      Need this test because the dm references the vectors that
350      reference the dm, so destroying the dm calls destroy on the
351      vectors that cause another destroy on the dm
352   */
353   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
354   ((PetscObject) (*dm))->refct = 0;
355   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
356     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
357     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
358   }
359   for (nlink=(*dm)->namedglobal; nlink; nlink=nnext) { /* Destroy the named vectors */
360     nnext = nlink->next;
361     if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name);
362     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
363     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
364     ierr = PetscFree(nlink);CHKERRQ(ierr);
365   }
366   (*dm)->namedglobal = NULL;
367 
368   for (nlink=(*dm)->namedlocal; nlink; nlink=nnext) { /* Destroy the named local vectors */
369     nnext = nlink->next;
370     if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name);
371     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
372     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
373     ierr = PetscFree(nlink);CHKERRQ(ierr);
374   }
375   (*dm)->namedlocal = NULL;
376 
377   /* Destroy the list of hooks */
378   {
379     DMCoarsenHookLink link,next;
380     for (link=(*dm)->coarsenhook; link; link=next) {
381       next = link->next;
382       ierr = PetscFree(link);CHKERRQ(ierr);
383     }
384     (*dm)->coarsenhook = NULL;
385   }
386   {
387     DMRefineHookLink link,next;
388     for (link=(*dm)->refinehook; link; link=next) {
389       next = link->next;
390       ierr = PetscFree(link);CHKERRQ(ierr);
391     }
392     (*dm)->refinehook = NULL;
393   }
394   {
395     DMSubDomainHookLink link,next;
396     for (link=(*dm)->subdomainhook; link; link=next) {
397       next = link->next;
398       ierr = PetscFree(link);CHKERRQ(ierr);
399     }
400     (*dm)->subdomainhook = NULL;
401   }
402   {
403     DMGlobalToLocalHookLink link,next;
404     for (link=(*dm)->gtolhook; link; link=next) {
405       next = link->next;
406       ierr = PetscFree(link);CHKERRQ(ierr);
407     }
408     (*dm)->gtolhook = NULL;
409   }
410   /* Destroy the work arrays */
411   {
412     DMWorkLink link,next;
413     if ((*dm)->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out");
414     for (link=(*dm)->workin; link; link=next) {
415       next = link->next;
416       ierr = PetscFree(link->mem);CHKERRQ(ierr);
417       ierr = PetscFree(link);CHKERRQ(ierr);
418     }
419     (*dm)->workin = NULL;
420   }
421 
422   ierr = PetscObjectDestroy(&(*dm)->dmksp);CHKERRQ(ierr);
423   ierr = PetscObjectDestroy(&(*dm)->dmsnes);CHKERRQ(ierr);
424   ierr = PetscObjectDestroy(&(*dm)->dmts);CHKERRQ(ierr);
425 
426   if ((*dm)->ctx && (*dm)->ctxdestroy) {
427     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
428   }
429   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
430   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
431   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
432   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
433   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
434   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
435   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
436 
437   ierr = PetscSectionDestroy(&(*dm)->defaultSection);CHKERRQ(ierr);
438   ierr = PetscSectionDestroy(&(*dm)->defaultGlobalSection);CHKERRQ(ierr);
439   ierr = PetscLayoutDestroy(&(*dm)->map);CHKERRQ(ierr);
440   ierr = PetscSFDestroy(&(*dm)->sf);CHKERRQ(ierr);
441   ierr = PetscSFDestroy(&(*dm)->defaultSF);CHKERRQ(ierr);
442 
443   ierr = DMDestroy(&(*dm)->coordinateDM);CHKERRQ(ierr);
444   ierr = VecDestroy(&(*dm)->coordinates);CHKERRQ(ierr);
445   ierr = VecDestroy(&(*dm)->coordinatesLocal);CHKERRQ(ierr);
446 
447   for (f = 0; f < (*dm)->numFields; ++f) {
448     ierr = PetscObjectDestroy(&(*dm)->fields[f]);CHKERRQ(ierr);
449   }
450   ierr = PetscFree((*dm)->fields);CHKERRQ(ierr);
451   /* if memory was published with AMS then destroy it */
452   ierr = PetscObjectAMSViewOff((PetscObject)*dm);CHKERRQ(ierr);
453 
454   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
455   /* We do not destroy (*dm)->data here so that we can reference count backend objects */
456   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
457   PetscFunctionReturn(0);
458 }
459 
460 #undef __FUNCT__
461 #define __FUNCT__ "DMSetUp"
462 /*@
463     DMSetUp - sets up the data structures inside a DM object
464 
465     Collective on DM
466 
467     Input Parameter:
468 .   dm - the DM object to setup
469 
470     Level: developer
471 
472 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
473 
474 @*/
475 PetscErrorCode  DMSetUp(DM dm)
476 {
477   PetscErrorCode ierr;
478 
479   PetscFunctionBegin;
480   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
481   if (dm->setupcalled) PetscFunctionReturn(0);
482   if (dm->ops->setup) {
483     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
484   }
485   dm->setupcalled = PETSC_TRUE;
486   PetscFunctionReturn(0);
487 }
488 
489 #undef __FUNCT__
490 #define __FUNCT__ "DMSetFromOptions"
491 /*@
492     DMSetFromOptions - sets parameters in a DM from the options database
493 
494     Collective on DM
495 
496     Input Parameter:
497 .   dm - the DM object to set options for
498 
499     Options Database:
500 +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
501 .   -dm_vec_type <type>  type of vector to create inside DM
502 .   -dm_mat_type <type>  type of matrix to create inside DM
503 -   -dm_coloring_type <global or ghosted>
504 
505     Level: developer
506 
507 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
508 
509 @*/
510 PetscErrorCode  DMSetFromOptions(DM dm)
511 {
512   char           typeName[256] = MATAIJ;
513   PetscBool      flg;
514   PetscErrorCode ierr;
515 
516   PetscFunctionBegin;
517   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
518   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
519   ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,NULL);CHKERRQ(ierr);
520   ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
521   if (flg) {
522     ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
523   }
524   ierr = PetscOptionsList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype ? dm->mattype : typeName,typeName,sizeof(typeName),&flg);CHKERRQ(ierr);
525   if (flg) {
526     ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr);
527   }
528   ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,NULL);CHKERRQ(ierr);
529   if (dm->ops->setfromoptions) {
530     ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
531   }
532   /* process any options handlers added with PetscObjectAddOptionsHandler() */
533   ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
534   ierr = PetscOptionsEnd();CHKERRQ(ierr);
535   ierr = DMViewFromOptions(dm,NULL,"-dm_view");CHKERRQ(ierr);
536   PetscFunctionReturn(0);
537 }
538 
539 #undef __FUNCT__
540 #define __FUNCT__ "DMView"
541 /*@C
542     DMView - Views a vector packer or DMDA.
543 
544     Collective on DM
545 
546     Input Parameter:
547 +   dm - the DM object to view
548 -   v - the viewer
549 
550     Level: developer
551 
552 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
553 
554 @*/
555 PetscErrorCode  DMView(DM dm,PetscViewer v)
556 {
557   PetscErrorCode ierr;
558   PetscBool      isbinary;
559 
560   PetscFunctionBegin;
561   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
562   if (!v) {
563     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm),&v);CHKERRQ(ierr);
564   }
565   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
566   if (isbinary) {
567     PetscInt classid = DM_FILE_CLASSID;
568     char     type[256];
569 
570     ierr = PetscViewerBinaryWrite(v,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
571     ierr = PetscStrncpy(type,((PetscObject)dm)->type_name,256);CHKERRQ(ierr);
572     ierr = PetscViewerBinaryWrite(v,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
573   }
574   if (dm->ops->view) {
575     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
576   }
577   PetscFunctionReturn(0);
578 }
579 
580 #undef __FUNCT__
581 #define __FUNCT__ "DMCreateGlobalVector"
582 /*@
583     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
584 
585     Collective on DM
586 
587     Input Parameter:
588 .   dm - the DM object
589 
590     Output Parameter:
591 .   vec - the global vector
592 
593     Level: beginner
594 
595 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
596 
597 @*/
598 PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
599 {
600   PetscErrorCode ierr;
601 
602   PetscFunctionBegin;
603   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
604   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
605   PetscFunctionReturn(0);
606 }
607 
608 #undef __FUNCT__
609 #define __FUNCT__ "DMCreateLocalVector"
610 /*@
611     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
612 
613     Not Collective
614 
615     Input Parameter:
616 .   dm - the DM object
617 
618     Output Parameter:
619 .   vec - the local vector
620 
621     Level: beginner
622 
623 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
624 
625 @*/
626 PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
627 {
628   PetscErrorCode ierr;
629 
630   PetscFunctionBegin;
631   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
632   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
633   PetscFunctionReturn(0);
634 }
635 
636 #undef __FUNCT__
637 #define __FUNCT__ "DMGetLocalToGlobalMapping"
638 /*@
639    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
640 
641    Collective on DM
642 
643    Input Parameter:
644 .  dm - the DM that provides the mapping
645 
646    Output Parameter:
647 .  ltog - the mapping
648 
649    Level: intermediate
650 
651    Notes:
652    This mapping can then be used by VecSetLocalToGlobalMapping() or
653    MatSetLocalToGlobalMapping().
654 
655 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
656 @*/
657 PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
658 {
659   PetscErrorCode ierr;
660 
661   PetscFunctionBegin;
662   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
663   PetscValidPointer(ltog,2);
664   if (!dm->ltogmap) {
665     PetscSection section, sectionGlobal;
666 
667     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
668     if (section) {
669       PetscInt *ltog;
670       PetscInt pStart, pEnd, size, p, l;
671 
672       ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
673       ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
674       ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
675       ierr = PetscMalloc(size * sizeof(PetscInt), &ltog);CHKERRQ(ierr); /* We want the local+overlap size */
676       for (p = pStart, l = 0; p < pEnd; ++p) {
677         PetscInt dof, off, c;
678 
679         /* Should probably use constrained dofs */
680         ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
681         ierr = PetscSectionGetOffset(sectionGlobal, p, &off);CHKERRQ(ierr);
682         for (c = 0; c < dof; ++c, ++l) {
683           ltog[l] = off+c;
684         }
685       }
686       ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, size, ltog, PETSC_OWN_POINTER, &dm->ltogmap);CHKERRQ(ierr);
687       ierr = PetscLogObjectParent(dm, dm->ltogmap);CHKERRQ(ierr);
688     } else {
689       if (!dm->ops->getlocaltoglobalmapping) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
690       ierr = (*dm->ops->getlocaltoglobalmapping)(dm);CHKERRQ(ierr);
691     }
692   }
693   *ltog = dm->ltogmap;
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
699 /*@
700    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
701 
702    Collective on DM
703 
704    Input Parameter:
705 .  da - the distributed array that provides the mapping
706 
707    Output Parameter:
708 .  ltog - the block mapping
709 
710    Level: intermediate
711 
712    Notes:
713    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
714    MatSetLocalToGlobalMappingBlock().
715 
716 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
717 @*/
718 PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
719 {
720   PetscErrorCode ierr;
721 
722   PetscFunctionBegin;
723   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
724   PetscValidPointer(ltog,2);
725   if (!dm->ltogmapb) {
726     PetscInt bs;
727     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
728     if (bs > 1) {
729       if (!dm->ops->getlocaltoglobalmappingblock) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
730       ierr = (*dm->ops->getlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
731     } else {
732       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
733       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
734     }
735   }
736   *ltog = dm->ltogmapb;
737   PetscFunctionReturn(0);
738 }
739 
740 #undef __FUNCT__
741 #define __FUNCT__ "DMGetBlockSize"
742 /*@
743    DMGetBlockSize - Gets the inherent block size associated with a DM
744 
745    Not Collective
746 
747    Input Parameter:
748 .  dm - the DM with block structure
749 
750    Output Parameter:
751 .  bs - the block size, 1 implies no exploitable block structure
752 
753    Level: intermediate
754 
755 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
756 @*/
757 PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
758 {
759   PetscFunctionBegin;
760   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
761   PetscValidPointer(bs,2);
762   if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet");
763   *bs = dm->bs;
764   PetscFunctionReturn(0);
765 }
766 
767 #undef __FUNCT__
768 #define __FUNCT__ "DMCreateInterpolation"
769 /*@
770     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
771 
772     Collective on DM
773 
774     Input Parameter:
775 +   dm1 - the DM object
776 -   dm2 - the second, finer DM object
777 
778     Output Parameter:
779 +  mat - the interpolation
780 -  vec - the scaling (optional)
781 
782     Level: developer
783 
784     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
785         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
786 
787         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMGetCoordinateDM()) to interpolate the mesh coordinate vectors
788         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
789 
790 
791 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
792 
793 @*/
794 PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
795 {
796   PetscErrorCode ierr;
797 
798   PetscFunctionBegin;
799   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
800   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
801   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
802   PetscFunctionReturn(0);
803 }
804 
805 #undef __FUNCT__
806 #define __FUNCT__ "DMCreateInjection"
807 /*@
808     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
809 
810     Collective on DM
811 
812     Input Parameter:
813 +   dm1 - the DM object
814 -   dm2 - the second, finer DM object
815 
816     Output Parameter:
817 .   ctx - the injection
818 
819     Level: developer
820 
821    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
822         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
823 
824 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
825 
826 @*/
827 PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
828 {
829   PetscErrorCode ierr;
830 
831   PetscFunctionBegin;
832   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
833   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
834   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
835   PetscFunctionReturn(0);
836 }
837 
838 #undef __FUNCT__
839 #define __FUNCT__ "DMCreateColoring"
840 /*@C
841     DMCreateColoring - Gets coloring for a DMDA or DMComposite
842 
843     Collective on DM
844 
845     Input Parameter:
846 +   dm - the DM object
847 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
848 -   matype - either MATAIJ or MATBAIJ
849 
850     Output Parameter:
851 .   coloring - the coloring
852 
853     Level: developer
854 
855 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
856 
857 @*/
858 PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,MatType mtype,ISColoring *coloring)
859 {
860   PetscErrorCode ierr;
861 
862   PetscFunctionBegin;
863   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
864   if (!dm->ops->getcoloring) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No coloring for this type of DM yet");
865   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
866   PetscFunctionReturn(0);
867 }
868 
869 #undef __FUNCT__
870 #define __FUNCT__ "DMCreateMatrix"
871 /*@C
872     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
873 
874     Collective on DM
875 
876     Input Parameter:
877 +   dm - the DM object
878 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
879             any type which inherits from one of these (such as MATAIJ)
880 
881     Output Parameter:
882 .   mat - the empty Jacobian
883 
884     Level: beginner
885 
886     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
887        do not need to do it yourself.
888 
889        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
890        the nonzero pattern call DMDASetMatPreallocateOnly()
891 
892        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
893        internally by PETSc.
894 
895        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
896        the indices for the global numbering for DMDAs which is complicated.
897 
898 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
899 
900 @*/
901 PetscErrorCode  DMCreateMatrix(DM dm,MatType mtype,Mat *mat)
902 {
903   PetscErrorCode ierr;
904 
905   PetscFunctionBegin;
906   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
907 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
908   ierr = MatInitializePackage();CHKERRQ(ierr);
909 #endif
910   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
911   PetscValidPointer(mat,3);
912   if (dm->mattype) {
913     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
914   } else {
915     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
916   }
917   PetscFunctionReturn(0);
918 }
919 
920 #undef __FUNCT__
921 #define __FUNCT__ "DMSetMatrixPreallocateOnly"
922 /*@
923   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
924     preallocated but the nonzero structure and zero values will not be set.
925 
926   Logically Collective on DMDA
927 
928   Input Parameter:
929 + dm - the DM
930 - only - PETSC_TRUE if only want preallocation
931 
932   Level: developer
933 .seealso DMCreateMatrix()
934 @*/
935 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
936 {
937   PetscFunctionBegin;
938   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
939   dm->prealloc_only = only;
940   PetscFunctionReturn(0);
941 }
942 
943 #undef __FUNCT__
944 #define __FUNCT__ "DMGetWorkArray"
945 /*@C
946   DMGetWorkArray - Gets a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray()
947 
948   Not Collective
949 
950   Input Parameters:
951 + dm - the DM object
952 . count - The minium size
953 - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT)
954 
955   Output Parameter:
956 . array - the work array
957 
958   Level: developer
959 
960 .seealso DMDestroy(), DMCreate()
961 @*/
962 PetscErrorCode DMGetWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem)
963 {
964   PetscErrorCode ierr;
965   DMWorkLink     link;
966   size_t         size;
967 
968   PetscFunctionBegin;
969   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
970   PetscValidPointer(mem,4);
971   if (dm->workin) {
972     link       = dm->workin;
973     dm->workin = dm->workin->next;
974   } else {
975     ierr = PetscNewLog(dm,struct _DMWorkLink,&link);CHKERRQ(ierr);
976   }
977   ierr = PetscDataTypeGetSize(dtype,&size);CHKERRQ(ierr);
978   if (size*count > link->bytes) {
979     ierr        = PetscFree(link->mem);CHKERRQ(ierr);
980     ierr        = PetscMalloc(size*count,&link->mem);CHKERRQ(ierr);
981     link->bytes = size*count;
982   }
983   link->next   = dm->workout;
984   dm->workout  = link;
985   *(void**)mem = link->mem;
986   PetscFunctionReturn(0);
987 }
988 
989 #undef __FUNCT__
990 #define __FUNCT__ "DMRestoreWorkArray"
991 /*@C
992   DMRestoreWorkArray - Restores a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray()
993 
994   Not Collective
995 
996   Input Parameters:
997 + dm - the DM object
998 . count - The minium size
999 - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT)
1000 
1001   Output Parameter:
1002 . array - the work array
1003 
1004   Level: developer
1005 
1006 .seealso DMDestroy(), DMCreate()
1007 @*/
1008 PetscErrorCode DMRestoreWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem)
1009 {
1010   DMWorkLink *p,link;
1011 
1012   PetscFunctionBegin;
1013   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1014   PetscValidPointer(mem,4);
1015   for (p=&dm->workout; (link=*p); p=&link->next) {
1016     if (link->mem == *(void**)mem) {
1017       *p           = link->next;
1018       link->next   = dm->workin;
1019       dm->workin   = link;
1020       *(void**)mem = NULL;
1021       PetscFunctionReturn(0);
1022     }
1023   }
1024   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out");
1025   PetscFunctionReturn(0);
1026 }
1027 
1028 #undef __FUNCT__
1029 #define __FUNCT__ "DMSetNullSpaceConstructor"
1030 PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt field, MatNullSpace *nullSpace))
1031 {
1032   PetscFunctionBegin;
1033   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1034   if (field >= 10) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %d >= 10 fields", field);
1035   dm->nullspaceConstructors[field] = nullsp;
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 #undef __FUNCT__
1040 #define __FUNCT__ "DMCreateFieldIS"
1041 /*@C
1042   DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field
1043 
1044   Not collective
1045 
1046   Input Parameter:
1047 . dm - the DM object
1048 
1049   Output Parameters:
1050 + numFields  - The number of fields (or NULL if not requested)
1051 . fieldNames - The name for each field (or NULL if not requested)
1052 - fields     - The global indices for each field (or NULL if not requested)
1053 
1054   Level: intermediate
1055 
1056   Notes:
1057   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
1058   PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with
1059   PetscFree().
1060 
1061 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
1062 @*/
1063 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1064 {
1065   PetscSection   section, sectionGlobal;
1066   PetscErrorCode ierr;
1067 
1068   PetscFunctionBegin;
1069   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1070   if (numFields) {
1071     PetscValidPointer(numFields,2);
1072     *numFields = 0;
1073   }
1074   if (fieldNames) {
1075     PetscValidPointer(fieldNames,3);
1076     *fieldNames = NULL;
1077   }
1078   if (fields) {
1079     PetscValidPointer(fields,4);
1080     *fields = NULL;
1081   }
1082   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1083   if (section) {
1084     PetscInt *fieldSizes, **fieldIndices;
1085     PetscInt nF, f, pStart, pEnd, p;
1086 
1087     ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
1088     ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
1089     ierr = PetscMalloc2(nF,PetscInt,&fieldSizes,nF,PetscInt*,&fieldIndices);CHKERRQ(ierr);
1090     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
1091     for (f = 0; f < nF; ++f) {
1092       fieldSizes[f] = 0;
1093     }
1094     for (p = pStart; p < pEnd; ++p) {
1095       PetscInt gdof;
1096 
1097       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1098       if (gdof > 0) {
1099         for (f = 0; f < nF; ++f) {
1100           PetscInt fdof, fcdof;
1101 
1102           ierr           = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
1103           ierr           = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
1104           fieldSizes[f] += fdof-fcdof;
1105         }
1106       }
1107     }
1108     for (f = 0; f < nF; ++f) {
1109       ierr          = PetscMalloc(fieldSizes[f] * sizeof(PetscInt), &fieldIndices[f]);CHKERRQ(ierr);
1110       fieldSizes[f] = 0;
1111     }
1112     for (p = pStart; p < pEnd; ++p) {
1113       PetscInt gdof, goff;
1114 
1115       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1116       if (gdof > 0) {
1117         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
1118         for (f = 0; f < nF; ++f) {
1119           PetscInt fdof, fcdof, fc;
1120 
1121           ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
1122           ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
1123           for (fc = 0; fc < fdof-fcdof; ++fc, ++fieldSizes[f]) {
1124             fieldIndices[f][fieldSizes[f]] = goff++;
1125           }
1126         }
1127       }
1128     }
1129     if (numFields) *numFields = nF;
1130     if (fieldNames) {
1131       ierr = PetscMalloc(nF * sizeof(char*), fieldNames);CHKERRQ(ierr);
1132       for (f = 0; f < nF; ++f) {
1133         const char *fieldName;
1134 
1135         ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr);
1136         ierr = PetscStrallocpy(fieldName, (char**) &(*fieldNames)[f]);CHKERRQ(ierr);
1137       }
1138     }
1139     if (fields) {
1140       ierr = PetscMalloc(nF * sizeof(IS), fields);CHKERRQ(ierr);
1141       for (f = 0; f < nF; ++f) {
1142         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]);CHKERRQ(ierr);
1143       }
1144     }
1145     ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr);
1146   } else if (dm->ops->createfieldis) {
1147     ierr = (*dm->ops->createfieldis)(dm, numFields, fieldNames, fields);CHKERRQ(ierr);
1148   }
1149   PetscFunctionReturn(0);
1150 }
1151 
1152 
1153 #undef __FUNCT__
1154 #define __FUNCT__ "DMCreateFieldDecomposition"
1155 /*@C
1156   DMCreateFieldDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems
1157                           corresponding to different fields: each IS contains the global indices of the dofs of the
1158                           corresponding field. The optional list of DMs define the DM for each subproblem.
1159                           Generalizes DMCreateFieldIS().
1160 
1161   Not collective
1162 
1163   Input Parameter:
1164 . dm - the DM object
1165 
1166   Output Parameters:
1167 + len       - The number of subproblems in the field decomposition (or NULL if not requested)
1168 . namelist  - The name for each field (or NULL if not requested)
1169 . islist    - The global indices for each field (or NULL if not requested)
1170 - dmlist    - The DMs for each field subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined)
1171 
1172   Level: intermediate
1173 
1174   Notes:
1175   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
1176   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
1177   and all of the arrays should be freed with PetscFree().
1178 
1179 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1180 @*/
1181 PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1182 {
1183   PetscErrorCode ierr;
1184 
1185   PetscFunctionBegin;
1186   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1187   if (len) {
1188     PetscValidPointer(len,2);
1189     *len = 0;
1190   }
1191   if (namelist) {
1192     PetscValidPointer(namelist,3);
1193     *namelist = 0;
1194   }
1195   if (islist) {
1196     PetscValidPointer(islist,4);
1197     *islist = 0;
1198   }
1199   if (dmlist) {
1200     PetscValidPointer(dmlist,5);
1201     *dmlist = 0;
1202   }
1203   /*
1204    Is it a good idea to apply the following check across all impls?
1205    Perhaps some impls can have a well-defined decomposition before DMSetUp?
1206    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
1207    */
1208   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
1209   if (!dm->ops->createfielddecomposition) {
1210     PetscSection section;
1211     PetscInt     numFields, f;
1212 
1213     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1214     if (section) {ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);}
1215     if (section && numFields && dm->ops->createsubdm) {
1216       *len = numFields;
1217       ierr = PetscMalloc3(numFields,char*,namelist,numFields,IS,islist,numFields,DM,dmlist);CHKERRQ(ierr);
1218       for (f = 0; f < numFields; ++f) {
1219         const char *fieldName;
1220 
1221         ierr = DMCreateSubDM(dm, 1, &f, &(*islist)[f], &(*dmlist)[f]);CHKERRQ(ierr);
1222         ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr);
1223         ierr = PetscStrallocpy(fieldName, (char**) &(*namelist)[f]);CHKERRQ(ierr);
1224       }
1225     } else {
1226       ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr);
1227       /* By default there are no DMs associated with subproblems. */
1228       if (dmlist) *dmlist = NULL;
1229     }
1230   } else {
1231     ierr = (*dm->ops->createfielddecomposition)(dm,len,namelist,islist,dmlist);CHKERRQ(ierr);
1232   }
1233   PetscFunctionReturn(0);
1234 }
1235 
1236 #undef __FUNCT__
1237 #define __FUNCT__ "DMCreateSubDM"
1238 /*@C
1239   DMCreateSubDM - Returns an IS and DM encapsulating a subproblem defined by the fields passed in.
1240                   The fields are defined by DMCreateFieldIS().
1241 
1242   Not collective
1243 
1244   Input Parameters:
1245 + dm - the DM object
1246 . numFields - number of fields in this subproblem
1247 - len       - The number of subproblems in the decomposition (or NULL if not requested)
1248 
1249   Output Parameters:
1250 . is - The global indices for the subproblem
1251 - dm - The DM for the subproblem
1252 
1253   Level: intermediate
1254 
1255 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1256 @*/
1257 PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1258 {
1259   PetscErrorCode ierr;
1260 
1261   PetscFunctionBegin;
1262   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1263   PetscValidPointer(fields,3);
1264   if (is) PetscValidPointer(is,4);
1265   if (subdm) PetscValidPointer(subdm,5);
1266   if (dm->ops->createsubdm) {
1267     ierr = (*dm->ops->createsubdm)(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1268   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateSubDM implementation defined");
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 
1273 #undef __FUNCT__
1274 #define __FUNCT__ "DMCreateDomainDecomposition"
1275 /*@C
1276   DMCreateDomainDecomposition - Returns lists of IS objects defining a decomposition of a problem into subproblems
1277                           corresponding to restrictions to pairs nested subdomains: each IS contains the global
1278                           indices of the dofs of the corresponding subdomains.  The inner subdomains conceptually
1279                           define a nonoverlapping covering, while outer subdomains can overlap.
1280                           The optional list of DMs define the DM for each subproblem.
1281 
1282   Not collective
1283 
1284   Input Parameter:
1285 . dm - the DM object
1286 
1287   Output Parameters:
1288 + len         - The number of subproblems in the domain decomposition (or NULL if not requested)
1289 . namelist    - The name for each subdomain (or NULL if not requested)
1290 . innerislist - The global indices for each inner subdomain (or NULL, if not requested)
1291 . outerislist - The global indices for each outer subdomain (or NULL, if not requested)
1292 - dmlist      - The DMs for each subdomain subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined)
1293 
1294   Level: intermediate
1295 
1296   Notes:
1297   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
1298   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
1299   and all of the arrays should be freed with PetscFree().
1300 
1301 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateDomainDecompositionDM(), DMCreateFieldDecomposition()
1302 @*/
1303 PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist)
1304 {
1305   PetscErrorCode      ierr;
1306   DMSubDomainHookLink link;
1307   PetscInt            i,l;
1308 
1309   PetscFunctionBegin;
1310   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1311   if (len)           {PetscValidPointer(len,2);            *len         = 0;}
1312   if (namelist)      {PetscValidPointer(namelist,3);       *namelist    = NULL;}
1313   if (innerislist)   {PetscValidPointer(innerislist,4);    *innerislist = NULL;}
1314   if (outerislist)   {PetscValidPointer(outerislist,5);    *outerislist = NULL;}
1315   if (dmlist)        {PetscValidPointer(dmlist,6);         *dmlist      = NULL;}
1316   /*
1317    Is it a good idea to apply the following check across all impls?
1318    Perhaps some impls can have a well-defined decomposition before DMSetUp?
1319    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
1320    */
1321   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
1322   if (dm->ops->createdomaindecomposition) {
1323     ierr = (*dm->ops->createdomaindecomposition)(dm,&l,namelist,innerislist,outerislist,dmlist);CHKERRQ(ierr);
1324     /* copy subdomain hooks and context over to the subdomain DMs */
1325     if (dmlist) {
1326       for (i = 0; i < l; i++) {
1327         for (link=dm->subdomainhook; link; link=link->next) {
1328           if (link->ddhook) {ierr = (*link->ddhook)(dm,(*dmlist)[i],link->ctx);CHKERRQ(ierr);}
1329         }
1330         (*dmlist)[i]->ctx = dm->ctx;
1331       }
1332     }
1333     if (len) *len = l;
1334   }
1335   PetscFunctionReturn(0);
1336 }
1337 
1338 
1339 #undef __FUNCT__
1340 #define __FUNCT__ "DMCreateDomainDecompositionScatters"
1341 /*@C
1342   DMCreateDomainDecompositionScatters - Returns scatters to the subdomain vectors from the global vector
1343 
1344   Not collective
1345 
1346   Input Parameters:
1347 + dm - the DM object
1348 . n  - the number of subdomain scatters
1349 - subdms - the local subdomains
1350 
1351   Output Parameters:
1352 + n     - the number of scatters returned
1353 . iscat - scatter from global vector to nonoverlapping global vector entries on subdomain
1354 . oscat - scatter from global vector to overlapping global vector entries on subdomain
1355 - gscat - scatter from global vector to local vector on subdomain (fills in ghosts)
1356 
1357   Notes: This is an alternative to the iis and ois arguments in DMCreateDomainDecomposition that allow for the solution
1358   of general nonlinear problems with overlapping subdomain methods.  While merely having index sets that enable subsets
1359   of the residual equations to be created is fine for linear problems, nonlinear problems require local assembly of
1360   solution and residual data.
1361 
1362   Level: developer
1363 
1364 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
1365 @*/
1366 PetscErrorCode DMCreateDomainDecompositionScatters(DM dm,PetscInt n,DM *subdms,VecScatter **iscat,VecScatter **oscat,VecScatter **gscat)
1367 {
1368   PetscErrorCode ierr;
1369 
1370   PetscFunctionBegin;
1371   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1372   PetscValidPointer(subdms,3);
1373   if (dm->ops->createddscatters) {
1374     ierr = (*dm->ops->createddscatters)(dm,n,subdms,iscat,oscat,gscat);CHKERRQ(ierr);
1375   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateDomainDecompositionLocalScatter implementation defined");
1376   PetscFunctionReturn(0);
1377 }
1378 
1379 #undef __FUNCT__
1380 #define __FUNCT__ "DMRefine"
1381 /*@
1382   DMRefine - Refines a DM object
1383 
1384   Collective on DM
1385 
1386   Input Parameter:
1387 + dm   - the DM object
1388 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
1389 
1390   Output Parameter:
1391 . dmf - the refined DM, or NULL
1392 
1393   Note: If no refinement was done, the return value is NULL
1394 
1395   Level: developer
1396 
1397 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1398 @*/
1399 PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
1400 {
1401   PetscErrorCode   ierr;
1402   DMRefineHookLink link;
1403 
1404   PetscFunctionBegin;
1405   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1406   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
1407   if (*dmf) {
1408     (*dmf)->ops->creatematrix = dm->ops->creatematrix;
1409 
1410     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
1411 
1412     (*dmf)->ctx       = dm->ctx;
1413     (*dmf)->leveldown = dm->leveldown;
1414     (*dmf)->levelup   = dm->levelup + 1;
1415 
1416     ierr = DMSetMatType(*dmf,dm->mattype);CHKERRQ(ierr);
1417     for (link=dm->refinehook; link; link=link->next) {
1418       if (link->refinehook) {
1419         ierr = (*link->refinehook)(dm,*dmf,link->ctx);CHKERRQ(ierr);
1420       }
1421     }
1422   }
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 #undef __FUNCT__
1427 #define __FUNCT__ "DMRefineHookAdd"
1428 /*@C
1429    DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid
1430 
1431    Logically Collective
1432 
1433    Input Arguments:
1434 +  coarse - nonlinear solver context on which to run a hook when restricting to a coarser level
1435 .  refinehook - function to run when setting up a coarser level
1436 .  interphook - function to run to update data on finer levels (once per SNESSolve())
1437 -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1438 
1439    Calling sequence of refinehook:
1440 $    refinehook(DM coarse,DM fine,void *ctx);
1441 
1442 +  coarse - coarse level DM
1443 .  fine - fine level DM to interpolate problem to
1444 -  ctx - optional user-defined function context
1445 
1446    Calling sequence for interphook:
1447 $    interphook(DM coarse,Mat interp,DM fine,void *ctx)
1448 
1449 +  coarse - coarse level DM
1450 .  interp - matrix interpolating a coarse-level solution to the finer grid
1451 .  fine - fine level DM to update
1452 -  ctx - optional user-defined function context
1453 
1454    Level: advanced
1455 
1456    Notes:
1457    This function is only needed if auxiliary data needs to be passed to fine grids while grid sequencing
1458 
1459    If this function is called multiple times, the hooks will be run in the order they are added.
1460 
1461    This function is currently not available from Fortran.
1462 
1463 .seealso: DMCoarsenHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1464 @*/
1465 PetscErrorCode DMRefineHookAdd(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx)
1466 {
1467   PetscErrorCode   ierr;
1468   DMRefineHookLink link,*p;
1469 
1470   PetscFunctionBegin;
1471   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
1472   for (p=&coarse->refinehook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1473   ierr             = PetscMalloc(sizeof(struct _DMRefineHookLink),&link);CHKERRQ(ierr);
1474   link->refinehook = refinehook;
1475   link->interphook = interphook;
1476   link->ctx        = ctx;
1477   link->next       = NULL;
1478   *p               = link;
1479   PetscFunctionReturn(0);
1480 }
1481 
1482 #undef __FUNCT__
1483 #define __FUNCT__ "DMInterpolate"
1484 /*@
1485    DMInterpolate - interpolates user-defined problem data to a finer DM by running hooks registered by DMRefineHookAdd()
1486 
1487    Collective if any hooks are
1488 
1489    Input Arguments:
1490 +  coarse - coarser DM to use as a base
1491 .  restrct - interpolation matrix, apply using MatInterpolate()
1492 -  fine - finer DM to update
1493 
1494    Level: developer
1495 
1496 .seealso: DMRefineHookAdd(), MatInterpolate()
1497 @*/
1498 PetscErrorCode DMInterpolate(DM coarse,Mat interp,DM fine)
1499 {
1500   PetscErrorCode   ierr;
1501   DMRefineHookLink link;
1502 
1503   PetscFunctionBegin;
1504   for (link=fine->refinehook; link; link=link->next) {
1505     if (link->interphook) {
1506       ierr = (*link->interphook)(coarse,interp,fine,link->ctx);CHKERRQ(ierr);
1507     }
1508   }
1509   PetscFunctionReturn(0);
1510 }
1511 
1512 #undef __FUNCT__
1513 #define __FUNCT__ "DMGetRefineLevel"
1514 /*@
1515     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
1516 
1517     Not Collective
1518 
1519     Input Parameter:
1520 .   dm - the DM object
1521 
1522     Output Parameter:
1523 .   level - number of refinements
1524 
1525     Level: developer
1526 
1527 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1528 
1529 @*/
1530 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
1531 {
1532   PetscFunctionBegin;
1533   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1534   *level = dm->levelup;
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 #undef __FUNCT__
1539 #define __FUNCT__ "DMGlobalToLocalHookAdd"
1540 /*@C
1541    DMGlobalToLocalHookAdd - adds a callback to be run when global to local is called
1542 
1543    Logically Collective
1544 
1545    Input Arguments:
1546 +  dm - the DM
1547 .  beginhook - function to run at the beginning of DMGlobalToLocalBegin()
1548 .  endhook - function to run after DMGlobalToLocalEnd() has completed
1549 -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1550 
1551    Calling sequence for beginhook:
1552 $    beginhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx)
1553 
1554 +  dm - global DM
1555 .  g - global vector
1556 .  mode - mode
1557 .  l - local vector
1558 -  ctx - optional user-defined function context
1559 
1560 
1561    Calling sequence for endhook:
1562 $    endhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx)
1563 
1564 +  global - global DM
1565 -  ctx - optional user-defined function context
1566 
1567    Level: advanced
1568 
1569 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1570 @*/
1571 PetscErrorCode DMGlobalToLocalHookAdd(DM dm,PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*),void *ctx)
1572 {
1573   PetscErrorCode          ierr;
1574   DMGlobalToLocalHookLink link,*p;
1575 
1576   PetscFunctionBegin;
1577   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1578   for (p=&dm->gtolhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1579   ierr            = PetscMalloc(sizeof(struct _DMGlobalToLocalHookLink),&link);CHKERRQ(ierr);
1580   link->beginhook = beginhook;
1581   link->endhook   = endhook;
1582   link->ctx       = ctx;
1583   link->next      = NULL;
1584   *p              = link;
1585   PetscFunctionReturn(0);
1586 }
1587 
1588 #undef __FUNCT__
1589 #define __FUNCT__ "DMGlobalToLocalBegin"
1590 /*@
1591     DMGlobalToLocalBegin - Begins updating local vectors from global vector
1592 
1593     Neighbor-wise Collective on DM
1594 
1595     Input Parameters:
1596 +   dm - the DM object
1597 .   g - the global vector
1598 .   mode - INSERT_VALUES or ADD_VALUES
1599 -   l - the local vector
1600 
1601 
1602     Level: beginner
1603 
1604 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1605 
1606 @*/
1607 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
1608 {
1609   PetscSF                 sf;
1610   PetscErrorCode          ierr;
1611   DMGlobalToLocalHookLink link;
1612 
1613   PetscFunctionBegin;
1614   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1615   for (link=dm->gtolhook; link; link=link->next) {
1616     if (link->beginhook) {
1617       ierr = (*link->beginhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);
1618     }
1619   }
1620   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
1621   if (sf) {
1622     PetscScalar *lArray, *gArray;
1623 
1624     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
1625     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
1626     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
1627     ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
1628     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
1629     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
1630   } else {
1631     ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1632   }
1633   PetscFunctionReturn(0);
1634 }
1635 
1636 #undef __FUNCT__
1637 #define __FUNCT__ "DMGlobalToLocalEnd"
1638 /*@
1639     DMGlobalToLocalEnd - Ends updating local vectors from global vector
1640 
1641     Neighbor-wise Collective on DM
1642 
1643     Input Parameters:
1644 +   dm - the DM object
1645 .   g - the global vector
1646 .   mode - INSERT_VALUES or ADD_VALUES
1647 -   l - the local vector
1648 
1649 
1650     Level: beginner
1651 
1652 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1653 
1654 @*/
1655 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
1656 {
1657   PetscSF                 sf;
1658   PetscErrorCode          ierr;
1659   PetscScalar             *lArray, *gArray;
1660   DMGlobalToLocalHookLink link;
1661 
1662   PetscFunctionBegin;
1663   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1664   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
1665   if (sf) {
1666     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
1667 
1668     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
1669     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
1670     ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
1671     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
1672     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
1673   } else {
1674     ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1675   }
1676   for (link=dm->gtolhook; link; link=link->next) {
1677     if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);}
1678   }
1679   PetscFunctionReturn(0);
1680 }
1681 
1682 #undef __FUNCT__
1683 #define __FUNCT__ "DMLocalToGlobalBegin"
1684 /*@
1685     DMLocalToGlobalBegin - updates global vectors from local vectors
1686 
1687     Neighbor-wise Collective on DM
1688 
1689     Input Parameters:
1690 +   dm - the DM object
1691 .   l - the local vector
1692 .   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
1693            base point.
1694 - - the global vector
1695 
1696     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
1697            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
1698            global array to the final global array with VecAXPY().
1699 
1700     Level: beginner
1701 
1702 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
1703 
1704 @*/
1705 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
1706 {
1707   PetscSF        sf;
1708   PetscErrorCode ierr;
1709 
1710   PetscFunctionBegin;
1711   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1712   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
1713   if (sf) {
1714     MPI_Op      op;
1715     PetscScalar *lArray, *gArray;
1716 
1717     switch (mode) {
1718     case INSERT_VALUES:
1719     case INSERT_ALL_VALUES:
1720       op = MPIU_REPLACE; break;
1721     case ADD_VALUES:
1722     case ADD_ALL_VALUES:
1723       op = MPI_SUM; break;
1724     default:
1725       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
1726     }
1727     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
1728     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
1729     ierr = PetscSFReduceBegin(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr);
1730     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
1731     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
1732   } else {
1733     ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
1734   }
1735   PetscFunctionReturn(0);
1736 }
1737 
1738 #undef __FUNCT__
1739 #define __FUNCT__ "DMLocalToGlobalEnd"
1740 /*@
1741     DMLocalToGlobalEnd - updates global vectors from local vectors
1742 
1743     Neighbor-wise Collective on DM
1744 
1745     Input Parameters:
1746 +   dm - the DM object
1747 .   l - the local vector
1748 .   mode - INSERT_VALUES or ADD_VALUES
1749 -   g - the global vector
1750 
1751 
1752     Level: beginner
1753 
1754 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
1755 
1756 @*/
1757 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
1758 {
1759   PetscSF        sf;
1760   PetscErrorCode ierr;
1761 
1762   PetscFunctionBegin;
1763   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1764   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
1765   if (sf) {
1766     MPI_Op      op;
1767     PetscScalar *lArray, *gArray;
1768 
1769     switch (mode) {
1770     case INSERT_VALUES:
1771     case INSERT_ALL_VALUES:
1772       op = MPIU_REPLACE; break;
1773     case ADD_VALUES:
1774     case ADD_ALL_VALUES:
1775       op = MPI_SUM; break;
1776     default:
1777       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
1778     }
1779     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
1780     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
1781     ierr = PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr);
1782     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
1783     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
1784   } else {
1785     ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
1786   }
1787   PetscFunctionReturn(0);
1788 }
1789 
1790 #undef __FUNCT__
1791 #define __FUNCT__ "DMLocalToLocalBegin"
1792 /* TODO: Fix this documentation!  Still a copy of DMGlobalToLocalEnd. */
1793 /*@
1794     DMLocalToLocalBegin - Begins updating local vectors from global vector
1795 
1796     Neighbor-wise Collective on DM
1797 
1798     Input Parameters:
1799 +   dm - the DM object
1800 .   g - the global vector
1801 .   mode - INSERT_VALUES or ADD_VALUES
1802 -   l - the local vector
1803 
1804 
1805     Level: intermediate
1806 
1807 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1808 
1809 @*/
1810 PetscErrorCode  DMLocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
1811 {
1812   PetscSF                 sf;
1813   PetscErrorCode          ierr;
1814 
1815   PetscFunctionBegin;
1816   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1817   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
1818   if (sf) {
1819     PetscScalar *lArray, *gArray;
1820 
1821     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
1822     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
1823     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
1824     ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
1825     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
1826     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
1827   } else {
1828     ierr = (*dm->ops->localtolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1829   }
1830   PetscFunctionReturn(0);
1831 }
1832 
1833 #undef __FUNCT__
1834 #define __FUNCT__ "DMLocalToLocalEnd"
1835 /* TODO: Fix this documentation!  Still a copy of DMGlobalToLocalEnd. */
1836 /*@
1837     DMLocalToLocalEnd - Ends updating local vectors from global vector
1838 
1839     Neighbor-wise Collective on DM
1840 
1841     Input Parameters:
1842 +   dm - the DM object
1843 .   g - the global vector
1844 .   mode - INSERT_VALUES or ADD_VALUES
1845 -   l - the local vector
1846 
1847 
1848     Level: intermediate
1849 
1850 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1851 
1852 @*/
1853 PetscErrorCode  DMLocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
1854 {
1855   PetscSF                 sf;
1856   PetscErrorCode          ierr;
1857   PetscScalar             *lArray, *gArray;
1858   DMGlobalToLocalHookLink link;
1859 
1860   PetscFunctionBegin;
1861   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1862   ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr);
1863   if (sf) {
1864     if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode);
1865 
1866     ierr = VecGetArray(l, &lArray);CHKERRQ(ierr);
1867     ierr = VecGetArray(g, &gArray);CHKERRQ(ierr);
1868     ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr);
1869     ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr);
1870     ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr);
1871   } else {
1872     ierr = (*dm->ops->localtolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1873   }
1874   for (link=dm->gtolhook; link; link=link->next) {
1875     if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);}
1876   }
1877   PetscFunctionReturn(0);
1878 }
1879 
1880 
1881 #undef __FUNCT__
1882 #define __FUNCT__ "DMCoarsen"
1883 /*@
1884     DMCoarsen - Coarsens a DM object
1885 
1886     Collective on DM
1887 
1888     Input Parameter:
1889 +   dm - the DM object
1890 -   comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
1891 
1892     Output Parameter:
1893 .   dmc - the coarsened DM
1894 
1895     Level: developer
1896 
1897 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1898 
1899 @*/
1900 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
1901 {
1902   PetscErrorCode    ierr;
1903   DMCoarsenHookLink link;
1904 
1905   PetscFunctionBegin;
1906   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1907   ierr                      = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
1908   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
1909   ierr                      = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1910   (*dmc)->ctx               = dm->ctx;
1911   (*dmc)->levelup           = dm->levelup;
1912   (*dmc)->leveldown         = dm->leveldown + 1;
1913   ierr                      = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr);
1914   for (link=dm->coarsenhook; link; link=link->next) {
1915     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1916   }
1917   PetscFunctionReturn(0);
1918 }
1919 
1920 #undef __FUNCT__
1921 #define __FUNCT__ "DMCoarsenHookAdd"
1922 /*@C
1923    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1924 
1925    Logically Collective
1926 
1927    Input Arguments:
1928 +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1929 .  coarsenhook - function to run when setting up a coarser level
1930 .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
1931 -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1932 
1933    Calling sequence of coarsenhook:
1934 $    coarsenhook(DM fine,DM coarse,void *ctx);
1935 
1936 +  fine - fine level DM
1937 .  coarse - coarse level DM to restrict problem to
1938 -  ctx - optional user-defined function context
1939 
1940    Calling sequence for restricthook:
1941 $    restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx)
1942 
1943 +  fine - fine level DM
1944 .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1945 .  rscale - scaling vector for restriction
1946 .  inject - matrix restricting by injection
1947 .  coarse - coarse level DM to update
1948 -  ctx - optional user-defined function context
1949 
1950    Level: advanced
1951 
1952    Notes:
1953    This function is only needed if auxiliary data needs to be set up on coarse grids.
1954 
1955    If this function is called multiple times, the hooks will be run in the order they are added.
1956 
1957    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1958    extract the finest level information from its context (instead of from the SNES).
1959 
1960    This function is currently not available from Fortran.
1961 
1962 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1963 @*/
1964 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1965 {
1966   PetscErrorCode    ierr;
1967   DMCoarsenHookLink link,*p;
1968 
1969   PetscFunctionBegin;
1970   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
1971   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1972   ierr               = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1973   link->coarsenhook  = coarsenhook;
1974   link->restricthook = restricthook;
1975   link->ctx          = ctx;
1976   link->next         = NULL;
1977   *p                 = link;
1978   PetscFunctionReturn(0);
1979 }
1980 
1981 #undef __FUNCT__
1982 #define __FUNCT__ "DMRestrict"
1983 /*@
1984    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1985 
1986    Collective if any hooks are
1987 
1988    Input Arguments:
1989 +  fine - finer DM to use as a base
1990 .  restrct - restriction matrix, apply using MatRestrict()
1991 .  inject - injection matrix, also use MatRestrict()
1992 -  coarse - coarer DM to update
1993 
1994    Level: developer
1995 
1996 .seealso: DMCoarsenHookAdd(), MatRestrict()
1997 @*/
1998 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1999 {
2000   PetscErrorCode    ierr;
2001   DMCoarsenHookLink link;
2002 
2003   PetscFunctionBegin;
2004   for (link=fine->coarsenhook; link; link=link->next) {
2005     if (link->restricthook) {
2006       ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);
2007     }
2008   }
2009   PetscFunctionReturn(0);
2010 }
2011 
2012 #undef __FUNCT__
2013 #define __FUNCT__ "DMSubDomainHookAdd"
2014 /*@C
2015    DMSubDomainHookAdd - adds a callback to be run when restricting a problem to the coarse grid
2016 
2017    Logically Collective
2018 
2019    Input Arguments:
2020 +  global - global DM
2021 .  ddhook - function to run to pass data to the decomposition DM upon its creation
2022 .  restricthook - function to run to update data on block solve (at the beginning of the block solve)
2023 -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
2024 
2025 
2026    Calling sequence for ddhook:
2027 $    ddhook(DM global,DM block,void *ctx)
2028 
2029 +  global - global DM
2030 .  block  - block DM
2031 -  ctx - optional user-defined function context
2032 
2033    Calling sequence for restricthook:
2034 $    restricthook(DM global,VecScatter out,VecScatter in,DM block,void *ctx)
2035 
2036 +  global - global DM
2037 .  out    - scatter to the outer (with ghost and overlap points) block vector
2038 .  in     - scatter to block vector values only owned locally
2039 .  block  - block DM
2040 -  ctx - optional user-defined function context
2041 
2042    Level: advanced
2043 
2044    Notes:
2045    This function is only needed if auxiliary data needs to be set up on subdomain DMs.
2046 
2047    If this function is called multiple times, the hooks will be run in the order they are added.
2048 
2049    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
2050    extract the global information from its context (instead of from the SNES).
2051 
2052    This function is currently not available from Fortran.
2053 
2054 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
2055 @*/
2056 PetscErrorCode DMSubDomainHookAdd(DM global,PetscErrorCode (*ddhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*),void *ctx)
2057 {
2058   PetscErrorCode      ierr;
2059   DMSubDomainHookLink link,*p;
2060 
2061   PetscFunctionBegin;
2062   PetscValidHeaderSpecific(global,DM_CLASSID,1);
2063   for (p=&global->subdomainhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
2064   ierr               = PetscMalloc(sizeof(struct _DMSubDomainHookLink),&link);CHKERRQ(ierr);
2065   link->restricthook = restricthook;
2066   link->ddhook       = ddhook;
2067   link->ctx          = ctx;
2068   link->next         = NULL;
2069   *p                 = link;
2070   PetscFunctionReturn(0);
2071 }
2072 
2073 #undef __FUNCT__
2074 #define __FUNCT__ "DMSubDomainRestrict"
2075 /*@
2076    DMSubDomainRestrict - restricts user-defined problem data to a block DM by running hooks registered by DMSubDomainHookAdd()
2077 
2078    Collective if any hooks are
2079 
2080    Input Arguments:
2081 +  fine - finer DM to use as a base
2082 .  oscatter - scatter from domain global vector filling subdomain global vector with overlap
2083 .  gscatter - scatter from domain global vector filling subdomain local vector with ghosts
2084 -  coarse - coarer DM to update
2085 
2086    Level: developer
2087 
2088 .seealso: DMCoarsenHookAdd(), MatRestrict()
2089 @*/
2090 PetscErrorCode DMSubDomainRestrict(DM global,VecScatter oscatter,VecScatter gscatter,DM subdm)
2091 {
2092   PetscErrorCode      ierr;
2093   DMSubDomainHookLink link;
2094 
2095   PetscFunctionBegin;
2096   for (link=global->subdomainhook; link; link=link->next) {
2097     if (link->restricthook) {
2098       ierr = (*link->restricthook)(global,oscatter,gscatter,subdm,link->ctx);CHKERRQ(ierr);
2099     }
2100   }
2101   PetscFunctionReturn(0);
2102 }
2103 
2104 #undef __FUNCT__
2105 #define __FUNCT__ "DMGetCoarsenLevel"
2106 /*@
2107     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
2108 
2109     Not Collective
2110 
2111     Input Parameter:
2112 .   dm - the DM object
2113 
2114     Output Parameter:
2115 .   level - number of coarsenings
2116 
2117     Level: developer
2118 
2119 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
2120 
2121 @*/
2122 PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
2123 {
2124   PetscFunctionBegin;
2125   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2126   *level = dm->leveldown;
2127   PetscFunctionReturn(0);
2128 }
2129 
2130 
2131 
2132 #undef __FUNCT__
2133 #define __FUNCT__ "DMRefineHierarchy"
2134 /*@C
2135     DMRefineHierarchy - Refines a DM object, all levels at once
2136 
2137     Collective on DM
2138 
2139     Input Parameter:
2140 +   dm - the DM object
2141 -   nlevels - the number of levels of refinement
2142 
2143     Output Parameter:
2144 .   dmf - the refined DM hierarchy
2145 
2146     Level: developer
2147 
2148 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
2149 
2150 @*/
2151 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
2152 {
2153   PetscErrorCode ierr;
2154 
2155   PetscFunctionBegin;
2156   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2157   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
2158   if (nlevels == 0) PetscFunctionReturn(0);
2159   if (dm->ops->refinehierarchy) {
2160     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
2161   } else if (dm->ops->refine) {
2162     PetscInt i;
2163 
2164     ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr);
2165     for (i=1; i<nlevels; i++) {
2166       ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr);
2167     }
2168   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
2169   PetscFunctionReturn(0);
2170 }
2171 
2172 #undef __FUNCT__
2173 #define __FUNCT__ "DMCoarsenHierarchy"
2174 /*@C
2175     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
2176 
2177     Collective on DM
2178 
2179     Input Parameter:
2180 +   dm - the DM object
2181 -   nlevels - the number of levels of coarsening
2182 
2183     Output Parameter:
2184 .   dmc - the coarsened DM hierarchy
2185 
2186     Level: developer
2187 
2188 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
2189 
2190 @*/
2191 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
2192 {
2193   PetscErrorCode ierr;
2194 
2195   PetscFunctionBegin;
2196   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2197   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
2198   if (nlevels == 0) PetscFunctionReturn(0);
2199   PetscValidPointer(dmc,3);
2200   if (dm->ops->coarsenhierarchy) {
2201     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
2202   } else if (dm->ops->coarsen) {
2203     PetscInt i;
2204 
2205     ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr);
2206     for (i=1; i<nlevels; i++) {
2207       ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr);
2208     }
2209   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
2210   PetscFunctionReturn(0);
2211 }
2212 
2213 #undef __FUNCT__
2214 #define __FUNCT__ "DMCreateAggregates"
2215 /*@
2216    DMCreateAggregates - Gets the aggregates that map between
2217    grids associated with two DMs.
2218 
2219    Collective on DM
2220 
2221    Input Parameters:
2222 +  dmc - the coarse grid DM
2223 -  dmf - the fine grid DM
2224 
2225    Output Parameters:
2226 .  rest - the restriction matrix (transpose of the projection matrix)
2227 
2228    Level: intermediate
2229 
2230 .keywords: interpolation, restriction, multigrid
2231 
2232 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
2233 @*/
2234 PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
2235 {
2236   PetscErrorCode ierr;
2237 
2238   PetscFunctionBegin;
2239   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
2240   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
2241   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
2242   PetscFunctionReturn(0);
2243 }
2244 
2245 #undef __FUNCT__
2246 #define __FUNCT__ "DMSetApplicationContextDestroy"
2247 /*@C
2248     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
2249 
2250     Not Collective
2251 
2252     Input Parameters:
2253 +   dm - the DM object
2254 -   destroy - the destroy function
2255 
2256     Level: intermediate
2257 
2258 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2259 
2260 @*/
2261 PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
2262 {
2263   PetscFunctionBegin;
2264   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2265   dm->ctxdestroy = destroy;
2266   PetscFunctionReturn(0);
2267 }
2268 
2269 #undef __FUNCT__
2270 #define __FUNCT__ "DMSetApplicationContext"
2271 /*@
2272     DMSetApplicationContext - Set a user context into a DM object
2273 
2274     Not Collective
2275 
2276     Input Parameters:
2277 +   dm - the DM object
2278 -   ctx - the user context
2279 
2280     Level: intermediate
2281 
2282 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2283 
2284 @*/
2285 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
2286 {
2287   PetscFunctionBegin;
2288   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2289   dm->ctx = ctx;
2290   PetscFunctionReturn(0);
2291 }
2292 
2293 #undef __FUNCT__
2294 #define __FUNCT__ "DMGetApplicationContext"
2295 /*@
2296     DMGetApplicationContext - Gets a user context from a DM object
2297 
2298     Not Collective
2299 
2300     Input Parameter:
2301 .   dm - the DM object
2302 
2303     Output Parameter:
2304 .   ctx - the user context
2305 
2306     Level: intermediate
2307 
2308 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2309 
2310 @*/
2311 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
2312 {
2313   PetscFunctionBegin;
2314   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2315   *(void**)ctx = dm->ctx;
2316   PetscFunctionReturn(0);
2317 }
2318 
2319 #undef __FUNCT__
2320 #define __FUNCT__ "DMSetVariableBounds"
2321 /*@C
2322     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
2323 
2324     Logically Collective on DM
2325 
2326     Input Parameter:
2327 +   dm - the DM object
2328 -   f - the function that computes variable bounds used by SNESVI (use NULL to cancel a previous function that was set)
2329 
2330     Level: intermediate
2331 
2332 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(),
2333          DMSetJacobian()
2334 
2335 @*/
2336 PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
2337 {
2338   PetscFunctionBegin;
2339   dm->ops->computevariablebounds = f;
2340   PetscFunctionReturn(0);
2341 }
2342 
2343 #undef __FUNCT__
2344 #define __FUNCT__ "DMHasVariableBounds"
2345 /*@
2346     DMHasVariableBounds - does the DM object have a variable bounds function?
2347 
2348     Not Collective
2349 
2350     Input Parameter:
2351 .   dm - the DM object to destroy
2352 
2353     Output Parameter:
2354 .   flg - PETSC_TRUE if the variable bounds function exists
2355 
2356     Level: developer
2357 
2358 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2359 
2360 @*/
2361 PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
2362 {
2363   PetscFunctionBegin;
2364   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
2365   PetscFunctionReturn(0);
2366 }
2367 
2368 #undef __FUNCT__
2369 #define __FUNCT__ "DMComputeVariableBounds"
2370 /*@C
2371     DMComputeVariableBounds - compute variable bounds used by SNESVI.
2372 
2373     Logically Collective on DM
2374 
2375     Input Parameters:
2376 +   dm - the DM object to destroy
2377 -   x  - current solution at which the bounds are computed
2378 
2379     Output parameters:
2380 +   xl - lower bound
2381 -   xu - upper bound
2382 
2383     Level: intermediate
2384 
2385 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2386 
2387 @*/
2388 PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
2389 {
2390   PetscErrorCode ierr;
2391 
2392   PetscFunctionBegin;
2393   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
2394   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
2395   if (dm->ops->computevariablebounds) {
2396     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu);CHKERRQ(ierr);
2397   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "This DM is incapable of computing variable bounds.");
2398   PetscFunctionReturn(0);
2399 }
2400 
2401 #undef __FUNCT__
2402 #define __FUNCT__ "DMHasColoring"
2403 /*@
2404     DMHasColoring - does the DM object have a method of providing a coloring?
2405 
2406     Not Collective
2407 
2408     Input Parameter:
2409 .   dm - the DM object
2410 
2411     Output Parameter:
2412 .   flg - PETSC_TRUE if the DM has facilities for DMCreateColoring().
2413 
2414     Level: developer
2415 
2416 .seealso DMHasFunction(), DMCreateColoring()
2417 
2418 @*/
2419 PetscErrorCode  DMHasColoring(DM dm,PetscBool  *flg)
2420 {
2421   PetscFunctionBegin;
2422   *flg =  (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE;
2423   PetscFunctionReturn(0);
2424 }
2425 
2426 #undef  __FUNCT__
2427 #define __FUNCT__ "DMSetVec"
2428 /*@C
2429     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
2430 
2431     Collective on DM
2432 
2433     Input Parameter:
2434 +   dm - the DM object
2435 -   x - location to compute residual and Jacobian, if NULL is passed to those routines; will be NULL for linear problems.
2436 
2437     Level: developer
2438 
2439 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2440 
2441 @*/
2442 PetscErrorCode  DMSetVec(DM dm,Vec x)
2443 {
2444   PetscErrorCode ierr;
2445 
2446   PetscFunctionBegin;
2447   if (x) {
2448     if (!dm->x) {
2449       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
2450     }
2451     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
2452   } else if (dm->x) {
2453     ierr = VecDestroy(&dm->x);CHKERRQ(ierr);
2454   }
2455   PetscFunctionReturn(0);
2456 }
2457 
2458 PetscFunctionList DMList              = NULL;
2459 PetscBool         DMRegisterAllCalled = PETSC_FALSE;
2460 
2461 #undef __FUNCT__
2462 #define __FUNCT__ "DMSetType"
2463 /*@C
2464   DMSetType - Builds a DM, for a particular DM implementation.
2465 
2466   Collective on DM
2467 
2468   Input Parameters:
2469 + dm     - The DM object
2470 - method - The name of the DM type
2471 
2472   Options Database Key:
2473 . -dm_type <type> - Sets the DM type; use -help for a list of available types
2474 
2475   Notes:
2476   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
2477 
2478   Level: intermediate
2479 
2480 .keywords: DM, set, type
2481 .seealso: DMGetType(), DMCreate()
2482 @*/
2483 PetscErrorCode  DMSetType(DM dm, DMType method)
2484 {
2485   PetscErrorCode (*r)(DM);
2486   PetscBool      match;
2487   PetscErrorCode ierr;
2488 
2489   PetscFunctionBegin;
2490   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2491   ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
2492   if (match) PetscFunctionReturn(0);
2493 
2494   if (!DMRegisterAllCalled) {ierr = DMRegisterAll();CHKERRQ(ierr);}
2495   ierr = PetscFunctionListFind(DMList,method,&r);CHKERRQ(ierr);
2496   if (!r) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
2497 
2498   if (dm->ops->destroy) {
2499     ierr             = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
2500     dm->ops->destroy = NULL;
2501   }
2502   ierr = (*r)(dm);CHKERRQ(ierr);
2503   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
2504   PetscFunctionReturn(0);
2505 }
2506 
2507 #undef __FUNCT__
2508 #define __FUNCT__ "DMGetType"
2509 /*@C
2510   DMGetType - Gets the DM type name (as a string) from the DM.
2511 
2512   Not Collective
2513 
2514   Input Parameter:
2515 . dm  - The DM
2516 
2517   Output Parameter:
2518 . type - The DM type name
2519 
2520   Level: intermediate
2521 
2522 .keywords: DM, get, type, name
2523 .seealso: DMSetType(), DMCreate()
2524 @*/
2525 PetscErrorCode  DMGetType(DM dm, DMType *type)
2526 {
2527   PetscErrorCode ierr;
2528 
2529   PetscFunctionBegin;
2530   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2531   PetscValidCharPointer(type,2);
2532   if (!DMRegisterAllCalled) {
2533     ierr = DMRegisterAll();CHKERRQ(ierr);
2534   }
2535   *type = ((PetscObject)dm)->type_name;
2536   PetscFunctionReturn(0);
2537 }
2538 
2539 #undef __FUNCT__
2540 #define __FUNCT__ "DMConvert"
2541 /*@C
2542   DMConvert - Converts a DM to another DM, either of the same or different type.
2543 
2544   Collective on DM
2545 
2546   Input Parameters:
2547 + dm - the DM
2548 - newtype - new DM type (use "same" for the same type)
2549 
2550   Output Parameter:
2551 . M - pointer to new DM
2552 
2553   Notes:
2554   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
2555   the MPI communicator of the generated DM is always the same as the communicator
2556   of the input DM.
2557 
2558   Level: intermediate
2559 
2560 .seealso: DMCreate()
2561 @*/
2562 PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M)
2563 {
2564   DM             B;
2565   char           convname[256];
2566   PetscBool      sametype, issame;
2567   PetscErrorCode ierr;
2568 
2569   PetscFunctionBegin;
2570   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2571   PetscValidType(dm,1);
2572   PetscValidPointer(M,3);
2573   ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
2574   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
2575   {
2576     PetscErrorCode (*conv)(DM, DMType, DM*) = NULL;
2577 
2578     /*
2579        Order of precedence:
2580        1) See if a specialized converter is known to the current DM.
2581        2) See if a specialized converter is known to the desired DM class.
2582        3) See if a good general converter is registered for the desired class
2583        4) See if a good general converter is known for the current matrix.
2584        5) Use a really basic converter.
2585     */
2586 
2587     /* 1) See if a specialized converter is known to the current DM and the desired class */
2588     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
2589     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
2590     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
2591     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2592     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2593     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,&conv);CHKERRQ(ierr);
2594     if (conv) goto foundconv;
2595 
2596     /* 2)  See if a specialized converter is known to the desired DM class. */
2597     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &B);CHKERRQ(ierr);
2598     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
2599     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
2600     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
2601     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
2602     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2603     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2604     ierr = PetscObjectQueryFunction((PetscObject)B,convname,&conv);CHKERRQ(ierr);
2605     if (conv) {
2606       ierr = DMDestroy(&B);CHKERRQ(ierr);
2607       goto foundconv;
2608     }
2609 
2610 #if 0
2611     /* 3) See if a good general converter is registered for the desired class */
2612     conv = B->ops->convertfrom;
2613     ierr = DMDestroy(&B);CHKERRQ(ierr);
2614     if (conv) goto foundconv;
2615 
2616     /* 4) See if a good general converter is known for the current matrix */
2617     if (dm->ops->convert) {
2618       conv = dm->ops->convert;
2619     }
2620     if (conv) goto foundconv;
2621 #endif
2622 
2623     /* 5) Use a really basic converter. */
2624     SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
2625 
2626 foundconv:
2627     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2628     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
2629     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2630   }
2631   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
2632   PetscFunctionReturn(0);
2633 }
2634 
2635 /*--------------------------------------------------------------------------------------------------------------------*/
2636 
2637 #undef __FUNCT__
2638 #define __FUNCT__ "DMRegister"
2639 /*@C
2640   DMRegister -  Adds a new DM component implementation
2641 
2642   Not Collective
2643 
2644   Input Parameters:
2645 + name        - The name of a new user-defined creation routine
2646 - create_func - The creation routine itself
2647 
2648   Notes:
2649   DMRegister() may be called multiple times to add several user-defined DMs
2650 
2651 
2652   Sample usage:
2653 .vb
2654     DMRegister("my_da", MyDMCreate);
2655 .ve
2656 
2657   Then, your DM type can be chosen with the procedural interface via
2658 .vb
2659     DMCreate(MPI_Comm, DM *);
2660     DMSetType(DM,"my_da");
2661 .ve
2662    or at runtime via the option
2663 .vb
2664     -da_type my_da
2665 .ve
2666 
2667   Level: advanced
2668 
2669 .keywords: DM, register
2670 .seealso: DMRegisterAll(), DMRegisterDestroy()
2671 
2672 @*/
2673 PetscErrorCode  DMRegister(const char sname[],PetscErrorCode (*function)(DM))
2674 {
2675   PetscErrorCode ierr;
2676 
2677   PetscFunctionBegin;
2678   ierr = PetscFunctionListAdd(&DMList,sname,function);CHKERRQ(ierr);
2679   PetscFunctionReturn(0);
2680 }
2681 
2682 #undef __FUNCT__
2683 #define __FUNCT__ "DMLoad"
2684 /*@C
2685   DMLoad - Loads a DM that has been stored in binary  with DMView().
2686 
2687   Collective on PetscViewer
2688 
2689   Input Parameters:
2690 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2691            some related function before a call to DMLoad().
2692 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2693            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2694 
2695    Level: intermediate
2696 
2697   Notes:
2698    The type is determined by the data in the file, any type set into the DM before this call is ignored.
2699 
2700   Notes for advanced users:
2701   Most users should not need to know the details of the binary storage
2702   format, since DMLoad() and DMView() completely hide these details.
2703   But for anyone who's interested, the standard binary matrix storage
2704   format is
2705 .vb
2706      has not yet been determined
2707 .ve
2708 
2709 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2710 @*/
2711 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2712 {
2713   PetscErrorCode ierr;
2714   PetscBool      isbinary;
2715   PetscInt       classid;
2716   char           type[256];
2717 
2718   PetscFunctionBegin;
2719   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2720   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2721   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
2722   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
2723 
2724   ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
2725   if (classid != DM_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not DM next in file");
2726   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
2727   ierr = DMSetType(newdm, type);CHKERRQ(ierr);
2728   if (newdm->ops->load) {
2729     ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2730   }
2731   PetscFunctionReturn(0);
2732 }
2733 
2734 /******************************** FEM Support **********************************/
2735 
2736 #undef __FUNCT__
2737 #define __FUNCT__ "DMPrintCellVector"
2738 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[])
2739 {
2740   PetscInt       f;
2741   PetscErrorCode ierr;
2742 
2743   PetscFunctionBegin;
2744   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2745   for (f = 0; f < len; ++f) {
2746     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
2747   }
2748   PetscFunctionReturn(0);
2749 }
2750 
2751 #undef __FUNCT__
2752 #define __FUNCT__ "DMPrintCellMatrix"
2753 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
2754 {
2755   PetscInt       f, g;
2756   PetscErrorCode ierr;
2757 
2758   PetscFunctionBegin;
2759   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2760   for (f = 0; f < rows; ++f) {
2761     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
2762     for (g = 0; g < cols; ++g) {
2763       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
2764     }
2765     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
2766   }
2767   PetscFunctionReturn(0);
2768 }
2769 
2770 #undef __FUNCT__
2771 #define __FUNCT__ "DMGetDefaultSection"
2772 /*@
2773   DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM.
2774 
2775   Input Parameter:
2776 . dm - The DM
2777 
2778   Output Parameter:
2779 . section - The PetscSection
2780 
2781   Level: intermediate
2782 
2783   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
2784 
2785 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
2786 @*/
2787 PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section)
2788 {
2789   PetscFunctionBegin;
2790   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2791   PetscValidPointer(section, 2);
2792   *section = dm->defaultSection;
2793   PetscFunctionReturn(0);
2794 }
2795 
2796 #undef __FUNCT__
2797 #define __FUNCT__ "DMSetDefaultSection"
2798 /*@
2799   DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM.
2800 
2801   Input Parameters:
2802 + dm - The DM
2803 - section - The PetscSection
2804 
2805   Level: intermediate
2806 
2807   Note: Any existing Section will be destroyed
2808 
2809 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
2810 @*/
2811 PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section)
2812 {
2813   PetscInt       numFields;
2814   PetscInt       f;
2815   PetscErrorCode ierr;
2816 
2817   PetscFunctionBegin;
2818   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2819   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
2820   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
2821   ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr);
2822   dm->defaultSection = section;
2823   ierr = PetscSectionGetNumFields(dm->defaultSection, &numFields);CHKERRQ(ierr);
2824   if (numFields) {
2825     ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr);
2826     for (f = 0; f < numFields; ++f) {
2827       const char *name;
2828 
2829       ierr = PetscSectionGetFieldName(dm->defaultSection, f, &name);CHKERRQ(ierr);
2830       ierr = PetscObjectSetName(dm->fields[f], name);CHKERRQ(ierr);
2831     }
2832   }
2833   /* The global section will be rebuilt in the next call to DMGetDefaultGlobalSection(). */
2834   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
2835   PetscFunctionReturn(0);
2836 }
2837 
2838 #undef __FUNCT__
2839 #define __FUNCT__ "DMGetDefaultGlobalSection"
2840 /*@
2841   DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM.
2842 
2843   Collective on DM
2844 
2845   Input Parameter:
2846 . dm - The DM
2847 
2848   Output Parameter:
2849 . section - The PetscSection
2850 
2851   Level: intermediate
2852 
2853   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
2854 
2855 .seealso: DMSetDefaultSection(), DMGetDefaultSection()
2856 @*/
2857 PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section)
2858 {
2859   PetscErrorCode ierr;
2860 
2861   PetscFunctionBegin;
2862   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2863   PetscValidPointer(section, 2);
2864   if (!dm->defaultGlobalSection) {
2865     if (!dm->defaultSection || !dm->sf) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSection and PetscSF in order to create a global PetscSection");
2866     ierr = PetscSectionCreateGlobalSection(dm->defaultSection, dm->sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr);
2867     ierr = PetscLayoutDestroy(&dm->map);CHKERRQ(ierr);
2868     ierr = PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm),dm->defaultGlobalSection,&dm->map);CHKERRQ(ierr);
2869   }
2870   *section = dm->defaultGlobalSection;
2871   PetscFunctionReturn(0);
2872 }
2873 
2874 #undef __FUNCT__
2875 #define __FUNCT__ "DMSetDefaultGlobalSection"
2876 /*@
2877   DMSetDefaultGlobalSection - Set the PetscSection encoding the global data layout for the DM.
2878 
2879   Input Parameters:
2880 + dm - The DM
2881 - section - The PetscSection, or NULL
2882 
2883   Level: intermediate
2884 
2885   Note: Any existing Section will be destroyed
2886 
2887 .seealso: DMGetDefaultGlobalSection(), DMSetDefaultSection()
2888 @*/
2889 PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection section)
2890 {
2891   PetscErrorCode ierr;
2892 
2893   PetscFunctionBegin;
2894   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2895   if (section) PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
2896   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
2897   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
2898   dm->defaultGlobalSection = section;
2899   PetscFunctionReturn(0);
2900 }
2901 
2902 #undef __FUNCT__
2903 #define __FUNCT__ "DMGetDefaultSF"
2904 /*@
2905   DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set,
2906   it is created from the default PetscSection layouts in the DM.
2907 
2908   Input Parameter:
2909 . dm - The DM
2910 
2911   Output Parameter:
2912 . sf - The PetscSF
2913 
2914   Level: intermediate
2915 
2916   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
2917 
2918 .seealso: DMSetDefaultSF(), DMCreateDefaultSF()
2919 @*/
2920 PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf)
2921 {
2922   PetscInt       nroots;
2923   PetscErrorCode ierr;
2924 
2925   PetscFunctionBegin;
2926   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2927   PetscValidPointer(sf, 2);
2928   ierr = PetscSFGetGraph(dm->defaultSF, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
2929   if (nroots < 0) {
2930     PetscSection section, gSection;
2931 
2932     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
2933     if (section) {
2934       ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
2935       ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr);
2936     } else {
2937       *sf = NULL;
2938       PetscFunctionReturn(0);
2939     }
2940   }
2941   *sf = dm->defaultSF;
2942   PetscFunctionReturn(0);
2943 }
2944 
2945 #undef __FUNCT__
2946 #define __FUNCT__ "DMSetDefaultSF"
2947 /*@
2948   DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM
2949 
2950   Input Parameters:
2951 + dm - The DM
2952 - sf - The PetscSF
2953 
2954   Level: intermediate
2955 
2956   Note: Any previous SF is destroyed
2957 
2958 .seealso: DMGetDefaultSF(), DMCreateDefaultSF()
2959 @*/
2960 PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf)
2961 {
2962   PetscErrorCode ierr;
2963 
2964   PetscFunctionBegin;
2965   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2966   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
2967   ierr          = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr);
2968   dm->defaultSF = sf;
2969   PetscFunctionReturn(0);
2970 }
2971 
2972 #undef __FUNCT__
2973 #define __FUNCT__ "DMCreateDefaultSF"
2974 /*@C
2975   DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections
2976   describing the data layout.
2977 
2978   Input Parameters:
2979 + dm - The DM
2980 . localSection - PetscSection describing the local data layout
2981 - globalSection - PetscSection describing the global data layout
2982 
2983   Level: intermediate
2984 
2985 .seealso: DMGetDefaultSF(), DMSetDefaultSF()
2986 @*/
2987 PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection)
2988 {
2989   MPI_Comm       comm;
2990   PetscLayout    layout;
2991   const PetscInt *ranges;
2992   PetscInt       *local;
2993   PetscSFNode    *remote;
2994   PetscInt       pStart, pEnd, p, nroots, nleaves, l;
2995   PetscMPIInt    size, rank;
2996   PetscErrorCode ierr;
2997 
2998   PetscFunctionBegin;
2999   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
3000   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3001   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
3002   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
3003   ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr);
3004   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr);
3005   ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
3006   ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
3007   ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr);
3008   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
3009   ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
3010   for (p = pStart, nleaves = 0; p < pEnd; ++p) {
3011     PetscInt gdof, gcdof;
3012 
3013     ierr     = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
3014     ierr     = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
3015     nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
3016   }
3017   ierr = PetscMalloc(nleaves * sizeof(PetscInt), &local);CHKERRQ(ierr);
3018   ierr = PetscMalloc(nleaves * sizeof(PetscSFNode), &remote);CHKERRQ(ierr);
3019   for (p = pStart, l = 0; p < pEnd; ++p) {
3020     const PetscInt *cind;
3021     PetscInt       dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
3022 
3023     ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr);
3024     ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr);
3025     ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr);
3026     ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr);
3027     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
3028     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
3029     ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
3030     if (!gdof) continue; /* Censored point */
3031     gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
3032     if (gsize != dof-cdof) {
3033       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);
3034       cdof = 0; /* Ignore constraints */
3035     }
3036     for (d = 0, c = 0; d < dof; ++d) {
3037       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
3038       local[l+d-c] = off+d;
3039     }
3040     if (gdof < 0) {
3041       for (d = 0; d < gsize; ++d, ++l) {
3042         PetscInt offset = -(goff+1) + d, r;
3043 
3044         ierr = PetscFindInt(offset,size,ranges,&r);CHKERRQ(ierr);
3045         if (r < 0) r = -(r+2);
3046         remote[l].rank  = r;
3047         remote[l].index = offset - ranges[r];
3048       }
3049     } else {
3050       for (d = 0; d < gsize; ++d, ++l) {
3051         remote[l].rank  = rank;
3052         remote[l].index = goff+d - ranges[rank];
3053       }
3054     }
3055   }
3056   if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %d != nleaves %d", l, nleaves);
3057   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
3058   ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
3059   PetscFunctionReturn(0);
3060 }
3061 
3062 #undef __FUNCT__
3063 #define __FUNCT__ "DMGetPointSF"
3064 /*@
3065   DMGetPointSF - Get the PetscSF encoding the parallel section point overlap for the DM.
3066 
3067   Input Parameter:
3068 . dm - The DM
3069 
3070   Output Parameter:
3071 . sf - The PetscSF
3072 
3073   Level: intermediate
3074 
3075   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
3076 
3077 .seealso: DMSetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3078 @*/
3079 PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf)
3080 {
3081   PetscFunctionBegin;
3082   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3083   PetscValidPointer(sf, 2);
3084   *sf = dm->sf;
3085   PetscFunctionReturn(0);
3086 }
3087 
3088 #undef __FUNCT__
3089 #define __FUNCT__ "DMSetPointSF"
3090 /*@
3091   DMSetPointSF - Set the PetscSF encoding the parallel section point overlap for the DM.
3092 
3093   Input Parameters:
3094 + dm - The DM
3095 - sf - The PetscSF
3096 
3097   Level: intermediate
3098 
3099 .seealso: DMGetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3100 @*/
3101 PetscErrorCode DMSetPointSF(DM dm, PetscSF sf)
3102 {
3103   PetscErrorCode ierr;
3104 
3105   PetscFunctionBegin;
3106   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3107   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
3108   ierr   = PetscSFDestroy(&dm->sf);CHKERRQ(ierr);
3109   ierr   = PetscObjectReference((PetscObject) sf);CHKERRQ(ierr);
3110   dm->sf = sf;
3111   PetscFunctionReturn(0);
3112 }
3113 
3114 #undef __FUNCT__
3115 #define __FUNCT__ "DMGetNumFields"
3116 PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields)
3117 {
3118   PetscFunctionBegin;
3119   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3120   PetscValidPointer(numFields, 2);
3121   *numFields = dm->numFields;
3122   PetscFunctionReturn(0);
3123 }
3124 
3125 #undef __FUNCT__
3126 #define __FUNCT__ "DMSetNumFields"
3127 PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields)
3128 {
3129   PetscInt       f;
3130   PetscErrorCode ierr;
3131 
3132   PetscFunctionBegin;
3133   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3134   for (f = 0; f < dm->numFields; ++f) {
3135     ierr = PetscObjectDestroy(&dm->fields[f]);CHKERRQ(ierr);
3136   }
3137   ierr          = PetscFree(dm->fields);CHKERRQ(ierr);
3138   dm->numFields = numFields;
3139   ierr          = PetscMalloc(dm->numFields * sizeof(PetscObject), &dm->fields);CHKERRQ(ierr);
3140   for (f = 0; f < dm->numFields; ++f) {
3141     ierr = PetscContainerCreate(PetscObjectComm((PetscObject)dm), (PetscContainer*) &dm->fields[f]);CHKERRQ(ierr);
3142   }
3143   PetscFunctionReturn(0);
3144 }
3145 
3146 #undef __FUNCT__
3147 #define __FUNCT__ "DMGetField"
3148 PetscErrorCode DMGetField(DM dm, PetscInt f, PetscObject *field)
3149 {
3150   PetscFunctionBegin;
3151   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3152   PetscValidPointer(field, 2);
3153   if (!dm->fields) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Fields have not been setup in this DM. Call DMSetNumFields()");
3154   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);
3155   *field = dm->fields[f];
3156   PetscFunctionReturn(0);
3157 }
3158 
3159 #undef __FUNCT__
3160 #define __FUNCT__ "DMSetCoordinates"
3161 /*@
3162   DMSetCoordinates - Sets into the DM a global vector that holds the coordinates
3163 
3164   Collective on DM
3165 
3166   Input Parameters:
3167 + dm - the DM
3168 - c - coordinate vector
3169 
3170   Note:
3171   The coordinates do include those for ghost points, which are in the local vector
3172 
3173   Level: intermediate
3174 
3175 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3176 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLoca(), DMGetCoordinateDM()
3177 @*/
3178 PetscErrorCode DMSetCoordinates(DM dm, Vec c)
3179 {
3180   PetscErrorCode ierr;
3181 
3182   PetscFunctionBegin;
3183   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3184   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
3185   ierr            = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
3186   ierr            = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
3187   dm->coordinates = c;
3188   ierr            = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
3189   PetscFunctionReturn(0);
3190 }
3191 
3192 #undef __FUNCT__
3193 #define __FUNCT__ "DMSetCoordinatesLocal"
3194 /*@
3195   DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates
3196 
3197   Collective on DM
3198 
3199    Input Parameters:
3200 +  dm - the DM
3201 -  c - coordinate vector
3202 
3203   Note:
3204   The coordinates of ghost points can be set using DMSetCoordinates()
3205   followed by DMGetCoordinatesLocal(). This is intended to enable the
3206   setting of ghost coordinates outside of the domain.
3207 
3208   Level: intermediate
3209 
3210 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3211 .seealso: DMGetCoordinatesLocal(), DMSetCoordinates(), DMGetCoordinates(), DMGetCoordinateDM()
3212 @*/
3213 PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
3214 {
3215   PetscErrorCode ierr;
3216 
3217   PetscFunctionBegin;
3218   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3219   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
3220   ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
3221   ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
3222 
3223   dm->coordinatesLocal = c;
3224 
3225   ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
3226   PetscFunctionReturn(0);
3227 }
3228 
3229 #undef __FUNCT__
3230 #define __FUNCT__ "DMGetCoordinates"
3231 /*@
3232   DMGetCoordinates - Gets a global vector with the coordinates associated with the DM.
3233 
3234   Not Collective
3235 
3236   Input Parameter:
3237 . dm - the DM
3238 
3239   Output Parameter:
3240 . c - global coordinate vector
3241 
3242   Note:
3243   This is a borrowed reference, so the user should NOT destroy this vector
3244 
3245   Each process has only the local coordinates (does NOT have the ghost coordinates).
3246 
3247   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
3248   and (x_0,y_0,z_0,x_1,y_1,z_1...)
3249 
3250   Level: intermediate
3251 
3252 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3253 .seealso: DMSetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM()
3254 @*/
3255 PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
3256 {
3257   PetscErrorCode ierr;
3258 
3259   PetscFunctionBegin;
3260   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3261   PetscValidPointer(c,2);
3262   if (!dm->coordinates && dm->coordinatesLocal) {
3263     DM cdm = NULL;
3264 
3265     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3266     ierr = DMCreateGlobalVector(cdm, &dm->coordinates);CHKERRQ(ierr);
3267     ierr = PetscObjectSetName((PetscObject) dm->coordinates, "coordinates");CHKERRQ(ierr);
3268     ierr = DMLocalToGlobalBegin(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
3269     ierr = DMLocalToGlobalEnd(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
3270   }
3271   *c = dm->coordinates;
3272   PetscFunctionReturn(0);
3273 }
3274 
3275 #undef __FUNCT__
3276 #define __FUNCT__ "DMGetCoordinatesLocal"
3277 /*@
3278   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM.
3279 
3280   Collective on DM
3281 
3282   Input Parameter:
3283 . dm - the DM
3284 
3285   Output Parameter:
3286 . c - coordinate vector
3287 
3288   Note:
3289   This is a borrowed reference, so the user should NOT destroy this vector
3290 
3291   Each process has the local and ghost coordinates
3292 
3293   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
3294   and (x_0,y_0,z_0,x_1,y_1,z_1...)
3295 
3296   Level: intermediate
3297 
3298 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3299 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMSetCoordinates(), DMGetCoordinateDM()
3300 @*/
3301 PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
3302 {
3303   PetscErrorCode ierr;
3304 
3305   PetscFunctionBegin;
3306   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3307   PetscValidPointer(c,2);
3308   if (!dm->coordinatesLocal && dm->coordinates) {
3309     DM cdm = NULL;
3310 
3311     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3312     ierr = DMCreateLocalVector(cdm, &dm->coordinatesLocal);CHKERRQ(ierr);
3313     ierr = PetscObjectSetName((PetscObject) dm->coordinatesLocal, "coordinates");CHKERRQ(ierr);
3314     ierr = DMGlobalToLocalBegin(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
3315     ierr = DMGlobalToLocalEnd(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
3316   }
3317   *c = dm->coordinatesLocal;
3318   PetscFunctionReturn(0);
3319 }
3320 
3321 #undef __FUNCT__
3322 #define __FUNCT__ "DMGetCoordinateDM"
3323 /*@
3324   DMGetCoordinateDM - Gets the DM that scatters between global and local coordinates
3325 
3326   Collective on DM
3327 
3328   Input Parameter:
3329 . dm - the DM
3330 
3331   Output Parameter:
3332 . cdm - coordinate DM
3333 
3334   Level: intermediate
3335 
3336 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3337 .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
3338 @*/
3339 PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
3340 {
3341   PetscErrorCode ierr;
3342 
3343   PetscFunctionBegin;
3344   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3345   PetscValidPointer(cdm,2);
3346   if (!dm->coordinateDM) {
3347     if (!dm->ops->createcoordinatedm) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unable to create coordinates for this DM");
3348     ierr = (*dm->ops->createcoordinatedm)(dm, &dm->coordinateDM);CHKERRQ(ierr);
3349   }
3350   *cdm = dm->coordinateDM;
3351   PetscFunctionReturn(0);
3352 }
3353 
3354 #undef __FUNCT__
3355 #define __FUNCT__ "DMLocatePoints"
3356 /*@
3357   DMLocatePoints - Locate the points in v in the mesh and return an IS of the containing cells
3358 
3359   Not collective
3360 
3361   Input Parameters:
3362 + dm - The DM
3363 - v - The Vec of points
3364 
3365   Output Parameter:
3366 . cells - The local cell numbers for cells which contain the points
3367 
3368   Level: developer
3369 
3370 .keywords: point location, mesh
3371 .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
3372 @*/
3373 PetscErrorCode DMLocatePoints(DM dm, Vec v, IS *cells)
3374 {
3375   PetscErrorCode ierr;
3376 
3377   PetscFunctionBegin;
3378   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3379   PetscValidHeaderSpecific(v,VEC_CLASSID,2);
3380   PetscValidPointer(cells,3);
3381   if (dm->ops->locatepoints) {
3382     ierr = (*dm->ops->locatepoints)(dm,v,cells);CHKERRQ(ierr);
3383   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Point location not available for this DM");
3384   PetscFunctionReturn(0);
3385 }
3386