xref: /petsc/src/dm/interface/dm.c (revision be7a6d03210faab0df8c3f3b720c976929326be8)
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 /*@
1793    DMLocalToLocalBegin - Maps from a local vector (including ghost points
1794    that contain irrelevant values) to another local vector where the ghost
1795    points in the second are set correctly. Must be followed by DMLocalToLocalEnd().
1796 
1797    Neighbor-wise Collective on DM and Vec
1798 
1799    Input Parameters:
1800 +  dm - the DM object
1801 .  g - the original local vector
1802 -  mode - one of INSERT_VALUES or ADD_VALUES
1803 
1804    Output Parameter:
1805 .  l  - the local vector with correct ghost values
1806 
1807    Level: intermediate
1808 
1809    Notes:
1810    The local vectors used here need not be the same as those
1811    obtained from DMCreateLocalVector(), BUT they
1812    must have the same parallel data layout; they could, for example, be
1813    obtained with VecDuplicate() from the DM originating vectors.
1814 
1815 .keywords: DM, local-to-local, begin
1816 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalEnd(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1817 
1818 @*/
1819 PetscErrorCode  DMLocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
1820 {
1821   PetscErrorCode          ierr;
1822 
1823   PetscFunctionBegin;
1824   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1825   ierr = (*dm->ops->localtolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1826   PetscFunctionReturn(0);
1827 }
1828 
1829 #undef __FUNCT__
1830 #define __FUNCT__ "DMLocalToLocalEnd"
1831 /*@
1832    DMLocalToLocalEnd - Maps from a local vector (including ghost points
1833    that contain irrelevant values) to another local vector where the ghost
1834    points in the second are set correctly. Must be preceded by DMLocalToLocalBegin().
1835 
1836    Neighbor-wise Collective on DM and Vec
1837 
1838    Input Parameters:
1839 +  da - the DM object
1840 .  g - the original local vector
1841 -  mode - one of INSERT_VALUES or ADD_VALUES
1842 
1843    Output Parameter:
1844 .  l  - the local vector with correct ghost values
1845 
1846    Level: intermediate
1847 
1848    Notes:
1849    The local vectors used here need not be the same as those
1850    obtained from DMCreateLocalVector(), BUT they
1851    must have the same parallel data layout; they could, for example, be
1852    obtained with VecDuplicate() from the DM originating vectors.
1853 
1854 .keywords: DM, local-to-local, end
1855 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1856 
1857 @*/
1858 PetscErrorCode  DMLocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
1859 {
1860   PetscErrorCode          ierr;
1861 
1862   PetscFunctionBegin;
1863   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1864   ierr = (*dm->ops->localtolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1865   PetscFunctionReturn(0);
1866 }
1867 
1868 
1869 #undef __FUNCT__
1870 #define __FUNCT__ "DMCoarsen"
1871 /*@
1872     DMCoarsen - Coarsens a DM object
1873 
1874     Collective on DM
1875 
1876     Input Parameter:
1877 +   dm - the DM object
1878 -   comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
1879 
1880     Output Parameter:
1881 .   dmc - the coarsened DM
1882 
1883     Level: developer
1884 
1885 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1886 
1887 @*/
1888 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
1889 {
1890   PetscErrorCode    ierr;
1891   DMCoarsenHookLink link;
1892 
1893   PetscFunctionBegin;
1894   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1895   ierr                      = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
1896   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
1897   ierr                      = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1898   (*dmc)->ctx               = dm->ctx;
1899   (*dmc)->levelup           = dm->levelup;
1900   (*dmc)->leveldown         = dm->leveldown + 1;
1901   ierr                      = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr);
1902   for (link=dm->coarsenhook; link; link=link->next) {
1903     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1904   }
1905   PetscFunctionReturn(0);
1906 }
1907 
1908 #undef __FUNCT__
1909 #define __FUNCT__ "DMCoarsenHookAdd"
1910 /*@C
1911    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1912 
1913    Logically Collective
1914 
1915    Input Arguments:
1916 +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1917 .  coarsenhook - function to run when setting up a coarser level
1918 .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
1919 -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
1920 
1921    Calling sequence of coarsenhook:
1922 $    coarsenhook(DM fine,DM coarse,void *ctx);
1923 
1924 +  fine - fine level DM
1925 .  coarse - coarse level DM to restrict problem to
1926 -  ctx - optional user-defined function context
1927 
1928    Calling sequence for restricthook:
1929 $    restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx)
1930 
1931 +  fine - fine level DM
1932 .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1933 .  rscale - scaling vector for restriction
1934 .  inject - matrix restricting by injection
1935 .  coarse - coarse level DM to update
1936 -  ctx - optional user-defined function context
1937 
1938    Level: advanced
1939 
1940    Notes:
1941    This function is only needed if auxiliary data needs to be set up on coarse grids.
1942 
1943    If this function is called multiple times, the hooks will be run in the order they are added.
1944 
1945    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1946    extract the finest level information from its context (instead of from the SNES).
1947 
1948    This function is currently not available from Fortran.
1949 
1950 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1951 @*/
1952 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1953 {
1954   PetscErrorCode    ierr;
1955   DMCoarsenHookLink link,*p;
1956 
1957   PetscFunctionBegin;
1958   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
1959   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1960   ierr               = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1961   link->coarsenhook  = coarsenhook;
1962   link->restricthook = restricthook;
1963   link->ctx          = ctx;
1964   link->next         = NULL;
1965   *p                 = link;
1966   PetscFunctionReturn(0);
1967 }
1968 
1969 #undef __FUNCT__
1970 #define __FUNCT__ "DMRestrict"
1971 /*@
1972    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1973 
1974    Collective if any hooks are
1975 
1976    Input Arguments:
1977 +  fine - finer DM to use as a base
1978 .  restrct - restriction matrix, apply using MatRestrict()
1979 .  inject - injection matrix, also use MatRestrict()
1980 -  coarse - coarer DM to update
1981 
1982    Level: developer
1983 
1984 .seealso: DMCoarsenHookAdd(), MatRestrict()
1985 @*/
1986 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1987 {
1988   PetscErrorCode    ierr;
1989   DMCoarsenHookLink link;
1990 
1991   PetscFunctionBegin;
1992   for (link=fine->coarsenhook; link; link=link->next) {
1993     if (link->restricthook) {
1994       ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);
1995     }
1996   }
1997   PetscFunctionReturn(0);
1998 }
1999 
2000 #undef __FUNCT__
2001 #define __FUNCT__ "DMSubDomainHookAdd"
2002 /*@C
2003    DMSubDomainHookAdd - adds a callback to be run when restricting a problem to the coarse grid
2004 
2005    Logically Collective
2006 
2007    Input Arguments:
2008 +  global - global DM
2009 .  ddhook - function to run to pass data to the decomposition DM upon its creation
2010 .  restricthook - function to run to update data on block solve (at the beginning of the block solve)
2011 -  ctx - [optional] user-defined context for provide data for the hooks (may be NULL)
2012 
2013 
2014    Calling sequence for ddhook:
2015 $    ddhook(DM global,DM block,void *ctx)
2016 
2017 +  global - global DM
2018 .  block  - block DM
2019 -  ctx - optional user-defined function context
2020 
2021    Calling sequence for restricthook:
2022 $    restricthook(DM global,VecScatter out,VecScatter in,DM block,void *ctx)
2023 
2024 +  global - global DM
2025 .  out    - scatter to the outer (with ghost and overlap points) block vector
2026 .  in     - scatter to block vector values only owned locally
2027 .  block  - block DM
2028 -  ctx - optional user-defined function context
2029 
2030    Level: advanced
2031 
2032    Notes:
2033    This function is only needed if auxiliary data needs to be set up on subdomain DMs.
2034 
2035    If this function is called multiple times, the hooks will be run in the order they are added.
2036 
2037    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
2038    extract the global information from its context (instead of from the SNES).
2039 
2040    This function is currently not available from Fortran.
2041 
2042 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
2043 @*/
2044 PetscErrorCode DMSubDomainHookAdd(DM global,PetscErrorCode (*ddhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*),void *ctx)
2045 {
2046   PetscErrorCode      ierr;
2047   DMSubDomainHookLink link,*p;
2048 
2049   PetscFunctionBegin;
2050   PetscValidHeaderSpecific(global,DM_CLASSID,1);
2051   for (p=&global->subdomainhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
2052   ierr               = PetscMalloc(sizeof(struct _DMSubDomainHookLink),&link);CHKERRQ(ierr);
2053   link->restricthook = restricthook;
2054   link->ddhook       = ddhook;
2055   link->ctx          = ctx;
2056   link->next         = NULL;
2057   *p                 = link;
2058   PetscFunctionReturn(0);
2059 }
2060 
2061 #undef __FUNCT__
2062 #define __FUNCT__ "DMSubDomainRestrict"
2063 /*@
2064    DMSubDomainRestrict - restricts user-defined problem data to a block DM by running hooks registered by DMSubDomainHookAdd()
2065 
2066    Collective if any hooks are
2067 
2068    Input Arguments:
2069 +  fine - finer DM to use as a base
2070 .  oscatter - scatter from domain global vector filling subdomain global vector with overlap
2071 .  gscatter - scatter from domain global vector filling subdomain local vector with ghosts
2072 -  coarse - coarer DM to update
2073 
2074    Level: developer
2075 
2076 .seealso: DMCoarsenHookAdd(), MatRestrict()
2077 @*/
2078 PetscErrorCode DMSubDomainRestrict(DM global,VecScatter oscatter,VecScatter gscatter,DM subdm)
2079 {
2080   PetscErrorCode      ierr;
2081   DMSubDomainHookLink link;
2082 
2083   PetscFunctionBegin;
2084   for (link=global->subdomainhook; link; link=link->next) {
2085     if (link->restricthook) {
2086       ierr = (*link->restricthook)(global,oscatter,gscatter,subdm,link->ctx);CHKERRQ(ierr);
2087     }
2088   }
2089   PetscFunctionReturn(0);
2090 }
2091 
2092 #undef __FUNCT__
2093 #define __FUNCT__ "DMGetCoarsenLevel"
2094 /*@
2095     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
2096 
2097     Not Collective
2098 
2099     Input Parameter:
2100 .   dm - the DM object
2101 
2102     Output Parameter:
2103 .   level - number of coarsenings
2104 
2105     Level: developer
2106 
2107 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
2108 
2109 @*/
2110 PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
2111 {
2112   PetscFunctionBegin;
2113   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2114   *level = dm->leveldown;
2115   PetscFunctionReturn(0);
2116 }
2117 
2118 
2119 
2120 #undef __FUNCT__
2121 #define __FUNCT__ "DMRefineHierarchy"
2122 /*@C
2123     DMRefineHierarchy - Refines a DM object, all levels at once
2124 
2125     Collective on DM
2126 
2127     Input Parameter:
2128 +   dm - the DM object
2129 -   nlevels - the number of levels of refinement
2130 
2131     Output Parameter:
2132 .   dmf - the refined DM hierarchy
2133 
2134     Level: developer
2135 
2136 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
2137 
2138 @*/
2139 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
2140 {
2141   PetscErrorCode ierr;
2142 
2143   PetscFunctionBegin;
2144   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2145   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
2146   if (nlevels == 0) PetscFunctionReturn(0);
2147   if (dm->ops->refinehierarchy) {
2148     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
2149   } else if (dm->ops->refine) {
2150     PetscInt i;
2151 
2152     ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr);
2153     for (i=1; i<nlevels; i++) {
2154       ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr);
2155     }
2156   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
2157   PetscFunctionReturn(0);
2158 }
2159 
2160 #undef __FUNCT__
2161 #define __FUNCT__ "DMCoarsenHierarchy"
2162 /*@C
2163     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
2164 
2165     Collective on DM
2166 
2167     Input Parameter:
2168 +   dm - the DM object
2169 -   nlevels - the number of levels of coarsening
2170 
2171     Output Parameter:
2172 .   dmc - the coarsened DM hierarchy
2173 
2174     Level: developer
2175 
2176 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
2177 
2178 @*/
2179 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
2180 {
2181   PetscErrorCode ierr;
2182 
2183   PetscFunctionBegin;
2184   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2185   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
2186   if (nlevels == 0) PetscFunctionReturn(0);
2187   PetscValidPointer(dmc,3);
2188   if (dm->ops->coarsenhierarchy) {
2189     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
2190   } else if (dm->ops->coarsen) {
2191     PetscInt i;
2192 
2193     ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr);
2194     for (i=1; i<nlevels; i++) {
2195       ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr);
2196     }
2197   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
2198   PetscFunctionReturn(0);
2199 }
2200 
2201 #undef __FUNCT__
2202 #define __FUNCT__ "DMCreateAggregates"
2203 /*@
2204    DMCreateAggregates - Gets the aggregates that map between
2205    grids associated with two DMs.
2206 
2207    Collective on DM
2208 
2209    Input Parameters:
2210 +  dmc - the coarse grid DM
2211 -  dmf - the fine grid DM
2212 
2213    Output Parameters:
2214 .  rest - the restriction matrix (transpose of the projection matrix)
2215 
2216    Level: intermediate
2217 
2218 .keywords: interpolation, restriction, multigrid
2219 
2220 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
2221 @*/
2222 PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
2223 {
2224   PetscErrorCode ierr;
2225 
2226   PetscFunctionBegin;
2227   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
2228   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
2229   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
2230   PetscFunctionReturn(0);
2231 }
2232 
2233 #undef __FUNCT__
2234 #define __FUNCT__ "DMSetApplicationContextDestroy"
2235 /*@C
2236     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
2237 
2238     Not Collective
2239 
2240     Input Parameters:
2241 +   dm - the DM object
2242 -   destroy - the destroy function
2243 
2244     Level: intermediate
2245 
2246 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2247 
2248 @*/
2249 PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
2250 {
2251   PetscFunctionBegin;
2252   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2253   dm->ctxdestroy = destroy;
2254   PetscFunctionReturn(0);
2255 }
2256 
2257 #undef __FUNCT__
2258 #define __FUNCT__ "DMSetApplicationContext"
2259 /*@
2260     DMSetApplicationContext - Set a user context into a DM object
2261 
2262     Not Collective
2263 
2264     Input Parameters:
2265 +   dm - the DM object
2266 -   ctx - the user context
2267 
2268     Level: intermediate
2269 
2270 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2271 
2272 @*/
2273 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
2274 {
2275   PetscFunctionBegin;
2276   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2277   dm->ctx = ctx;
2278   PetscFunctionReturn(0);
2279 }
2280 
2281 #undef __FUNCT__
2282 #define __FUNCT__ "DMGetApplicationContext"
2283 /*@
2284     DMGetApplicationContext - Gets a user context from a DM object
2285 
2286     Not Collective
2287 
2288     Input Parameter:
2289 .   dm - the DM object
2290 
2291     Output Parameter:
2292 .   ctx - the user context
2293 
2294     Level: intermediate
2295 
2296 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2297 
2298 @*/
2299 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
2300 {
2301   PetscFunctionBegin;
2302   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2303   *(void**)ctx = dm->ctx;
2304   PetscFunctionReturn(0);
2305 }
2306 
2307 #undef __FUNCT__
2308 #define __FUNCT__ "DMSetVariableBounds"
2309 /*@C
2310     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
2311 
2312     Logically Collective on DM
2313 
2314     Input Parameter:
2315 +   dm - the DM object
2316 -   f - the function that computes variable bounds used by SNESVI (use NULL to cancel a previous function that was set)
2317 
2318     Level: intermediate
2319 
2320 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(),
2321          DMSetJacobian()
2322 
2323 @*/
2324 PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
2325 {
2326   PetscFunctionBegin;
2327   dm->ops->computevariablebounds = f;
2328   PetscFunctionReturn(0);
2329 }
2330 
2331 #undef __FUNCT__
2332 #define __FUNCT__ "DMHasVariableBounds"
2333 /*@
2334     DMHasVariableBounds - does the DM object have a variable bounds function?
2335 
2336     Not Collective
2337 
2338     Input Parameter:
2339 .   dm - the DM object to destroy
2340 
2341     Output Parameter:
2342 .   flg - PETSC_TRUE if the variable bounds function exists
2343 
2344     Level: developer
2345 
2346 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2347 
2348 @*/
2349 PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
2350 {
2351   PetscFunctionBegin;
2352   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
2353   PetscFunctionReturn(0);
2354 }
2355 
2356 #undef __FUNCT__
2357 #define __FUNCT__ "DMComputeVariableBounds"
2358 /*@C
2359     DMComputeVariableBounds - compute variable bounds used by SNESVI.
2360 
2361     Logically Collective on DM
2362 
2363     Input Parameters:
2364 +   dm - the DM object to destroy
2365 -   x  - current solution at which the bounds are computed
2366 
2367     Output parameters:
2368 +   xl - lower bound
2369 -   xu - upper bound
2370 
2371     Level: intermediate
2372 
2373 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2374 
2375 @*/
2376 PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
2377 {
2378   PetscErrorCode ierr;
2379 
2380   PetscFunctionBegin;
2381   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
2382   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
2383   if (dm->ops->computevariablebounds) {
2384     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu);CHKERRQ(ierr);
2385   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "This DM is incapable of computing variable bounds.");
2386   PetscFunctionReturn(0);
2387 }
2388 
2389 #undef __FUNCT__
2390 #define __FUNCT__ "DMHasColoring"
2391 /*@
2392     DMHasColoring - does the DM object have a method of providing a coloring?
2393 
2394     Not Collective
2395 
2396     Input Parameter:
2397 .   dm - the DM object
2398 
2399     Output Parameter:
2400 .   flg - PETSC_TRUE if the DM has facilities for DMCreateColoring().
2401 
2402     Level: developer
2403 
2404 .seealso DMHasFunction(), DMCreateColoring()
2405 
2406 @*/
2407 PetscErrorCode  DMHasColoring(DM dm,PetscBool  *flg)
2408 {
2409   PetscFunctionBegin;
2410   *flg =  (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE;
2411   PetscFunctionReturn(0);
2412 }
2413 
2414 #undef  __FUNCT__
2415 #define __FUNCT__ "DMSetVec"
2416 /*@C
2417     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
2418 
2419     Collective on DM
2420 
2421     Input Parameter:
2422 +   dm - the DM object
2423 -   x - location to compute residual and Jacobian, if NULL is passed to those routines; will be NULL for linear problems.
2424 
2425     Level: developer
2426 
2427 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
2428 
2429 @*/
2430 PetscErrorCode  DMSetVec(DM dm,Vec x)
2431 {
2432   PetscErrorCode ierr;
2433 
2434   PetscFunctionBegin;
2435   if (x) {
2436     if (!dm->x) {
2437       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
2438     }
2439     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
2440   } else if (dm->x) {
2441     ierr = VecDestroy(&dm->x);CHKERRQ(ierr);
2442   }
2443   PetscFunctionReturn(0);
2444 }
2445 
2446 PetscFunctionList DMList              = NULL;
2447 PetscBool         DMRegisterAllCalled = PETSC_FALSE;
2448 
2449 #undef __FUNCT__
2450 #define __FUNCT__ "DMSetType"
2451 /*@C
2452   DMSetType - Builds a DM, for a particular DM implementation.
2453 
2454   Collective on DM
2455 
2456   Input Parameters:
2457 + dm     - The DM object
2458 - method - The name of the DM type
2459 
2460   Options Database Key:
2461 . -dm_type <type> - Sets the DM type; use -help for a list of available types
2462 
2463   Notes:
2464   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
2465 
2466   Level: intermediate
2467 
2468 .keywords: DM, set, type
2469 .seealso: DMGetType(), DMCreate()
2470 @*/
2471 PetscErrorCode  DMSetType(DM dm, DMType method)
2472 {
2473   PetscErrorCode (*r)(DM);
2474   PetscBool      match;
2475   PetscErrorCode ierr;
2476 
2477   PetscFunctionBegin;
2478   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2479   ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
2480   if (match) PetscFunctionReturn(0);
2481 
2482   if (!DMRegisterAllCalled) {ierr = DMRegisterAll();CHKERRQ(ierr);}
2483   ierr = PetscFunctionListFind(DMList,method,&r);CHKERRQ(ierr);
2484   if (!r) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
2485 
2486   if (dm->ops->destroy) {
2487     ierr             = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
2488     dm->ops->destroy = NULL;
2489   }
2490   ierr = (*r)(dm);CHKERRQ(ierr);
2491   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
2492   PetscFunctionReturn(0);
2493 }
2494 
2495 #undef __FUNCT__
2496 #define __FUNCT__ "DMGetType"
2497 /*@C
2498   DMGetType - Gets the DM type name (as a string) from the DM.
2499 
2500   Not Collective
2501 
2502   Input Parameter:
2503 . dm  - The DM
2504 
2505   Output Parameter:
2506 . type - The DM type name
2507 
2508   Level: intermediate
2509 
2510 .keywords: DM, get, type, name
2511 .seealso: DMSetType(), DMCreate()
2512 @*/
2513 PetscErrorCode  DMGetType(DM dm, DMType *type)
2514 {
2515   PetscErrorCode ierr;
2516 
2517   PetscFunctionBegin;
2518   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
2519   PetscValidCharPointer(type,2);
2520   if (!DMRegisterAllCalled) {
2521     ierr = DMRegisterAll();CHKERRQ(ierr);
2522   }
2523   *type = ((PetscObject)dm)->type_name;
2524   PetscFunctionReturn(0);
2525 }
2526 
2527 #undef __FUNCT__
2528 #define __FUNCT__ "DMConvert"
2529 /*@C
2530   DMConvert - Converts a DM to another DM, either of the same or different type.
2531 
2532   Collective on DM
2533 
2534   Input Parameters:
2535 + dm - the DM
2536 - newtype - new DM type (use "same" for the same type)
2537 
2538   Output Parameter:
2539 . M - pointer to new DM
2540 
2541   Notes:
2542   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
2543   the MPI communicator of the generated DM is always the same as the communicator
2544   of the input DM.
2545 
2546   Level: intermediate
2547 
2548 .seealso: DMCreate()
2549 @*/
2550 PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M)
2551 {
2552   DM             B;
2553   char           convname[256];
2554   PetscBool      sametype, issame;
2555   PetscErrorCode ierr;
2556 
2557   PetscFunctionBegin;
2558   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2559   PetscValidType(dm,1);
2560   PetscValidPointer(M,3);
2561   ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
2562   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
2563   {
2564     PetscErrorCode (*conv)(DM, DMType, DM*) = NULL;
2565 
2566     /*
2567        Order of precedence:
2568        1) See if a specialized converter is known to the current DM.
2569        2) See if a specialized converter is known to the desired DM class.
2570        3) See if a good general converter is registered for the desired class
2571        4) See if a good general converter is known for the current matrix.
2572        5) Use a really basic converter.
2573     */
2574 
2575     /* 1) See if a specialized converter is known to the current DM and the desired class */
2576     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
2577     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
2578     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
2579     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2580     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2581     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,&conv);CHKERRQ(ierr);
2582     if (conv) goto foundconv;
2583 
2584     /* 2)  See if a specialized converter is known to the desired DM class. */
2585     ierr = DMCreate(PetscObjectComm((PetscObject)dm), &B);CHKERRQ(ierr);
2586     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
2587     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
2588     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
2589     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
2590     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2591     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2592     ierr = PetscObjectQueryFunction((PetscObject)B,convname,&conv);CHKERRQ(ierr);
2593     if (conv) {
2594       ierr = DMDestroy(&B);CHKERRQ(ierr);
2595       goto foundconv;
2596     }
2597 
2598 #if 0
2599     /* 3) See if a good general converter is registered for the desired class */
2600     conv = B->ops->convertfrom;
2601     ierr = DMDestroy(&B);CHKERRQ(ierr);
2602     if (conv) goto foundconv;
2603 
2604     /* 4) See if a good general converter is known for the current matrix */
2605     if (dm->ops->convert) {
2606       conv = dm->ops->convert;
2607     }
2608     if (conv) goto foundconv;
2609 #endif
2610 
2611     /* 5) Use a really basic converter. */
2612     SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
2613 
2614 foundconv:
2615     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2616     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
2617     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2618   }
2619   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
2620   PetscFunctionReturn(0);
2621 }
2622 
2623 /*--------------------------------------------------------------------------------------------------------------------*/
2624 
2625 #undef __FUNCT__
2626 #define __FUNCT__ "DMRegister"
2627 /*@C
2628   DMRegister -  Adds a new DM component implementation
2629 
2630   Not Collective
2631 
2632   Input Parameters:
2633 + name        - The name of a new user-defined creation routine
2634 - create_func - The creation routine itself
2635 
2636   Notes:
2637   DMRegister() may be called multiple times to add several user-defined DMs
2638 
2639 
2640   Sample usage:
2641 .vb
2642     DMRegister("my_da", MyDMCreate);
2643 .ve
2644 
2645   Then, your DM type can be chosen with the procedural interface via
2646 .vb
2647     DMCreate(MPI_Comm, DM *);
2648     DMSetType(DM,"my_da");
2649 .ve
2650    or at runtime via the option
2651 .vb
2652     -da_type my_da
2653 .ve
2654 
2655   Level: advanced
2656 
2657 .keywords: DM, register
2658 .seealso: DMRegisterAll(), DMRegisterDestroy()
2659 
2660 @*/
2661 PetscErrorCode  DMRegister(const char sname[],PetscErrorCode (*function)(DM))
2662 {
2663   PetscErrorCode ierr;
2664 
2665   PetscFunctionBegin;
2666   ierr = PetscFunctionListAdd(&DMList,sname,function);CHKERRQ(ierr);
2667   PetscFunctionReturn(0);
2668 }
2669 
2670 #undef __FUNCT__
2671 #define __FUNCT__ "DMLoad"
2672 /*@C
2673   DMLoad - Loads a DM that has been stored in binary  with DMView().
2674 
2675   Collective on PetscViewer
2676 
2677   Input Parameters:
2678 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2679            some related function before a call to DMLoad().
2680 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2681            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2682 
2683    Level: intermediate
2684 
2685   Notes:
2686    The type is determined by the data in the file, any type set into the DM before this call is ignored.
2687 
2688   Notes for advanced users:
2689   Most users should not need to know the details of the binary storage
2690   format, since DMLoad() and DMView() completely hide these details.
2691   But for anyone who's interested, the standard binary matrix storage
2692   format is
2693 .vb
2694      has not yet been determined
2695 .ve
2696 
2697 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2698 @*/
2699 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2700 {
2701   PetscErrorCode ierr;
2702   PetscBool      isbinary;
2703   PetscInt       classid;
2704   char           type[256];
2705 
2706   PetscFunctionBegin;
2707   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2708   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2709   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
2710   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
2711 
2712   ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
2713   if (classid != DM_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not DM next in file");
2714   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
2715   ierr = DMSetType(newdm, type);CHKERRQ(ierr);
2716   if (newdm->ops->load) {
2717     ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2718   }
2719   PetscFunctionReturn(0);
2720 }
2721 
2722 /******************************** FEM Support **********************************/
2723 
2724 #undef __FUNCT__
2725 #define __FUNCT__ "DMPrintCellVector"
2726 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[])
2727 {
2728   PetscInt       f;
2729   PetscErrorCode ierr;
2730 
2731   PetscFunctionBegin;
2732   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2733   for (f = 0; f < len; ++f) {
2734     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
2735   }
2736   PetscFunctionReturn(0);
2737 }
2738 
2739 #undef __FUNCT__
2740 #define __FUNCT__ "DMPrintCellMatrix"
2741 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
2742 {
2743   PetscInt       f, g;
2744   PetscErrorCode ierr;
2745 
2746   PetscFunctionBegin;
2747   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2748   for (f = 0; f < rows; ++f) {
2749     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
2750     for (g = 0; g < cols; ++g) {
2751       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
2752     }
2753     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
2754   }
2755   PetscFunctionReturn(0);
2756 }
2757 
2758 #undef __FUNCT__
2759 #define __FUNCT__ "DMGetDefaultSection"
2760 /*@
2761   DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM.
2762 
2763   Input Parameter:
2764 . dm - The DM
2765 
2766   Output Parameter:
2767 . section - The PetscSection
2768 
2769   Level: intermediate
2770 
2771   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
2772 
2773 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
2774 @*/
2775 PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section)
2776 {
2777   PetscFunctionBegin;
2778   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2779   PetscValidPointer(section, 2);
2780   *section = dm->defaultSection;
2781   PetscFunctionReturn(0);
2782 }
2783 
2784 #undef __FUNCT__
2785 #define __FUNCT__ "DMSetDefaultSection"
2786 /*@
2787   DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM.
2788 
2789   Input Parameters:
2790 + dm - The DM
2791 - section - The PetscSection
2792 
2793   Level: intermediate
2794 
2795   Note: Any existing Section will be destroyed
2796 
2797 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
2798 @*/
2799 PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section)
2800 {
2801   PetscInt       numFields;
2802   PetscInt       f;
2803   PetscErrorCode ierr;
2804 
2805   PetscFunctionBegin;
2806   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2807   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
2808   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
2809   ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr);
2810   dm->defaultSection = section;
2811   ierr = PetscSectionGetNumFields(dm->defaultSection, &numFields);CHKERRQ(ierr);
2812   if (numFields) {
2813     ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr);
2814     for (f = 0; f < numFields; ++f) {
2815       const char *name;
2816 
2817       ierr = PetscSectionGetFieldName(dm->defaultSection, f, &name);CHKERRQ(ierr);
2818       ierr = PetscObjectSetName(dm->fields[f], name);CHKERRQ(ierr);
2819     }
2820   }
2821   /* The global section will be rebuilt in the next call to DMGetDefaultGlobalSection(). */
2822   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
2823   PetscFunctionReturn(0);
2824 }
2825 
2826 #undef __FUNCT__
2827 #define __FUNCT__ "DMGetDefaultGlobalSection"
2828 /*@
2829   DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM.
2830 
2831   Collective on DM
2832 
2833   Input Parameter:
2834 . dm - The DM
2835 
2836   Output Parameter:
2837 . section - The PetscSection
2838 
2839   Level: intermediate
2840 
2841   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
2842 
2843 .seealso: DMSetDefaultSection(), DMGetDefaultSection()
2844 @*/
2845 PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section)
2846 {
2847   PetscErrorCode ierr;
2848 
2849   PetscFunctionBegin;
2850   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2851   PetscValidPointer(section, 2);
2852   if (!dm->defaultGlobalSection) {
2853     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");
2854     ierr = PetscSectionCreateGlobalSection(dm->defaultSection, dm->sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr);
2855     ierr = PetscLayoutDestroy(&dm->map);CHKERRQ(ierr);
2856     ierr = PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm),dm->defaultGlobalSection,&dm->map);CHKERRQ(ierr);
2857   }
2858   *section = dm->defaultGlobalSection;
2859   PetscFunctionReturn(0);
2860 }
2861 
2862 #undef __FUNCT__
2863 #define __FUNCT__ "DMSetDefaultGlobalSection"
2864 /*@
2865   DMSetDefaultGlobalSection - Set the PetscSection encoding the global data layout for the DM.
2866 
2867   Input Parameters:
2868 + dm - The DM
2869 - section - The PetscSection, or NULL
2870 
2871   Level: intermediate
2872 
2873   Note: Any existing Section will be destroyed
2874 
2875 .seealso: DMGetDefaultGlobalSection(), DMSetDefaultSection()
2876 @*/
2877 PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection section)
2878 {
2879   PetscErrorCode ierr;
2880 
2881   PetscFunctionBegin;
2882   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2883   if (section) PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2);
2884   ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr);
2885   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
2886   dm->defaultGlobalSection = section;
2887   PetscFunctionReturn(0);
2888 }
2889 
2890 #undef __FUNCT__
2891 #define __FUNCT__ "DMGetDefaultSF"
2892 /*@
2893   DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set,
2894   it is created from the default PetscSection layouts in the DM.
2895 
2896   Input Parameter:
2897 . dm - The DM
2898 
2899   Output Parameter:
2900 . sf - The PetscSF
2901 
2902   Level: intermediate
2903 
2904   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
2905 
2906 .seealso: DMSetDefaultSF(), DMCreateDefaultSF()
2907 @*/
2908 PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf)
2909 {
2910   PetscInt       nroots;
2911   PetscErrorCode ierr;
2912 
2913   PetscFunctionBegin;
2914   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2915   PetscValidPointer(sf, 2);
2916   ierr = PetscSFGetGraph(dm->defaultSF, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
2917   if (nroots < 0) {
2918     PetscSection section, gSection;
2919 
2920     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
2921     if (section) {
2922       ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
2923       ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr);
2924     } else {
2925       *sf = NULL;
2926       PetscFunctionReturn(0);
2927     }
2928   }
2929   *sf = dm->defaultSF;
2930   PetscFunctionReturn(0);
2931 }
2932 
2933 #undef __FUNCT__
2934 #define __FUNCT__ "DMSetDefaultSF"
2935 /*@
2936   DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM
2937 
2938   Input Parameters:
2939 + dm - The DM
2940 - sf - The PetscSF
2941 
2942   Level: intermediate
2943 
2944   Note: Any previous SF is destroyed
2945 
2946 .seealso: DMGetDefaultSF(), DMCreateDefaultSF()
2947 @*/
2948 PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf)
2949 {
2950   PetscErrorCode ierr;
2951 
2952   PetscFunctionBegin;
2953   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2954   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
2955   ierr          = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr);
2956   dm->defaultSF = sf;
2957   PetscFunctionReturn(0);
2958 }
2959 
2960 #undef __FUNCT__
2961 #define __FUNCT__ "DMCreateDefaultSF"
2962 /*@C
2963   DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections
2964   describing the data layout.
2965 
2966   Input Parameters:
2967 + dm - The DM
2968 . localSection - PetscSection describing the local data layout
2969 - globalSection - PetscSection describing the global data layout
2970 
2971   Level: intermediate
2972 
2973 .seealso: DMGetDefaultSF(), DMSetDefaultSF()
2974 @*/
2975 PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection)
2976 {
2977   MPI_Comm       comm;
2978   PetscLayout    layout;
2979   const PetscInt *ranges;
2980   PetscInt       *local;
2981   PetscSFNode    *remote;
2982   PetscInt       pStart, pEnd, p, nroots, nleaves, l;
2983   PetscMPIInt    size, rank;
2984   PetscErrorCode ierr;
2985 
2986   PetscFunctionBegin;
2987   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
2988   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2989   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
2990   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2991   ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr);
2992   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr);
2993   ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
2994   ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
2995   ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr);
2996   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
2997   ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
2998   for (p = pStart, nleaves = 0; p < pEnd; ++p) {
2999     PetscInt gdof, gcdof;
3000 
3001     ierr     = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
3002     ierr     = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
3003     nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
3004   }
3005   ierr = PetscMalloc(nleaves * sizeof(PetscInt), &local);CHKERRQ(ierr);
3006   ierr = PetscMalloc(nleaves * sizeof(PetscSFNode), &remote);CHKERRQ(ierr);
3007   for (p = pStart, l = 0; p < pEnd; ++p) {
3008     const PetscInt *cind;
3009     PetscInt       dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
3010 
3011     ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr);
3012     ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr);
3013     ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr);
3014     ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr);
3015     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
3016     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
3017     ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
3018     if (!gdof) continue; /* Censored point */
3019     gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
3020     if (gsize != dof-cdof) {
3021       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);
3022       cdof = 0; /* Ignore constraints */
3023     }
3024     for (d = 0, c = 0; d < dof; ++d) {
3025       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
3026       local[l+d-c] = off+d;
3027     }
3028     if (gdof < 0) {
3029       for (d = 0; d < gsize; ++d, ++l) {
3030         PetscInt offset = -(goff+1) + d, r;
3031 
3032         ierr = PetscFindInt(offset,size,ranges,&r);CHKERRQ(ierr);
3033         if (r < 0) r = -(r+2);
3034         remote[l].rank  = r;
3035         remote[l].index = offset - ranges[r];
3036       }
3037     } else {
3038       for (d = 0; d < gsize; ++d, ++l) {
3039         remote[l].rank  = rank;
3040         remote[l].index = goff+d - ranges[rank];
3041       }
3042     }
3043   }
3044   if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %d != nleaves %d", l, nleaves);
3045   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
3046   ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
3047   PetscFunctionReturn(0);
3048 }
3049 
3050 #undef __FUNCT__
3051 #define __FUNCT__ "DMGetPointSF"
3052 /*@
3053   DMGetPointSF - Get the PetscSF encoding the parallel section point overlap for the DM.
3054 
3055   Input Parameter:
3056 . dm - The DM
3057 
3058   Output Parameter:
3059 . sf - The PetscSF
3060 
3061   Level: intermediate
3062 
3063   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
3064 
3065 .seealso: DMSetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3066 @*/
3067 PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf)
3068 {
3069   PetscFunctionBegin;
3070   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3071   PetscValidPointer(sf, 2);
3072   *sf = dm->sf;
3073   PetscFunctionReturn(0);
3074 }
3075 
3076 #undef __FUNCT__
3077 #define __FUNCT__ "DMSetPointSF"
3078 /*@
3079   DMSetPointSF - Set the PetscSF encoding the parallel section point overlap for the DM.
3080 
3081   Input Parameters:
3082 + dm - The DM
3083 - sf - The PetscSF
3084 
3085   Level: intermediate
3086 
3087 .seealso: DMGetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF()
3088 @*/
3089 PetscErrorCode DMSetPointSF(DM dm, PetscSF sf)
3090 {
3091   PetscErrorCode ierr;
3092 
3093   PetscFunctionBegin;
3094   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3095   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
3096   ierr   = PetscSFDestroy(&dm->sf);CHKERRQ(ierr);
3097   ierr   = PetscObjectReference((PetscObject) sf);CHKERRQ(ierr);
3098   dm->sf = sf;
3099   PetscFunctionReturn(0);
3100 }
3101 
3102 #undef __FUNCT__
3103 #define __FUNCT__ "DMGetNumFields"
3104 PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields)
3105 {
3106   PetscFunctionBegin;
3107   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3108   PetscValidPointer(numFields, 2);
3109   *numFields = dm->numFields;
3110   PetscFunctionReturn(0);
3111 }
3112 
3113 #undef __FUNCT__
3114 #define __FUNCT__ "DMSetNumFields"
3115 PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields)
3116 {
3117   PetscInt       f;
3118   PetscErrorCode ierr;
3119 
3120   PetscFunctionBegin;
3121   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3122   for (f = 0; f < dm->numFields; ++f) {
3123     ierr = PetscObjectDestroy(&dm->fields[f]);CHKERRQ(ierr);
3124   }
3125   ierr          = PetscFree(dm->fields);CHKERRQ(ierr);
3126   dm->numFields = numFields;
3127   ierr          = PetscMalloc(dm->numFields * sizeof(PetscObject), &dm->fields);CHKERRQ(ierr);
3128   for (f = 0; f < dm->numFields; ++f) {
3129     ierr = PetscContainerCreate(PetscObjectComm((PetscObject)dm), (PetscContainer*) &dm->fields[f]);CHKERRQ(ierr);
3130   }
3131   PetscFunctionReturn(0);
3132 }
3133 
3134 #undef __FUNCT__
3135 #define __FUNCT__ "DMGetField"
3136 PetscErrorCode DMGetField(DM dm, PetscInt f, PetscObject *field)
3137 {
3138   PetscFunctionBegin;
3139   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3140   PetscValidPointer(field, 2);
3141   if (!dm->fields) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Fields have not been setup in this DM. Call DMSetNumFields()");
3142   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);
3143   *field = dm->fields[f];
3144   PetscFunctionReturn(0);
3145 }
3146 
3147 #undef __FUNCT__
3148 #define __FUNCT__ "DMRestrictHook_Coordinates"
3149 PetscErrorCode DMRestrictHook_Coordinates(DM dm,DM dmc,void *ctx)
3150 {
3151   DM dm_coord,dmc_coord;
3152   PetscErrorCode ierr;
3153   Vec coords,ccoords;
3154   VecScatter scat;
3155   PetscFunctionBegin;
3156   ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr);
3157   ierr = DMGetCoordinateDM(dmc,&dmc_coord);CHKERRQ(ierr);
3158   ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr);
3159   ierr = DMGetCoordinates(dmc,&ccoords);CHKERRQ(ierr);
3160   if (coords && !ccoords) {
3161     ierr = DMCreateGlobalVector(dmc_coord,&ccoords);CHKERRQ(ierr);
3162     ierr = DMCreateInjection(dmc_coord,dm_coord,&scat);CHKERRQ(ierr);
3163     ierr = VecScatterBegin(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3164     ierr = VecScatterEnd(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3165     ierr = DMSetCoordinates(dmc,ccoords);CHKERRQ(ierr);
3166     ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
3167     ierr = VecDestroy(&ccoords);CHKERRQ(ierr);
3168   }
3169   PetscFunctionReturn(0);
3170 }
3171 
3172 #undef __FUNCT__
3173 #define __FUNCT__ "DMSubDomainHook_Coordinates"
3174 static PetscErrorCode DMSubDomainHook_Coordinates(DM dm,DM subdm,void *ctx)
3175 {
3176   DM dm_coord,subdm_coord;
3177   PetscErrorCode ierr;
3178   Vec coords,ccoords,clcoords;
3179   VecScatter *scat_i,*scat_g;
3180   PetscFunctionBegin;
3181   ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr);
3182   ierr = DMGetCoordinateDM(subdm,&subdm_coord);CHKERRQ(ierr);
3183   ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr);
3184   ierr = DMGetCoordinates(subdm,&ccoords);CHKERRQ(ierr);
3185   if (coords && !ccoords) {
3186     ierr = DMCreateGlobalVector(subdm_coord,&ccoords);CHKERRQ(ierr);
3187     ierr = DMCreateLocalVector(subdm_coord,&clcoords);CHKERRQ(ierr);
3188     ierr = DMCreateDomainDecompositionScatters(dm_coord,1,&subdm_coord,NULL,&scat_i,&scat_g);CHKERRQ(ierr);
3189     ierr = VecScatterBegin(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3190     ierr = VecScatterBegin(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3191     ierr = VecScatterEnd(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3192     ierr = VecScatterEnd(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3193     ierr = DMSetCoordinates(subdm,ccoords);CHKERRQ(ierr);
3194     ierr = DMSetCoordinatesLocal(subdm,clcoords);CHKERRQ(ierr);
3195     ierr = VecScatterDestroy(&scat_i[0]);CHKERRQ(ierr);
3196     ierr = VecScatterDestroy(&scat_g[0]);CHKERRQ(ierr);
3197     ierr = VecDestroy(&ccoords);CHKERRQ(ierr);
3198     ierr = VecDestroy(&clcoords);CHKERRQ(ierr);
3199     ierr = PetscFree(scat_i);CHKERRQ(ierr);
3200     ierr = PetscFree(scat_g);CHKERRQ(ierr);
3201   }
3202   PetscFunctionReturn(0);
3203 }
3204 
3205 #undef __FUNCT__
3206 #define __FUNCT__ "DMSetCoordinates"
3207 /*@
3208   DMSetCoordinates - Sets into the DM a global vector that holds the coordinates
3209 
3210   Collective on DM
3211 
3212   Input Parameters:
3213 + dm - the DM
3214 - c - coordinate vector
3215 
3216   Note:
3217   The coordinates do include those for ghost points, which are in the local vector
3218 
3219   Level: intermediate
3220 
3221 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3222 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLoca(), DMGetCoordinateDM()
3223 @*/
3224 PetscErrorCode DMSetCoordinates(DM dm, Vec c)
3225 {
3226   PetscErrorCode ierr;
3227 
3228   PetscFunctionBegin;
3229   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3230   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
3231   ierr            = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
3232   ierr            = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
3233   dm->coordinates = c;
3234   ierr            = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
3235   ierr            = DMCoarsenHookAdd(dm,DMRestrictHook_Coordinates,NULL,NULL);CHKERRQ(ierr);
3236   ierr            = DMSubDomainHookAdd(dm,DMSubDomainHook_Coordinates,NULL,NULL);CHKERRQ(ierr);
3237   PetscFunctionReturn(0);
3238 }
3239 
3240 #undef __FUNCT__
3241 #define __FUNCT__ "DMSetCoordinatesLocal"
3242 /*@
3243   DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates
3244 
3245   Collective on DM
3246 
3247    Input Parameters:
3248 +  dm - the DM
3249 -  c - coordinate vector
3250 
3251   Note:
3252   The coordinates of ghost points can be set using DMSetCoordinates()
3253   followed by DMGetCoordinatesLocal(). This is intended to enable the
3254   setting of ghost coordinates outside of the domain.
3255 
3256   Level: intermediate
3257 
3258 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3259 .seealso: DMGetCoordinatesLocal(), DMSetCoordinates(), DMGetCoordinates(), DMGetCoordinateDM()
3260 @*/
3261 PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
3262 {
3263   PetscErrorCode ierr;
3264 
3265   PetscFunctionBegin;
3266   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3267   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
3268   ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr);
3269   ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr);
3270 
3271   dm->coordinatesLocal = c;
3272 
3273   ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
3274   PetscFunctionReturn(0);
3275 }
3276 
3277 #undef __FUNCT__
3278 #define __FUNCT__ "DMGetCoordinates"
3279 /*@
3280   DMGetCoordinates - Gets a global vector with the coordinates associated with the DM.
3281 
3282   Not Collective
3283 
3284   Input Parameter:
3285 . dm - the DM
3286 
3287   Output Parameter:
3288 . c - global coordinate vector
3289 
3290   Note:
3291   This is a borrowed reference, so the user should NOT destroy this vector
3292 
3293   Each process has only the local coordinates (does NOT have the ghost coordinates).
3294 
3295   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
3296   and (x_0,y_0,z_0,x_1,y_1,z_1...)
3297 
3298   Level: intermediate
3299 
3300 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3301 .seealso: DMSetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM()
3302 @*/
3303 PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
3304 {
3305   PetscErrorCode ierr;
3306 
3307   PetscFunctionBegin;
3308   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3309   PetscValidPointer(c,2);
3310   if (!dm->coordinates && dm->coordinatesLocal) {
3311     DM cdm = NULL;
3312 
3313     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3314     ierr = DMCreateGlobalVector(cdm, &dm->coordinates);CHKERRQ(ierr);
3315     ierr = PetscObjectSetName((PetscObject) dm->coordinates, "coordinates");CHKERRQ(ierr);
3316     ierr = DMLocalToGlobalBegin(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
3317     ierr = DMLocalToGlobalEnd(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr);
3318   }
3319   *c = dm->coordinates;
3320   PetscFunctionReturn(0);
3321 }
3322 
3323 #undef __FUNCT__
3324 #define __FUNCT__ "DMGetCoordinatesLocal"
3325 /*@
3326   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM.
3327 
3328   Collective on DM
3329 
3330   Input Parameter:
3331 . dm - the DM
3332 
3333   Output Parameter:
3334 . c - coordinate vector
3335 
3336   Note:
3337   This is a borrowed reference, so the user should NOT destroy this vector
3338 
3339   Each process has the local and ghost coordinates
3340 
3341   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
3342   and (x_0,y_0,z_0,x_1,y_1,z_1...)
3343 
3344   Level: intermediate
3345 
3346 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3347 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMSetCoordinates(), DMGetCoordinateDM()
3348 @*/
3349 PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
3350 {
3351   PetscErrorCode ierr;
3352 
3353   PetscFunctionBegin;
3354   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3355   PetscValidPointer(c,2);
3356   if (!dm->coordinatesLocal && dm->coordinates) {
3357     DM cdm = NULL;
3358 
3359     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
3360     ierr = DMCreateLocalVector(cdm, &dm->coordinatesLocal);CHKERRQ(ierr);
3361     ierr = PetscObjectSetName((PetscObject) dm->coordinatesLocal, "coordinates");CHKERRQ(ierr);
3362     ierr = DMGlobalToLocalBegin(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
3363     ierr = DMGlobalToLocalEnd(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr);
3364   }
3365   *c = dm->coordinatesLocal;
3366   PetscFunctionReturn(0);
3367 }
3368 
3369 #undef __FUNCT__
3370 #define __FUNCT__ "DMGetCoordinateDM"
3371 /*@
3372   DMGetCoordinateDM - Gets the DM that scatters between global and local coordinates
3373 
3374   Collective on DM
3375 
3376   Input Parameter:
3377 . dm - the DM
3378 
3379   Output Parameter:
3380 . cdm - coordinate DM
3381 
3382   Level: intermediate
3383 
3384 .keywords: distributed array, get, corners, nodes, local indices, coordinates
3385 .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
3386 @*/
3387 PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
3388 {
3389   PetscErrorCode ierr;
3390 
3391   PetscFunctionBegin;
3392   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3393   PetscValidPointer(cdm,2);
3394   if (!dm->coordinateDM) {
3395     if (!dm->ops->createcoordinatedm) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unable to create coordinates for this DM");
3396     ierr = (*dm->ops->createcoordinatedm)(dm, &dm->coordinateDM);CHKERRQ(ierr);
3397   }
3398   *cdm = dm->coordinateDM;
3399   PetscFunctionReturn(0);
3400 }
3401 
3402 #undef __FUNCT__
3403 #define __FUNCT__ "DMLocatePoints"
3404 /*@
3405   DMLocatePoints - Locate the points in v in the mesh and return an IS of the containing cells
3406 
3407   Not collective
3408 
3409   Input Parameters:
3410 + dm - The DM
3411 - v - The Vec of points
3412 
3413   Output Parameter:
3414 . cells - The local cell numbers for cells which contain the points
3415 
3416   Level: developer
3417 
3418 .keywords: point location, mesh
3419 .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal()
3420 @*/
3421 PetscErrorCode DMLocatePoints(DM dm, Vec v, IS *cells)
3422 {
3423   PetscErrorCode ierr;
3424 
3425   PetscFunctionBegin;
3426   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3427   PetscValidHeaderSpecific(v,VEC_CLASSID,2);
3428   PetscValidPointer(cells,3);
3429   if (dm->ops->locatepoints) {
3430     ierr = (*dm->ops->locatepoints)(dm,v,cells);CHKERRQ(ierr);
3431   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Point location not available for this DM");
3432   PetscFunctionReturn(0);
3433 }
3434