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