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