xref: /petsc/src/dm/interface/dm.c (revision 88ed4acec49ae3a81b39199da1871f2558231dad)
1 #include <petscsnes.h>
2 #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/
3 
4 PetscClassId  DM_CLASSID;
5 PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal;
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "DMCreate"
9 /*@
10   DMCreate - Creates an empty DM object. The type can then be set with DMSetType().
11 
12    If you never  call DMSetType()  it will generate an
13    error when you try to use the vector.
14 
15   Collective on MPI_Comm
16 
17   Input Parameter:
18 . comm - The communicator for the DM object
19 
20   Output Parameter:
21 . dm - The DM object
22 
23   Level: beginner
24 
25 .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE
26 @*/
27 PetscErrorCode  DMCreate(MPI_Comm comm,DM *dm)
28 {
29   DM             v;
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   PetscValidPointer(dm,2);
34   *dm = PETSC_NULL;
35 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
36   ierr = VecInitializePackage(PETSC_NULL);CHKERRQ(ierr);
37   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
38   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
39 #endif
40 
41   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, -1, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
42   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
43 
44 
45   v->workSize     = 0;
46   v->workArray    = PETSC_NULL;
47   v->ltogmap      = PETSC_NULL;
48   v->ltogmapb     = PETSC_NULL;
49   v->bs           = 1;
50   v->coloringtype = IS_COLORING_GLOBAL;
51   v->lf           = PETSC_NULL;
52   v->lj           = PETSC_NULL;
53   ierr = PetscSFCreate(comm, &v->sf);CHKERRQ(ierr);
54   ierr = PetscSFCreate(comm, &v->defaultSF);CHKERRQ(ierr);
55   v->defaultSection       = PETSC_NULL;
56   v->defaultGlobalSection = PETSC_NULL;
57 
58   *dm = v;
59   PetscFunctionReturn(0);
60 }
61 
62 
63 #undef __FUNCT__
64 #define __FUNCT__ "DMSetVecType"
65 /*@C
66        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
67 
68    Logically Collective on DMDA
69 
70    Input Parameter:
71 +  da - initial distributed array
72 .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
73 
74    Options Database:
75 .   -dm_vec_type ctype
76 
77    Level: intermediate
78 
79 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
80 @*/
81 PetscErrorCode  DMSetVecType(DM da,const VecType ctype)
82 {
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   PetscValidHeaderSpecific(da,DM_CLASSID,1);
87   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
88   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "DMSetMatType"
94 /*@C
95        DMSetMatType - Sets the type of matrix created with DMCreateMatrix()
96 
97    Logically Collective on DM
98 
99    Input Parameter:
100 +  dm - the DM context
101 .  ctype - the matrix type
102 
103    Options Database:
104 .   -dm_mat_type ctype
105 
106    Level: intermediate
107 
108 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType
109 @*/
110 PetscErrorCode  DMSetMatType(DM dm,const MatType ctype)
111 {
112   PetscErrorCode ierr;
113   PetscFunctionBegin;
114   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
115   ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
116   ierr = PetscStrallocpy(ctype,&dm->mattype);CHKERRQ(ierr);
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "DMSetOptionsPrefix"
122 /*@C
123    DMSetOptionsPrefix - Sets the prefix used for searching for all
124    DMDA options in the database.
125 
126    Logically Collective on DMDA
127 
128    Input Parameter:
129 +  da - the DMDA context
130 -  prefix - the prefix to prepend to all option names
131 
132    Notes:
133    A hyphen (-) must NOT be given at the beginning of the prefix name.
134    The first character of all runtime options is AUTOMATICALLY the hyphen.
135 
136    Level: advanced
137 
138 .keywords: DMDA, set, options, prefix, database
139 
140 .seealso: DMSetFromOptions()
141 @*/
142 PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
143 {
144   PetscErrorCode ierr;
145 
146   PetscFunctionBegin;
147   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
148   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
149   PetscFunctionReturn(0);
150 }
151 
152 #undef __FUNCT__
153 #define __FUNCT__ "DMDestroy"
154 /*@
155     DMDestroy - Destroys a vector packer or DMDA.
156 
157     Collective on DM
158 
159     Input Parameter:
160 .   dm - the DM object to destroy
161 
162     Level: developer
163 
164 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
165 
166 @*/
167 PetscErrorCode  DMDestroy(DM *dm)
168 {
169   PetscInt       i, cnt = 0;
170   DMCoarsenHookLink link,next;
171   DMNamedVecLink nlink,nnext;
172   PetscErrorCode ierr;
173 
174   PetscFunctionBegin;
175   if (!*dm) PetscFunctionReturn(0);
176   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
177 
178   /* count all the circular references of DM and its contained Vecs */
179   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
180     if ((*dm)->localin[i])  {cnt++;}
181     if ((*dm)->globalin[i]) {cnt++;}
182   }
183   for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++;
184   if ((*dm)->x) {
185     PetscObject obj;
186     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
187     if (obj == (PetscObject)*dm) cnt++;
188   }
189 
190   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
191   /*
192      Need this test because the dm references the vectors that
193      reference the dm, so destroying the dm calls destroy on the
194      vectors that cause another destroy on the dm
195   */
196   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
197   ((PetscObject) (*dm))->refct = 0;
198   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
199     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
200     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
201   }
202   for (nlink=(*dm)->namedglobal; nlink; nlink=nnext) { /* Destroy the named vectors */
203     nnext = nlink->next;
204     if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name);
205     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
206     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
207     ierr = PetscFree(nlink);CHKERRQ(ierr);
208   }
209   (*dm)->namedglobal = PETSC_NULL;
210 
211   /* Destroy the list of hooks */
212   for (link=(*dm)->coarsenhook; link; link=next) {
213     next = link->next;
214     ierr = PetscFree(link);CHKERRQ(ierr);
215   }
216   (*dm)->coarsenhook = PETSC_NULL;
217 
218   if ((*dm)->ctx && (*dm)->ctxdestroy) {
219     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
220   }
221   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
222   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
223   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
224   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
225   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
226   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
227   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
228   ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr);
229 
230   ierr = PetscSectionDestroy(&(*dm)->defaultSection);CHKERRQ(ierr);
231   ierr = PetscSectionDestroy(&(*dm)->defaultGlobalSection);CHKERRQ(ierr);
232   ierr = PetscSFDestroy(&(*dm)->sf);CHKERRQ(ierr);
233   ierr = PetscSFDestroy(&(*dm)->defaultSF);CHKERRQ(ierr);
234   /* if memory was published with AMS then destroy it */
235   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
236 
237   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
238   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
239   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
240   PetscFunctionReturn(0);
241 }
242 
243 #undef __FUNCT__
244 #define __FUNCT__ "DMSetUp"
245 /*@
246     DMSetUp - sets up the data structures inside a DM object
247 
248     Collective on DM
249 
250     Input Parameter:
251 .   dm - the DM object to setup
252 
253     Level: developer
254 
255 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
256 
257 @*/
258 PetscErrorCode  DMSetUp(DM dm)
259 {
260   PetscErrorCode ierr;
261 
262   PetscFunctionBegin;
263   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
264   if (dm->setupcalled) PetscFunctionReturn(0);
265   if (dm->ops->setup) {
266     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
267   }
268   dm->setupcalled = PETSC_TRUE;
269   PetscFunctionReturn(0);
270 }
271 
272 #undef __FUNCT__
273 #define __FUNCT__ "DMSetFromOptions"
274 /*@
275     DMSetFromOptions - sets parameters in a DM from the options database
276 
277     Collective on DM
278 
279     Input Parameter:
280 .   dm - the DM object to set options for
281 
282     Options Database:
283 +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
284 .   -dm_vec_type <type>  type of vector to create inside DM
285 .   -dm_mat_type <type>  type of matrix to create inside DM
286 -   -dm_coloring_type <global or ghosted>
287 
288     Level: developer
289 
290 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
291 
292 @*/
293 PetscErrorCode  DMSetFromOptions(DM dm)
294 {
295   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg;
296   PetscErrorCode ierr;
297   char           typeName[256] = MATAIJ;
298 
299   PetscFunctionBegin;
300   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
301   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
302     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
303     ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr);
304     ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr);
305     ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr);
306     ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,PETSC_NULL);CHKERRQ(ierr);
307     ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
308     if (flg) {
309       ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
310     }
311     ierr = PetscOptionsList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype?dm->mattype:typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr);
312     if (flg) {
313       ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr);
314     }
315     ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr);
316     if (dm->ops->setfromoptions) {
317       ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
318     }
319     /* process any options handlers added with PetscObjectAddOptionsHandler() */
320     ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
321   ierr = PetscOptionsEnd();CHKERRQ(ierr);
322   if (flg1) {
323     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
324   }
325   if (flg2) {
326     PetscViewer viewer;
327 
328     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
329     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
330     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
331     ierr = DMView(dm, viewer);CHKERRQ(ierr);
332     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
333   }
334   if (flg3) {
335     PetscViewer viewer;
336 
337     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
338     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
339     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
340     ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr);
341     ierr = DMView(dm, viewer);CHKERRQ(ierr);
342     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
343   }
344   if (flg4) {
345     PetscViewer viewer;
346 
347     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
348     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
349     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr);
350     ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr);
351     ierr = DMView(dm, viewer);CHKERRQ(ierr);
352     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
353   }
354   PetscFunctionReturn(0);
355 }
356 
357 #undef __FUNCT__
358 #define __FUNCT__ "DMView"
359 /*@C
360     DMView - Views a vector packer or DMDA.
361 
362     Collective on DM
363 
364     Input Parameter:
365 +   dm - the DM object to view
366 -   v - the viewer
367 
368     Level: developer
369 
370 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
371 
372 @*/
373 PetscErrorCode  DMView(DM dm,PetscViewer v)
374 {
375   PetscErrorCode ierr;
376 
377   PetscFunctionBegin;
378   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
379  if (!v) {
380     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
381   }
382   if (dm->ops->view) {
383     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
384   }
385   PetscFunctionReturn(0);
386 }
387 
388 #undef __FUNCT__
389 #define __FUNCT__ "DMCreateGlobalVector"
390 /*@
391     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
392 
393     Collective on DM
394 
395     Input Parameter:
396 .   dm - the DM object
397 
398     Output Parameter:
399 .   vec - the global vector
400 
401     Level: beginner
402 
403 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
404 
405 @*/
406 PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
407 {
408   PetscErrorCode ierr;
409 
410   PetscFunctionBegin;
411   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
412   if (dm->defaultSection) {
413     PetscSection gSection;
414     PetscInt     localSize;
415 
416     ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
417     ierr = PetscSectionGetConstrainedStorageSize(dm->defaultGlobalSection, &localSize);CHKERRQ(ierr);
418     ierr = VecCreate(((PetscObject) dm)->comm, vec);CHKERRQ(ierr);
419     ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
420     /* ierr = VecSetType(*vec, dm->vectype);CHKERRQ(ierr); */
421     ierr = VecSetFromOptions(*vec);CHKERRQ(ierr);
422     ierr = PetscObjectCompose((PetscObject) *vec, "DM", (PetscObject) dm);CHKERRQ(ierr);
423     /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
424     /* ierr = VecSetLocalToGlobalMappingBlock(*vec, dm->ltogmapb);CHKERRQ(ierr); */
425     /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */
426   } else {
427     ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
428   }
429   PetscFunctionReturn(0);
430 }
431 
432 #undef __FUNCT__
433 #define __FUNCT__ "DMCreateLocalVector"
434 /*@
435     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
436 
437     Not Collective
438 
439     Input Parameter:
440 .   dm - the DM object
441 
442     Output Parameter:
443 .   vec - the local vector
444 
445     Level: beginner
446 
447 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
448 
449 @*/
450 PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
451 {
452   PetscErrorCode ierr;
453 
454   PetscFunctionBegin;
455   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
456   if (dm->defaultSection) {
457     PetscInt localSize;
458 
459     ierr = PetscSectionGetStorageSize(dm->defaultSection, &localSize);CHKERRQ(ierr);
460     ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
461     ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
462     ierr = VecSetFromOptions(*vec);CHKERRQ(ierr);
463     ierr = PetscObjectCompose((PetscObject) *vec, "DM", (PetscObject) dm);CHKERRQ(ierr);
464   } else {
465     ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
466   }
467   PetscFunctionReturn(0);
468 }
469 
470 #undef __FUNCT__
471 #define __FUNCT__ "DMGetLocalToGlobalMapping"
472 /*@
473    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
474 
475    Collective on DM
476 
477    Input Parameter:
478 .  dm - the DM that provides the mapping
479 
480    Output Parameter:
481 .  ltog - the mapping
482 
483    Level: intermediate
484 
485    Notes:
486    This mapping can then be used by VecSetLocalToGlobalMapping() or
487    MatSetLocalToGlobalMapping().
488 
489 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
490 @*/
491 PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
492 {
493   PetscErrorCode ierr;
494 
495   PetscFunctionBegin;
496   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
497   PetscValidPointer(ltog,2);
498   if (!dm->ltogmap) {
499     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
500     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
501   }
502   *ltog = dm->ltogmap;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
508 /*@
509    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
510 
511    Collective on DM
512 
513    Input Parameter:
514 .  da - the distributed array that provides the mapping
515 
516    Output Parameter:
517 .  ltog - the block mapping
518 
519    Level: intermediate
520 
521    Notes:
522    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
523    MatSetLocalToGlobalMappingBlock().
524 
525 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
526 @*/
527 PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
528 {
529   PetscErrorCode ierr;
530 
531   PetscFunctionBegin;
532   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
533   PetscValidPointer(ltog,2);
534   if (!dm->ltogmapb) {
535     PetscInt bs;
536     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
537     if (bs > 1) {
538       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
539       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
540     } else {
541       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
542       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
543     }
544   }
545   *ltog = dm->ltogmapb;
546   PetscFunctionReturn(0);
547 }
548 
549 #undef __FUNCT__
550 #define __FUNCT__ "DMGetBlockSize"
551 /*@
552    DMGetBlockSize - Gets the inherent block size associated with a DM
553 
554    Not Collective
555 
556    Input Parameter:
557 .  dm - the DM with block structure
558 
559    Output Parameter:
560 .  bs - the block size, 1 implies no exploitable block structure
561 
562    Level: intermediate
563 
564 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
565 @*/
566 PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
567 {
568   PetscFunctionBegin;
569   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
570   PetscValidPointer(bs,2);
571   if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet");
572   *bs = dm->bs;
573   PetscFunctionReturn(0);
574 }
575 
576 #undef __FUNCT__
577 #define __FUNCT__ "DMCreateInterpolation"
578 /*@
579     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
580 
581     Collective on DM
582 
583     Input Parameter:
584 +   dm1 - the DM object
585 -   dm2 - the second, finer DM object
586 
587     Output Parameter:
588 +  mat - the interpolation
589 -  vec - the scaling (optional)
590 
591     Level: developer
592 
593     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
594         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
595 
596         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
597         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
598 
599 
600 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
601 
602 @*/
603 PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
604 {
605   PetscErrorCode ierr;
606 
607   PetscFunctionBegin;
608   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
609   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
610   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
611   PetscFunctionReturn(0);
612 }
613 
614 #undef __FUNCT__
615 #define __FUNCT__ "DMCreateInjection"
616 /*@
617     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
618 
619     Collective on DM
620 
621     Input Parameter:
622 +   dm1 - the DM object
623 -   dm2 - the second, finer DM object
624 
625     Output Parameter:
626 .   ctx - the injection
627 
628     Level: developer
629 
630    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
631         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
632 
633 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
634 
635 @*/
636 PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
637 {
638   PetscErrorCode ierr;
639 
640   PetscFunctionBegin;
641   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
642   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
643   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
644   PetscFunctionReturn(0);
645 }
646 
647 #undef __FUNCT__
648 #define __FUNCT__ "DMCreateColoring"
649 /*@C
650     DMCreateColoring - Gets coloring for a DMDA or DMComposite
651 
652     Collective on DM
653 
654     Input Parameter:
655 +   dm - the DM object
656 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
657 -   matype - either MATAIJ or MATBAIJ
658 
659     Output Parameter:
660 .   coloring - the coloring
661 
662     Level: developer
663 
664 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
665 
666 @*/
667 PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
668 {
669   PetscErrorCode ierr;
670 
671   PetscFunctionBegin;
672   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
673   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
674   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
675   PetscFunctionReturn(0);
676 }
677 
678 #undef __FUNCT__
679 #define __FUNCT__ "DMCreateMatrix"
680 /*@C
681     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
682 
683     Collective on DM
684 
685     Input Parameter:
686 +   dm - the DM object
687 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
688             any type which inherits from one of these (such as MATAIJ)
689 
690     Output Parameter:
691 .   mat - the empty Jacobian
692 
693     Level: beginner
694 
695     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
696        do not need to do it yourself.
697 
698        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
699        the nonzero pattern call DMDASetMatPreallocateOnly()
700 
701        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
702        internally by PETSc.
703 
704        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
705        the indices for the global numbering for DMDAs which is complicated.
706 
707 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
708 
709 @*/
710 PetscErrorCode  DMCreateMatrix(DM dm,const MatType mtype,Mat *mat)
711 {
712   PetscErrorCode ierr;
713 
714   PetscFunctionBegin;
715   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
716 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
717   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
718 #endif
719   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
720   PetscValidPointer(mat,3);
721   if (dm->mattype) {
722     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
723   } else {
724     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
725   }
726   PetscFunctionReturn(0);
727 }
728 
729 #undef __FUNCT__
730 #define __FUNCT__ "DMSetMatrixPreallocateOnly"
731 /*@
732   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
733     preallocated but the nonzero structure and zero values will not be set.
734 
735   Logically Collective on DMDA
736 
737   Input Parameter:
738 + dm - the DM
739 - only - PETSC_TRUE if only want preallocation
740 
741   Level: developer
742 .seealso DMCreateMatrix()
743 @*/
744 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
745 {
746   PetscFunctionBegin;
747   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
748   dm->prealloc_only = only;
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "DMGetWorkArray"
754 /*@C
755   DMGetWorkArray - Gets a work array guaranteed to be at least the input size
756 
757   Not Collective
758 
759   Input Parameters:
760 + dm - the DM object
761 - size - The minium size
762 
763   Output Parameter:
764 . array - the work array
765 
766   Level: developer
767 
768 .seealso DMDestroy(), DMCreate()
769 @*/
770 PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array)
771 {
772   PetscErrorCode ierr;
773 
774   PetscFunctionBegin;
775   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
776   PetscValidPointer(array,3);
777   if (size > dm->workSize) {
778     dm->workSize = size;
779     ierr = PetscFree(dm->workArray);CHKERRQ(ierr);
780     ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr);
781   }
782   *array = dm->workArray;
783   PetscFunctionReturn(0);
784 }
785 
786 #undef __FUNCT__
787 #define __FUNCT__ "DMCreateDecompositionDM"
788 /*@C
789   DMCreateDecompositionDM - creates a DM that encapsulates a decomposition of the original DM.
790 
791   Not Collective
792 
793   Input Parameters:
794 + dm   - the DM object
795 - name - the name of the decomposition
796 
797   Output Parameter:
798 . ddm  - the decomposition DM (PETSC_NULL, if no such decomposition is known)
799 
800   Level: advanced
801 
802 .seealso DMDestroy(), DMCreate(), DMCreateDecomposition()
803 @*/
804 PetscErrorCode DMCreateDecompositionDM(DM dm, const char* name, DM *ddm)
805 {
806   PetscErrorCode ierr;
807 
808   PetscFunctionBegin;
809   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
810   PetscValidCharPointer(name,2);
811   PetscValidPointer(ddm,3);
812   if(!dm->ops->createdecompositiondm) {
813     *ddm = PETSC_NULL;
814   }
815   else {
816     ierr = (*dm->ops->createdecompositiondm)(dm,name,ddm); CHKERRQ(ierr);
817   }
818   PetscFunctionReturn(0);
819 }
820 
821 #undef __FUNCT__
822 #define __FUNCT__ "DMCreateFieldIS"
823 /*@C
824   DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field
825 
826   Not collective
827 
828   Input Parameter:
829 . dm - the DM object
830 
831   Output Parameters:
832 + numFields - The number of fields (or PETSC_NULL if not requested)
833 . names     - The name for each field (or PETSC_NULL if not requested)
834 - fields    - The global indices for each field (or PETSC_NULL if not requested)
835 
836   Level: intermediate
837 
838   Notes:
839   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
840   PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with
841   PetscFree().
842 
843 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
844 @*/
845 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***names, IS **fields)
846 {
847   PetscErrorCode ierr;
848 
849   PetscFunctionBegin;
850   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
851   if (numFields) {
852     PetscValidPointer(numFields,2);
853     *numFields = 0;
854   }
855   if (names) {
856     PetscValidPointer(names,3);
857     *names = PETSC_NULL;
858   }
859   if (fields) {
860     PetscValidPointer(fields,4);
861     *fields = PETSC_NULL;
862   }
863   if(dm->ops->createfieldis) {
864     ierr = (*dm->ops->createfieldis)(dm, numFields, names, fields);CHKERRQ(ierr);
865   }
866   PetscFunctionReturn(0);
867 }
868 
869 
870 #undef __FUNCT__
871 #define __FUNCT__ "DMCreateDecomposition"
872 /*@C
873   DMCreateDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems:
874                           each IS contains the global indices of the dofs of the corresponding subproblem.
875                           The optional list of DMs define the DM for each subproblem.
876                           Generalizes DMCreateFieldIS().
877 
878   Not collective
879 
880   Input Parameter:
881 . dm - the DM object
882 
883   Output Parameters:
884 + len       - The number of subproblems in the decomposition (or PETSC_NULL if not requested)
885 . namelist  - The name for each subproblem (or PETSC_NULL if not requested)
886 . islist    - The global indices for each subproblem (or PETSC_NULL if not requested)
887 - dmlist    - The DMs for each subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined)
888 
889   Level: intermediate
890 
891   Notes:
892   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
893   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
894   and all of the arrays should be freed with PetscFree().
895 
896 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
897 @*/
898 PetscErrorCode DMCreateDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
899 {
900   PetscErrorCode ierr;
901 
902   PetscFunctionBegin;
903   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
904   if (len) {PetscValidPointer(len,2);}
905   if (namelist) {PetscValidPointer(namelist,3);}
906   if (islist) {PetscValidPointer(islist,4);}
907   if (dmlist) {PetscValidPointer(dmlist,5);}
908   if(!dm->ops->createdecomposition) {
909     ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr);
910     /* By default there are no DMs associated with subproblems. */
911     if(dmlist) *dmlist = PETSC_NULL;
912   }
913   else {
914     ierr = (*dm->ops->createdecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr);
915   }
916   PetscFunctionReturn(0);
917 }
918 
919 
920 #undef __FUNCT__
921 #define __FUNCT__ "DMRefine"
922 /*@
923   DMRefine - Refines a DM object
924 
925   Collective on DM
926 
927   Input Parameter:
928 + dm   - the DM object
929 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
930 
931   Output Parameter:
932 . dmf - the refined DM, or PETSC_NULL
933 
934   Note: If no refinement was done, the return value is PETSC_NULL
935 
936   Level: developer
937 
938 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
939 @*/
940 PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
941 {
942   PetscErrorCode ierr;
943 
944   PetscFunctionBegin;
945   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
946   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
947   if (*dmf) {
948     (*dmf)->ops->creatematrix = dm->ops->creatematrix;
949     (*dmf)->ops->initialguess = dm->ops->initialguess;
950     (*dmf)->ops->function     = dm->ops->function;
951     (*dmf)->ops->functionj    = dm->ops->functionj;
952     if (dm->ops->jacobian != DMComputeJacobianDefault) {
953       (*dmf)->ops->jacobian     = dm->ops->jacobian;
954     }
955     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
956     (*dmf)->ctx     = dm->ctx;
957     (*dmf)->levelup = dm->levelup + 1;
958   }
959   PetscFunctionReturn(0);
960 }
961 
962 #undef __FUNCT__
963 #define __FUNCT__ "DMGetRefineLevel"
964 /*@
965     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
966 
967     Not Collective
968 
969     Input Parameter:
970 .   dm - the DM object
971 
972     Output Parameter:
973 .   level - number of refinements
974 
975     Level: developer
976 
977 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
978 
979 @*/
980 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
981 {
982   PetscFunctionBegin;
983   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
984   *level = dm->levelup;
985   PetscFunctionReturn(0);
986 }
987 
988 #undef __FUNCT__
989 #define __FUNCT__ "DMGlobalToLocalBegin"
990 /*@
991     DMGlobalToLocalBegin - Begins updating local vectors from global vector
992 
993     Neighbor-wise Collective on DM
994 
995     Input Parameters:
996 +   dm - the DM object
997 .   g - the global vector
998 .   mode - INSERT_VALUES or ADD_VALUES
999 -   l - the local vector
1000 
1001 
1002     Level: beginner
1003 
1004 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1005 
1006 @*/
1007 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
1008 {
1009   PetscErrorCode ierr;
1010 
1011   PetscFunctionBegin;
1012   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1013   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 #undef __FUNCT__
1018 #define __FUNCT__ "DMGlobalToLocalEnd"
1019 /*@
1020     DMGlobalToLocalEnd - Ends updating local vectors from global vector
1021 
1022     Neighbor-wise Collective on DM
1023 
1024     Input Parameters:
1025 +   dm - the DM object
1026 .   g - the global vector
1027 .   mode - INSERT_VALUES or ADD_VALUES
1028 -   l - the local vector
1029 
1030 
1031     Level: beginner
1032 
1033 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
1034 
1035 @*/
1036 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
1037 {
1038   PetscErrorCode ierr;
1039 
1040   PetscFunctionBegin;
1041   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1042   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1043   PetscFunctionReturn(0);
1044 }
1045 
1046 #undef __FUNCT__
1047 #define __FUNCT__ "DMLocalToGlobalBegin"
1048 /*@
1049     DMLocalToGlobalBegin - updates global vectors from local vectors
1050 
1051     Neighbor-wise Collective on DM
1052 
1053     Input Parameters:
1054 +   dm - the DM object
1055 .   l - the local vector
1056 .   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
1057            base point.
1058 - - the global vector
1059 
1060     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
1061            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
1062            global array to the final global array with VecAXPY().
1063 
1064     Level: beginner
1065 
1066 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
1067 
1068 @*/
1069 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
1070 {
1071   PetscErrorCode ierr;
1072 
1073   PetscFunctionBegin;
1074   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1075   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 #undef __FUNCT__
1080 #define __FUNCT__ "DMLocalToGlobalEnd"
1081 /*@
1082     DMLocalToGlobalEnd - updates global vectors from local vectors
1083 
1084     Neighbor-wise Collective on DM
1085 
1086     Input Parameters:
1087 +   dm - the DM object
1088 .   l - the local vector
1089 .   mode - INSERT_VALUES or ADD_VALUES
1090 -   g - the global vector
1091 
1092 
1093     Level: beginner
1094 
1095 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
1096 
1097 @*/
1098 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
1099 {
1100   PetscErrorCode ierr;
1101 
1102   PetscFunctionBegin;
1103   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1104   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
1105   PetscFunctionReturn(0);
1106 }
1107 
1108 #undef __FUNCT__
1109 #define __FUNCT__ "DMComputeJacobianDefault"
1110 /*@
1111     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
1112 
1113     Collective on DM
1114 
1115     Input Parameter:
1116 +   dm - the DM object
1117 .   x - location to compute Jacobian at; may be ignored for linear problems
1118 .   A - matrix that defines the operator for the linear solve
1119 -   B - the matrix used to construct the preconditioner
1120 
1121     Level: developer
1122 
1123 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1124          DMSetFunction()
1125 
1126 @*/
1127 PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1128 {
1129   PetscErrorCode ierr;
1130 
1131   PetscFunctionBegin;
1132   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1133   *stflag = SAME_NONZERO_PATTERN;
1134   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
1135   if (A != B) {
1136     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1137     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1138   }
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 #undef __FUNCT__
1143 #define __FUNCT__ "DMCoarsen"
1144 /*@
1145     DMCoarsen - Coarsens a DM object
1146 
1147     Collective on DM
1148 
1149     Input Parameter:
1150 +   dm - the DM object
1151 -   comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
1152 
1153     Output Parameter:
1154 .   dmc - the coarsened DM
1155 
1156     Level: developer
1157 
1158 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1159 
1160 @*/
1161 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
1162 {
1163   PetscErrorCode ierr;
1164   DMCoarsenHookLink link;
1165 
1166   PetscFunctionBegin;
1167   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1168   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
1169   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
1170   (*dmc)->ops->initialguess = dm->ops->initialguess;
1171   (*dmc)->ops->function     = dm->ops->function;
1172   (*dmc)->ops->functionj    = dm->ops->functionj;
1173   if (dm->ops->jacobian != DMComputeJacobianDefault) {
1174     (*dmc)->ops->jacobian     = dm->ops->jacobian;
1175   }
1176   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1177   (*dmc)->ctx       = dm->ctx;
1178   (*dmc)->leveldown = dm->leveldown + 1;
1179   for (link=dm->coarsenhook; link; link=link->next) {
1180     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1181   }
1182   PetscFunctionReturn(0);
1183 }
1184 
1185 #undef __FUNCT__
1186 #define __FUNCT__ "DMCoarsenHookAdd"
1187 /*@
1188    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1189 
1190    Logically Collective
1191 
1192    Input Arguments:
1193 +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1194 .  coarsenhook - function to run when setting up a coarser level
1195 .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
1196 -  ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL)
1197 
1198    Calling sequence of coarsenhook:
1199 $    coarsenhook(DM fine,DM coarse,void *ctx);
1200 
1201 +  fine - fine level DM
1202 .  coarse - coarse level DM to restrict problem to
1203 -  ctx - optional user-defined function context
1204 
1205    Calling sequence for restricthook:
1206 $    restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx)
1207 
1208 +  fine - fine level DM
1209 .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1210 .  inject - matrix restricting by applying the transpose of injection
1211 .  coarse - coarse level DM to update
1212 -  ctx - optional user-defined function context
1213 
1214    Level: advanced
1215 
1216    Notes:
1217    This function is only needed if auxiliary data needs to be set up on coarse grids.
1218 
1219    If this function is called multiple times, the hooks will be run in the order they are added.
1220 
1221    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1222    extract the finest level information from its context (instead of from the SNES).
1223 
1224 .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1225 @*/
1226 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1227 {
1228   PetscErrorCode ierr;
1229   DMCoarsenHookLink link,*p;
1230 
1231   PetscFunctionBegin;
1232   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
1233   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1234   ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1235   link->coarsenhook = coarsenhook;
1236   link->restricthook = restricthook;
1237   link->ctx = ctx;
1238   link->next = PETSC_NULL;
1239   *p = link;
1240   PetscFunctionReturn(0);
1241 }
1242 
1243 #undef __FUNCT__
1244 #define __FUNCT__ "DMRestrict"
1245 /*@
1246    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1247 
1248    Collective if any hooks are
1249 
1250    Input Arguments:
1251 +  fine - finer DM to use as a base
1252 .  restrct - restriction matrix, apply using MatRestrict()
1253 .  inject - injection matrix, also use MatRestrict()
1254 -  coarse - coarer DM to update
1255 
1256    Level: developer
1257 
1258 .seealso: DMCoarsenHookAdd(), MatRestrict()
1259 @*/
1260 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1261 {
1262   PetscErrorCode ierr;
1263   DMCoarsenHookLink link;
1264 
1265   PetscFunctionBegin;
1266   for (link=fine->coarsenhook; link; link=link->next) {
1267     if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);}
1268   }
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 #undef __FUNCT__
1273 #define __FUNCT__ "DMGetCoarsenLevel"
1274 /*@
1275     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
1276 
1277     Not Collective
1278 
1279     Input Parameter:
1280 .   dm - the DM object
1281 
1282     Output Parameter:
1283 .   level - number of coarsenings
1284 
1285     Level: developer
1286 
1287 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1288 
1289 @*/
1290 PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
1291 {
1292   PetscFunctionBegin;
1293   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1294   *level = dm->leveldown;
1295   PetscFunctionReturn(0);
1296 }
1297 
1298 
1299 
1300 #undef __FUNCT__
1301 #define __FUNCT__ "DMRefineHierarchy"
1302 /*@C
1303     DMRefineHierarchy - Refines a DM object, all levels at once
1304 
1305     Collective on DM
1306 
1307     Input Parameter:
1308 +   dm - the DM object
1309 -   nlevels - the number of levels of refinement
1310 
1311     Output Parameter:
1312 .   dmf - the refined DM hierarchy
1313 
1314     Level: developer
1315 
1316 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1317 
1318 @*/
1319 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
1320 {
1321   PetscErrorCode ierr;
1322 
1323   PetscFunctionBegin;
1324   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1325   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1326   if (nlevels == 0) PetscFunctionReturn(0);
1327   if (dm->ops->refinehierarchy) {
1328     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
1329   } else if (dm->ops->refine) {
1330     PetscInt i;
1331 
1332     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
1333     for (i=1; i<nlevels; i++) {
1334       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
1335     }
1336   } else {
1337     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
1338   }
1339   PetscFunctionReturn(0);
1340 }
1341 
1342 #undef __FUNCT__
1343 #define __FUNCT__ "DMCoarsenHierarchy"
1344 /*@C
1345     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
1346 
1347     Collective on DM
1348 
1349     Input Parameter:
1350 +   dm - the DM object
1351 -   nlevels - the number of levels of coarsening
1352 
1353     Output Parameter:
1354 .   dmc - the coarsened DM hierarchy
1355 
1356     Level: developer
1357 
1358 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1359 
1360 @*/
1361 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
1362 {
1363   PetscErrorCode ierr;
1364 
1365   PetscFunctionBegin;
1366   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1367   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1368   if (nlevels == 0) PetscFunctionReturn(0);
1369   PetscValidPointer(dmc,3);
1370   if (dm->ops->coarsenhierarchy) {
1371     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
1372   } else if (dm->ops->coarsen) {
1373     PetscInt i;
1374 
1375     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
1376     for (i=1; i<nlevels; i++) {
1377       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
1378     }
1379   } else {
1380     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
1381   }
1382   PetscFunctionReturn(0);
1383 }
1384 
1385 #undef __FUNCT__
1386 #define __FUNCT__ "DMCreateAggregates"
1387 /*@
1388    DMCreateAggregates - Gets the aggregates that map between
1389    grids associated with two DMs.
1390 
1391    Collective on DM
1392 
1393    Input Parameters:
1394 +  dmc - the coarse grid DM
1395 -  dmf - the fine grid DM
1396 
1397    Output Parameters:
1398 .  rest - the restriction matrix (transpose of the projection matrix)
1399 
1400    Level: intermediate
1401 
1402 .keywords: interpolation, restriction, multigrid
1403 
1404 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
1405 @*/
1406 PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
1407 {
1408   PetscErrorCode ierr;
1409 
1410   PetscFunctionBegin;
1411   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1412   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
1413   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
1414   PetscFunctionReturn(0);
1415 }
1416 
1417 #undef __FUNCT__
1418 #define __FUNCT__ "DMSetApplicationContextDestroy"
1419 /*@C
1420     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
1421 
1422     Not Collective
1423 
1424     Input Parameters:
1425 +   dm - the DM object
1426 -   destroy - the destroy function
1427 
1428     Level: intermediate
1429 
1430 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1431 
1432 @*/
1433 PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
1434 {
1435   PetscFunctionBegin;
1436   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1437   dm->ctxdestroy = destroy;
1438   PetscFunctionReturn(0);
1439 }
1440 
1441 #undef __FUNCT__
1442 #define __FUNCT__ "DMSetApplicationContext"
1443 /*@
1444     DMSetApplicationContext - Set a user context into a DM object
1445 
1446     Not Collective
1447 
1448     Input Parameters:
1449 +   dm - the DM object
1450 -   ctx - the user context
1451 
1452     Level: intermediate
1453 
1454 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1455 
1456 @*/
1457 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
1458 {
1459   PetscFunctionBegin;
1460   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1461   dm->ctx = ctx;
1462   PetscFunctionReturn(0);
1463 }
1464 
1465 #undef __FUNCT__
1466 #define __FUNCT__ "DMGetApplicationContext"
1467 /*@
1468     DMGetApplicationContext - Gets a user context from a DM object
1469 
1470     Not Collective
1471 
1472     Input Parameter:
1473 .   dm - the DM object
1474 
1475     Output Parameter:
1476 .   ctx - the user context
1477 
1478     Level: intermediate
1479 
1480 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1481 
1482 @*/
1483 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
1484 {
1485   PetscFunctionBegin;
1486   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1487   *(void**)ctx = dm->ctx;
1488   PetscFunctionReturn(0);
1489 }
1490 
1491 #undef __FUNCT__
1492 #define __FUNCT__ "DMSetInitialGuess"
1493 /*@C
1494     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
1495 
1496     Logically Collective on DM
1497 
1498     Input Parameter:
1499 +   dm - the DM object to destroy
1500 -   f - the function to compute the initial guess
1501 
1502     Level: intermediate
1503 
1504 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1505 
1506 @*/
1507 PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
1508 {
1509   PetscFunctionBegin;
1510   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1511   dm->ops->initialguess = f;
1512   PetscFunctionReturn(0);
1513 }
1514 
1515 #undef __FUNCT__
1516 #define __FUNCT__ "DMSetFunction"
1517 /*@C
1518     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
1519 
1520     Logically Collective on DM
1521 
1522     Input Parameter:
1523 +   dm - the DM object
1524 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
1525 
1526     Level: intermediate
1527 
1528     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
1529            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
1530 
1531 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1532          DMSetJacobian()
1533 
1534 @*/
1535 PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1536 {
1537   PetscFunctionBegin;
1538   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1539   dm->ops->function = f;
1540   if (f) {
1541     dm->ops->functionj = f;
1542   }
1543   PetscFunctionReturn(0);
1544 }
1545 
1546 #undef __FUNCT__
1547 #define __FUNCT__ "DMSetJacobian"
1548 /*@C
1549     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
1550 
1551     Logically Collective on DM
1552 
1553     Input Parameter:
1554 +   dm - the DM object to destroy
1555 -   f - the function to compute the matrix entries
1556 
1557     Level: intermediate
1558 
1559 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1560          DMSetFunction()
1561 
1562 @*/
1563 PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
1564 {
1565   PetscFunctionBegin;
1566   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1567   dm->ops->jacobian = f;
1568   PetscFunctionReturn(0);
1569 }
1570 
1571 #undef __FUNCT__
1572 #define __FUNCT__ "DMSetVariableBounds"
1573 /*@C
1574     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
1575 
1576     Logically Collective on DM
1577 
1578     Input Parameter:
1579 +   dm - the DM object
1580 -   f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set)
1581 
1582     Level: intermediate
1583 
1584 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1585          DMSetJacobian()
1586 
1587 @*/
1588 PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1589 {
1590   PetscFunctionBegin;
1591   dm->ops->computevariablebounds = f;
1592   PetscFunctionReturn(0);
1593 }
1594 
1595 #undef __FUNCT__
1596 #define __FUNCT__ "DMHasVariableBounds"
1597 /*@
1598     DMHasVariableBounds - does the DM object have a variable bounds function?
1599 
1600     Not Collective
1601 
1602     Input Parameter:
1603 .   dm - the DM object to destroy
1604 
1605     Output Parameter:
1606 .   flg - PETSC_TRUE if the variable bounds function exists
1607 
1608     Level: developer
1609 
1610 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1611 
1612 @*/
1613 PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
1614 {
1615   PetscFunctionBegin;
1616   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
1617   PetscFunctionReturn(0);
1618 }
1619 
1620 #undef __FUNCT__
1621 #define __FUNCT__ "DMComputeVariableBounds"
1622 /*@C
1623     DMComputeVariableBounds - compute variable bounds used by SNESVI.
1624 
1625     Logically Collective on DM
1626 
1627     Input Parameters:
1628 +   dm - the DM object to destroy
1629 -   x  - current solution at which the bounds are computed
1630 
1631     Output parameters:
1632 +   xl - lower bound
1633 -   xu - upper bound
1634 
1635     Level: intermediate
1636 
1637 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1638          DMSetFunction(), DMSetVariableBounds()
1639 
1640 @*/
1641 PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
1642 {
1643   PetscErrorCode ierr;
1644   PetscFunctionBegin;
1645   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1646   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
1647   if(dm->ops->computevariablebounds) {
1648     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr);
1649   }
1650   else {
1651     ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr);
1652     ierr = VecSet(xu,SNES_VI_INF);  CHKERRQ(ierr);
1653   }
1654   PetscFunctionReturn(0);
1655 }
1656 
1657 #undef __FUNCT__
1658 #define __FUNCT__ "DMComputeInitialGuess"
1659 /*@
1660     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
1661 
1662     Collective on DM
1663 
1664     Input Parameter:
1665 +   dm - the DM object to destroy
1666 -   x - the vector to hold the initial guess values
1667 
1668     Level: developer
1669 
1670 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
1671 
1672 @*/
1673 PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
1674 {
1675   PetscErrorCode ierr;
1676   PetscFunctionBegin;
1677   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1678   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
1679   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
1680   PetscFunctionReturn(0);
1681 }
1682 
1683 #undef __FUNCT__
1684 #define __FUNCT__ "DMHasInitialGuess"
1685 /*@
1686     DMHasInitialGuess - does the DM object have an initial guess function
1687 
1688     Not Collective
1689 
1690     Input Parameter:
1691 .   dm - the DM object to destroy
1692 
1693     Output Parameter:
1694 .   flg - PETSC_TRUE if function exists
1695 
1696     Level: developer
1697 
1698 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1699 
1700 @*/
1701 PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
1702 {
1703   PetscFunctionBegin;
1704   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1705   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
1706   PetscFunctionReturn(0);
1707 }
1708 
1709 #undef __FUNCT__
1710 #define __FUNCT__ "DMHasFunction"
1711 /*@
1712     DMHasFunction - does the DM object have a function
1713 
1714     Not Collective
1715 
1716     Input Parameter:
1717 .   dm - the DM object to destroy
1718 
1719     Output Parameter:
1720 .   flg - PETSC_TRUE if function exists
1721 
1722     Level: developer
1723 
1724 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1725 
1726 @*/
1727 PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
1728 {
1729   PetscFunctionBegin;
1730   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1731   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
1732   PetscFunctionReturn(0);
1733 }
1734 
1735 #undef __FUNCT__
1736 #define __FUNCT__ "DMHasJacobian"
1737 /*@
1738     DMHasJacobian - does the DM object have a matrix function
1739 
1740     Not Collective
1741 
1742     Input Parameter:
1743 .   dm - the DM object to destroy
1744 
1745     Output Parameter:
1746 .   flg - PETSC_TRUE if function exists
1747 
1748     Level: developer
1749 
1750 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1751 
1752 @*/
1753 PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
1754 {
1755   PetscFunctionBegin;
1756   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1757   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
1758   PetscFunctionReturn(0);
1759 }
1760 
1761 #undef  __FUNCT__
1762 #define __FUNCT__ "DMSetVec"
1763 /*@C
1764     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
1765 
1766     Collective on DM
1767 
1768     Input Parameter:
1769 +   dm - the DM object
1770 -   x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems.
1771 
1772     Level: developer
1773 
1774 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1775          DMSetFunction(), DMSetJacobian(), DMSetVariableBounds()
1776 
1777 @*/
1778 PetscErrorCode  DMSetVec(DM dm,Vec x)
1779 {
1780   PetscErrorCode ierr;
1781   PetscFunctionBegin;
1782   if (x) {
1783     if (!dm->x) {
1784       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1785     }
1786     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1787   }
1788   else if(dm->x) {
1789     ierr = VecDestroy(&dm->x);  CHKERRQ(ierr);
1790   }
1791   PetscFunctionReturn(0);
1792 }
1793 
1794 
1795 #undef __FUNCT__
1796 #define __FUNCT__ "DMComputeFunction"
1797 /*@
1798     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
1799 
1800     Collective on DM
1801 
1802     Input Parameter:
1803 +   dm - the DM object to destroy
1804 .   x - the location where the function is evaluationed, may be ignored for linear problems
1805 -   b - the vector to hold the right hand side entries
1806 
1807     Level: developer
1808 
1809 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1810          DMSetJacobian()
1811 
1812 @*/
1813 PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
1814 {
1815   PetscErrorCode ierr;
1816   PetscFunctionBegin;
1817   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1818   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1819   PetscStackPush("DM user function");
1820   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1821   PetscStackPop;
1822   PetscFunctionReturn(0);
1823 }
1824 
1825 
1826 
1827 #undef __FUNCT__
1828 #define __FUNCT__ "DMComputeJacobian"
1829 /*@
1830     DMComputeJacobian - compute the matrix entries for the solver
1831 
1832     Collective on DM
1833 
1834     Input Parameter:
1835 +   dm - the DM object
1836 .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
1837 .   A - matrix that defines the operator for the linear solve
1838 -   B - the matrix used to construct the preconditioner
1839 
1840     Level: developer
1841 
1842 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1843          DMSetFunction()
1844 
1845 @*/
1846 PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1847 {
1848   PetscErrorCode ierr;
1849 
1850   PetscFunctionBegin;
1851   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1852   if (!dm->ops->jacobian) {
1853     ISColoring     coloring;
1854     MatFDColoring  fd;
1855     const MatType  mtype;
1856 
1857     ierr = PetscObjectGetType((PetscObject)B,&mtype);CHKERRQ(ierr);
1858     ierr = DMCreateColoring(dm,dm->coloringtype,mtype,&coloring);CHKERRQ(ierr);
1859     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1860     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1861     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1862     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1863     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
1864 
1865     dm->fd = fd;
1866     dm->ops->jacobian = DMComputeJacobianDefault;
1867 
1868     /* don't know why this is needed */
1869     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1870   }
1871   if (!x) x = dm->x;
1872   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1873 
1874   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1875   if (x) {
1876     if (!dm->x) {
1877       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1878     }
1879     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1880   }
1881   if (A != B) {
1882     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1883     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1884   }
1885   PetscFunctionReturn(0);
1886 }
1887 
1888 
1889 PetscFList DMList                       = PETSC_NULL;
1890 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1891 
1892 #undef __FUNCT__
1893 #define __FUNCT__ "DMSetType"
1894 /*@C
1895   DMSetType - Builds a DM, for a particular DM implementation.
1896 
1897   Collective on DM
1898 
1899   Input Parameters:
1900 + dm     - The DM object
1901 - method - The name of the DM type
1902 
1903   Options Database Key:
1904 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1905 
1906   Notes:
1907   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1908 
1909   Level: intermediate
1910 
1911 .keywords: DM, set, type
1912 .seealso: DMGetType(), DMCreate()
1913 @*/
1914 PetscErrorCode  DMSetType(DM dm, const DMType method)
1915 {
1916   PetscErrorCode (*r)(DM);
1917   PetscBool      match;
1918   PetscErrorCode ierr;
1919 
1920   PetscFunctionBegin;
1921   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1922   ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1923   if (match) PetscFunctionReturn(0);
1924 
1925   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1926   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1927   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1928 
1929   if (dm->ops->destroy) {
1930     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1931     dm->ops->destroy = PETSC_NULL;
1932   }
1933   ierr = (*r)(dm);CHKERRQ(ierr);
1934   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1935   PetscFunctionReturn(0);
1936 }
1937 
1938 #undef __FUNCT__
1939 #define __FUNCT__ "DMGetType"
1940 /*@C
1941   DMGetType - Gets the DM type name (as a string) from the DM.
1942 
1943   Not Collective
1944 
1945   Input Parameter:
1946 . dm  - The DM
1947 
1948   Output Parameter:
1949 . type - The DM type name
1950 
1951   Level: intermediate
1952 
1953 .keywords: DM, get, type, name
1954 .seealso: DMSetType(), DMCreate()
1955 @*/
1956 PetscErrorCode  DMGetType(DM dm, const DMType *type)
1957 {
1958   PetscErrorCode ierr;
1959 
1960   PetscFunctionBegin;
1961   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1962   PetscValidCharPointer(type,2);
1963   if (!DMRegisterAllCalled) {
1964     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1965   }
1966   *type = ((PetscObject)dm)->type_name;
1967   PetscFunctionReturn(0);
1968 }
1969 
1970 #undef __FUNCT__
1971 #define __FUNCT__ "DMConvert"
1972 /*@C
1973   DMConvert - Converts a DM to another DM, either of the same or different type.
1974 
1975   Collective on DM
1976 
1977   Input Parameters:
1978 + dm - the DM
1979 - newtype - new DM type (use "same" for the same type)
1980 
1981   Output Parameter:
1982 . M - pointer to new DM
1983 
1984   Notes:
1985   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
1986   the MPI communicator of the generated DM is always the same as the communicator
1987   of the input DM.
1988 
1989   Level: intermediate
1990 
1991 .seealso: DMCreate()
1992 @*/
1993 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
1994 {
1995   DM             B;
1996   char           convname[256];
1997   PetscBool      sametype, issame;
1998   PetscErrorCode ierr;
1999 
2000   PetscFunctionBegin;
2001   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2002   PetscValidType(dm,1);
2003   PetscValidPointer(M,3);
2004   ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
2005   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
2006   {
2007     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
2008 
2009     /*
2010        Order of precedence:
2011        1) See if a specialized converter is known to the current DM.
2012        2) See if a specialized converter is known to the desired DM class.
2013        3) See if a good general converter is registered for the desired class
2014        4) See if a good general converter is known for the current matrix.
2015        5) Use a really basic converter.
2016     */
2017 
2018     /* 1) See if a specialized converter is known to the current DM and the desired class */
2019     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
2020     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
2021     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
2022     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2023     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2024     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
2025     if (conv) goto foundconv;
2026 
2027     /* 2)  See if a specialized converter is known to the desired DM class. */
2028     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
2029     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
2030     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
2031     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
2032     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
2033     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2034     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2035     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
2036     if (conv) {
2037       ierr = DMDestroy(&B);CHKERRQ(ierr);
2038       goto foundconv;
2039     }
2040 
2041 #if 0
2042     /* 3) See if a good general converter is registered for the desired class */
2043     conv = B->ops->convertfrom;
2044     ierr = DMDestroy(&B);CHKERRQ(ierr);
2045     if (conv) goto foundconv;
2046 
2047     /* 4) See if a good general converter is known for the current matrix */
2048     if (dm->ops->convert) {
2049       conv = dm->ops->convert;
2050     }
2051     if (conv) goto foundconv;
2052 #endif
2053 
2054     /* 5) Use a really basic converter. */
2055     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
2056 
2057     foundconv:
2058     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2059     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
2060     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2061   }
2062   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
2063   PetscFunctionReturn(0);
2064 }
2065 
2066 /*--------------------------------------------------------------------------------------------------------------------*/
2067 
2068 #undef __FUNCT__
2069 #define __FUNCT__ "DMRegister"
2070 /*@C
2071   DMRegister - See DMRegisterDynamic()
2072 
2073   Level: advanced
2074 @*/
2075 PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
2076 {
2077   char fullname[PETSC_MAX_PATH_LEN];
2078   PetscErrorCode ierr;
2079 
2080   PetscFunctionBegin;
2081   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
2082   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
2083   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
2084   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
2085   PetscFunctionReturn(0);
2086 }
2087 
2088 
2089 /*--------------------------------------------------------------------------------------------------------------------*/
2090 #undef __FUNCT__
2091 #define __FUNCT__ "DMRegisterDestroy"
2092 /*@C
2093    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
2094 
2095    Not Collective
2096 
2097    Level: advanced
2098 
2099 .keywords: DM, register, destroy
2100 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
2101 @*/
2102 PetscErrorCode  DMRegisterDestroy(void)
2103 {
2104   PetscErrorCode ierr;
2105 
2106   PetscFunctionBegin;
2107   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
2108   DMRegisterAllCalled = PETSC_FALSE;
2109   PetscFunctionReturn(0);
2110 }
2111 
2112 #if defined(PETSC_HAVE_MATLAB_ENGINE)
2113 #include <mex.h>
2114 
2115 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
2116 
2117 #undef __FUNCT__
2118 #define __FUNCT__ "DMComputeFunction_Matlab"
2119 /*
2120    DMComputeFunction_Matlab - Calls the function that has been set with
2121                          DMSetFunctionMatlab().
2122 
2123    For linear problems x is null
2124 
2125 .seealso: DMSetFunction(), DMGetFunction()
2126 */
2127 PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
2128 {
2129   PetscErrorCode    ierr;
2130   DMMatlabContext   *sctx;
2131   int               nlhs = 1,nrhs = 4;
2132   mxArray	    *plhs[1],*prhs[4];
2133   long long int     lx = 0,ly = 0,ls = 0;
2134 
2135   PetscFunctionBegin;
2136   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2137   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
2138   PetscCheckSameComm(dm,1,y,3);
2139 
2140   /* call Matlab function in ctx with arguments x and y */
2141   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2142   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
2143   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2144   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
2145   prhs[0] =  mxCreateDoubleScalar((double)ls);
2146   prhs[1] =  mxCreateDoubleScalar((double)lx);
2147   prhs[2] =  mxCreateDoubleScalar((double)ly);
2148   prhs[3] =  mxCreateString(sctx->funcname);
2149   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
2150   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2151   mxDestroyArray(prhs[0]);
2152   mxDestroyArray(prhs[1]);
2153   mxDestroyArray(prhs[2]);
2154   mxDestroyArray(prhs[3]);
2155   mxDestroyArray(plhs[0]);
2156   PetscFunctionReturn(0);
2157 }
2158 
2159 
2160 #undef __FUNCT__
2161 #define __FUNCT__ "DMSetFunctionMatlab"
2162 /*
2163    DMSetFunctionMatlab - Sets the function evaluation routine
2164 
2165 */
2166 PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
2167 {
2168   PetscErrorCode    ierr;
2169   DMMatlabContext   *sctx;
2170 
2171   PetscFunctionBegin;
2172   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2173   /* currently sctx is memory bleed */
2174   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2175   if (!sctx) {
2176     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
2177   }
2178   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2179   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
2180   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
2181   PetscFunctionReturn(0);
2182 }
2183 
2184 #undef __FUNCT__
2185 #define __FUNCT__ "DMComputeJacobian_Matlab"
2186 /*
2187    DMComputeJacobian_Matlab - Calls the function that has been set with
2188                          DMSetJacobianMatlab().
2189 
2190    For linear problems x is null
2191 
2192 .seealso: DMSetFunction(), DMGetFunction()
2193 */
2194 PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
2195 {
2196   PetscErrorCode    ierr;
2197   DMMatlabContext   *sctx;
2198   int               nlhs = 2,nrhs = 5;
2199   mxArray	    *plhs[2],*prhs[5];
2200   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
2201 
2202   PetscFunctionBegin;
2203   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2204   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
2205 
2206   /* call MATLAB function in ctx with arguments x, A, and B */
2207   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2208   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
2209   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2210   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
2211   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
2212   prhs[0] =  mxCreateDoubleScalar((double)ls);
2213   prhs[1] =  mxCreateDoubleScalar((double)lx);
2214   prhs[2] =  mxCreateDoubleScalar((double)lA);
2215   prhs[3] =  mxCreateDoubleScalar((double)lB);
2216   prhs[4] =  mxCreateString(sctx->jacname);
2217   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
2218   *str    =  (MatStructure) mxGetScalar(plhs[0]);
2219   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2220   mxDestroyArray(prhs[0]);
2221   mxDestroyArray(prhs[1]);
2222   mxDestroyArray(prhs[2]);
2223   mxDestroyArray(prhs[3]);
2224   mxDestroyArray(prhs[4]);
2225   mxDestroyArray(plhs[0]);
2226   mxDestroyArray(plhs[1]);
2227   PetscFunctionReturn(0);
2228 }
2229 
2230 
2231 #undef __FUNCT__
2232 #define __FUNCT__ "DMSetJacobianMatlab"
2233 /*
2234    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
2235 
2236 */
2237 PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
2238 {
2239   PetscErrorCode    ierr;
2240   DMMatlabContext   *sctx;
2241 
2242   PetscFunctionBegin;
2243   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2244   /* currently sctx is memory bleed */
2245   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2246   if (!sctx) {
2247     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
2248   }
2249   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
2250   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
2251   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
2252   PetscFunctionReturn(0);
2253 }
2254 #endif
2255 
2256 #undef __FUNCT__
2257 #define __FUNCT__ "DMLoad"
2258 /*@C
2259   DMLoad - Loads a DM that has been stored in binary or HDF5 format
2260   with DMView().
2261 
2262   Collective on PetscViewer
2263 
2264   Input Parameters:
2265 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2266            some related function before a call to DMLoad().
2267 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2268            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2269 
2270    Level: intermediate
2271 
2272   Notes:
2273   Defaults to the DM DA.
2274 
2275   Notes for advanced users:
2276   Most users should not need to know the details of the binary storage
2277   format, since DMLoad() and DMView() completely hide these details.
2278   But for anyone who's interested, the standard binary matrix storage
2279   format is
2280 .vb
2281      has not yet been determined
2282 .ve
2283 
2284    In addition, PETSc automatically does the byte swapping for
2285 machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
2286 linux, Windows and the paragon; thus if you write your own binary
2287 read/write routines you have to swap the bytes; see PetscBinaryRead()
2288 and PetscBinaryWrite() to see how this may be done.
2289 
2290   Concepts: vector^loading from file
2291 
2292 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2293 @*/
2294 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2295 {
2296   PetscErrorCode ierr;
2297 
2298   PetscFunctionBegin;
2299   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2300   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2301 
2302   if (!((PetscObject)newdm)->type_name) {
2303     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
2304   }
2305   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2306   PetscFunctionReturn(0);
2307 }
2308 
2309 /******************************** FEM Support **********************************/
2310 
2311 #undef __FUNCT__
2312 #define __FUNCT__ "DMPrintCellVector"
2313 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) {
2314   PetscInt       f;
2315   PetscErrorCode ierr;
2316 
2317   PetscFunctionBegin;
2318   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2319   for(f = 0; f < len; ++f) {
2320     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
2321   }
2322   PetscFunctionReturn(0);
2323 }
2324 
2325 #undef __FUNCT__
2326 #define __FUNCT__ "DMPrintCellMatrix"
2327 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) {
2328   PetscInt       f, g;
2329   PetscErrorCode ierr;
2330 
2331   PetscFunctionBegin;
2332   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2333   for(f = 0; f < rows; ++f) {
2334     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
2335     for(g = 0; g < cols; ++g) {
2336       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
2337     }
2338     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
2339   }
2340   PetscFunctionReturn(0);
2341 }
2342 
2343 #undef __FUNCT__
2344 #define __FUNCT__ "DMGetLocalFunction"
2345 /*@C
2346   DMGetLocalFunction - Get the local residual function from this DM
2347 
2348   Not collective
2349 
2350   Input Parameter:
2351 . dm - The DM
2352 
2353   Output Parameter:
2354 . lf - The local residual function
2355 
2356    Calling sequence of lf:
2357 $    lf (SNES snes, Vec x, Vec f, void *ctx);
2358 
2359 +  snes - the SNES context
2360 .  x - local vector with the state at which to evaluate residual
2361 .  f - local vector to put residual in
2362 -  ctx - optional user-defined function context
2363 
2364   Level: intermediate
2365 
2366 .seealso DMSetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian()
2367 @*/
2368 PetscErrorCode DMGetLocalFunction(DM dm, PetscErrorCode (**lf)(DM, Vec, Vec, void *))
2369 {
2370   PetscFunctionBegin;
2371   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2372   if (lf) *lf = dm->lf;
2373   PetscFunctionReturn(0);
2374 }
2375 
2376 #undef __FUNCT__
2377 #define __FUNCT__ "DMSetLocalFunction"
2378 /*@C
2379   DMSetLocalFunction - Set the local residual function from this DM
2380 
2381   Not collective
2382 
2383   Input Parameters:
2384 + dm - The DM
2385 - lf - The local residual function
2386 
2387    Calling sequence of lf:
2388 $    lf (SNES snes, Vec x, Vec f, void *ctx);
2389 
2390 +  snes - the SNES context
2391 .  x - local vector with the state at which to evaluate residual
2392 .  f - local vector to put residual in
2393 -  ctx - optional user-defined function context
2394 
2395   Level: intermediate
2396 
2397 .seealso DMGetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian()
2398 @*/
2399 PetscErrorCode DMSetLocalFunction(DM dm, PetscErrorCode (*lf)(DM, Vec, Vec, void *))
2400 {
2401   PetscFunctionBegin;
2402   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2403   dm->lf = lf;
2404   PetscFunctionReturn(0);
2405 }
2406 
2407 #undef __FUNCT__
2408 #define __FUNCT__ "DMGetLocalJacobian"
2409 /*@C
2410   DMGetLocalJacobian - Get the local Jacobian function from this DM
2411 
2412   Not collective
2413 
2414   Input Parameter:
2415 . dm - The DM
2416 
2417   Output Parameter:
2418 . lj - The local Jacobian function
2419 
2420    Calling sequence of lj:
2421 $    lj (SNES snes, Vec x, Mat J, Mat M, void *ctx);
2422 
2423 +  snes - the SNES context
2424 .  x - local vector with the state at which to evaluate residual
2425 .  J - matrix to put Jacobian in
2426 .  M - matrix to use for defining Jacobian preconditioner
2427 -  ctx - optional user-defined function context
2428 
2429   Level: intermediate
2430 
2431 .seealso DMSetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction()
2432 @*/
2433 PetscErrorCode DMGetLocalJacobian(DM dm, PetscErrorCode (**lj)(DM, Vec, Mat, Mat, void *))
2434 {
2435   PetscFunctionBegin;
2436   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2437   if (lj) *lj = dm->lj;
2438   PetscFunctionReturn(0);
2439 }
2440 
2441 #undef __FUNCT__
2442 #define __FUNCT__ "DMSetLocalJacobian"
2443 /*@C
2444   DMSetLocalJacobian - Set the local Jacobian function from this DM
2445 
2446   Not collective
2447 
2448   Input Parameters:
2449 + dm - The DM
2450 - lj - The local Jacobian function
2451 
2452    Calling sequence of lj:
2453 $    lj (SNES snes, Vec x, Mat J, Mat M, void *ctx);
2454 
2455 +  snes - the SNES context
2456 .  x - local vector with the state at which to evaluate residual
2457 .  J - matrix to put Jacobian in
2458 .  M - matrix to use for defining Jacobian preconditioner
2459 -  ctx - optional user-defined function context
2460 
2461   Level: intermediate
2462 
2463 .seealso DMGetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction()
2464 @*/
2465 PetscErrorCode DMSetLocalJacobian(DM dm, PetscErrorCode (*lj)(DM, Vec, Mat,  Mat, void *))
2466 {
2467   PetscFunctionBegin;
2468   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2469   dm->lj = lj;
2470   PetscFunctionReturn(0);
2471 }
2472 
2473 #undef __FUNCT__
2474 #define __FUNCT__ "DMGetDefaultSection"
2475 /*@
2476   DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM.
2477 
2478   Input Parameter:
2479 . dm - The DM
2480 
2481   Output Parameter:
2482 . section - The PetscSection
2483 
2484   Level: intermediate
2485 
2486   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
2487 
2488 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
2489 @*/
2490 PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section) {
2491   PetscFunctionBegin;
2492   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2493   PetscValidPointer(section, 2);
2494   *section = dm->defaultSection;
2495   PetscFunctionReturn(0);
2496 }
2497 
2498 #undef __FUNCT__
2499 #define __FUNCT__ "DMSetDefaultSection"
2500 /*@
2501   DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM.
2502 
2503   Input Parameters:
2504 + dm - The DM
2505 - section - The PetscSection
2506 
2507   Level: intermediate
2508 
2509   Note: Any existing Section will be destroyed
2510 
2511 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection()
2512 @*/
2513 PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section) {
2514   PetscErrorCode ierr;
2515 
2516   PetscFunctionBegin;
2517   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2518   ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr);
2519   ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr);
2520   dm->defaultSection = section;
2521   PetscFunctionReturn(0);
2522 }
2523 
2524 #undef __FUNCT__
2525 #define __FUNCT__ "DMGetDefaultGlobalSection"
2526 /*@
2527   DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM.
2528 
2529   Input Parameter:
2530 . dm - The DM
2531 
2532   Output Parameter:
2533 . section - The PetscSection
2534 
2535   Level: intermediate
2536 
2537   Note: This gets a borrowed reference, so the user should not destroy this PetscSection.
2538 
2539 .seealso: DMSetDefaultSection(), DMGetDefaultSection()
2540 @*/
2541 PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section) {
2542   PetscErrorCode ierr;
2543 
2544   PetscFunctionBegin;
2545   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2546   PetscValidPointer(section, 2);
2547   if (!dm->defaultGlobalSection) {
2548     ierr = PetscSectionCreateGlobalSection(dm->defaultSection, dm->sf, &dm->defaultGlobalSection);CHKERRQ(ierr);
2549   }
2550   *section = dm->defaultGlobalSection;
2551   PetscFunctionReturn(0);
2552 }
2553 
2554 #undef __FUNCT__
2555 #define __FUNCT__ "DMGetDefaultSF"
2556 /*@
2557   DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set,
2558   it is created from the default PetscSection layouts in the DM.
2559 
2560   Input Parameter:
2561 . dm - The DM
2562 
2563   Output Parameter:
2564 . sf - The PetscSF
2565 
2566   Level: intermediate
2567 
2568   Note: This gets a borrowed reference, so the user should not destroy this PetscSF.
2569 
2570 .seealso: DMSetDefaultSF(), DMCreateDefaultSF()
2571 @*/
2572 PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf) {
2573   PetscInt       nroots;
2574   PetscErrorCode ierr;
2575 
2576   PetscFunctionBegin;
2577   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2578   PetscValidPointer(sf, 2);
2579   ierr = PetscSFGetGraph(dm->defaultSF, &nroots, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
2580   if (nroots < 0) {
2581     PetscSection section, gSection;
2582 
2583     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
2584     ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
2585     ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr);
2586   }
2587   *sf = dm->defaultSF;
2588   PetscFunctionReturn(0);
2589 }
2590 
2591 #undef __FUNCT__
2592 #define __FUNCT__ "DMSetDefaultSF"
2593 /*@
2594   DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM
2595 
2596   Input Parameters:
2597 + dm - The DM
2598 - sf - The PetscSF
2599 
2600   Level: intermediate
2601 
2602   Note: Any previous SF is destroyed
2603 
2604 .seealso: DMGetDefaultSF(), DMCreateDefaultSF()
2605 @*/
2606 PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf) {
2607   PetscErrorCode ierr;
2608 
2609   PetscFunctionBegin;
2610   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2611   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
2612   ierr = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr);
2613   dm->defaultSF = sf;
2614   PetscFunctionReturn(0);
2615 }
2616 
2617 #undef __FUNCT__
2618 #define __FUNCT__ "DMCreateDefaultSF"
2619 /*@C
2620   DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections
2621   describing the data layout.
2622 
2623   Input Parameters:
2624 + dm - The DM
2625 . localSection - PetscSection describing the local data layout
2626 - globalSection - PetscSection describing the global data layout
2627 
2628   Level: intermediate
2629 
2630 .seealso: DMGetDefaultSF(), DMSetDefaultSF()
2631 @*/
2632 PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection)
2633 {
2634   MPI_Comm        comm = ((PetscObject) dm)->comm;
2635   PetscLayout     layout;
2636   const PetscInt *ranges;
2637   PetscInt       *local;
2638   PetscSFNode    *remote;
2639   PetscInt        pStart, pEnd, p, nroots, nleaves, l;
2640   PetscMPIInt     size, rank;
2641   PetscErrorCode  ierr;
2642 
2643   PetscFunctionBegin;
2644   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2645   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
2646   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2647   ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr);
2648   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr);
2649   ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
2650   ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
2651   ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr);
2652   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
2653   ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
2654   for(p = pStart, nleaves = 0; p < pEnd; ++p) {
2655     PetscInt dof, cdof;
2656 
2657     ierr = PetscSectionGetDof(globalSection, p, &dof);CHKERRQ(ierr);
2658     ierr = PetscSectionGetConstraintDof(globalSection, p, &cdof);CHKERRQ(ierr);
2659     nleaves += dof < 0 ? -(dof+1)-cdof : dof-cdof;
2660   }
2661   ierr = PetscMalloc(nleaves * sizeof(PetscInt), &local);CHKERRQ(ierr);
2662   ierr = PetscMalloc(nleaves * sizeof(PetscSFNode), &remote);CHKERRQ(ierr);
2663   for(p = pStart, l = 0; p < pEnd; ++p) {
2664     PetscInt *cind;
2665     PetscInt  dof, gdof, cdof, dim, off, goff, d, c;
2666 
2667     ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr);
2668     ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr);
2669     ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr);
2670     ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr);
2671     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
2672     ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
2673     dim  = dof-cdof;
2674     for(d = 0, c = 0; d < dof; ++d) {
2675       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
2676       local[l+d-c] = off+d;
2677     }
2678     if (gdof < 0) {
2679       for(d = 0; d < dim; ++d, ++l) {
2680         PetscInt offset = -(goff+1) + d, r;
2681 
2682         for(r = 0; r < size; ++r) {
2683           if ((offset >= ranges[r]) && (offset < ranges[r+1])) break;
2684         }
2685         remote[l].rank  = r;
2686         remote[l].index = offset - ranges[r];
2687       }
2688     } else {
2689       for(d = 0; d < dim; ++d, ++l) {
2690         remote[l].rank  = rank;
2691         remote[l].index = goff+d - ranges[rank];
2692       }
2693     }
2694   }
2695   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
2696   ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
2697   PetscFunctionReturn(0);
2698 }
2699