xref: /petsc/src/dm/interface/dm.c (revision c4721b0ef1399f9b4324b73a33db78f8d329fdf8)
1 
2 #include <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   v->workSize     = 0;
45   v->workArray    = PETSC_NULL;
46   v->ltogmap      = PETSC_NULL;
47   v->ltogmapb     = PETSC_NULL;
48   v->bs           = 1;
49   v->coloringtype = IS_COLORING_GLOBAL;
50 
51   *dm = v;
52   PetscFunctionReturn(0);
53 }
54 
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "DMSetVecType"
58 /*@C
59        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
60 
61    Logically Collective on DMDA
62 
63    Input Parameter:
64 +  da - initial distributed array
65 .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
66 
67    Options Database:
68 .   -dm_vec_type ctype
69 
70    Level: intermediate
71 
72 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
73 @*/
74 PetscErrorCode  DMSetVecType(DM da,const VecType ctype)
75 {
76   PetscErrorCode ierr;
77 
78   PetscFunctionBegin;
79   PetscValidHeaderSpecific(da,DM_CLASSID,1);
80   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
81   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
82   PetscFunctionReturn(0);
83 }
84 
85 #undef __FUNCT__
86 #define __FUNCT__ "DMSetOptionsPrefix"
87 /*@C
88    DMSetOptionsPrefix - Sets the prefix used for searching for all
89    DMDA options in the database.
90 
91    Logically Collective on DMDA
92 
93    Input Parameter:
94 +  da - the DMDA context
95 -  prefix - the prefix to prepend to all option names
96 
97    Notes:
98    A hyphen (-) must NOT be given at the beginning of the prefix name.
99    The first character of all runtime options is AUTOMATICALLY the hyphen.
100 
101    Level: advanced
102 
103 .keywords: DMDA, set, options, prefix, database
104 
105 .seealso: DMSetFromOptions()
106 @*/
107 PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
108 {
109   PetscErrorCode ierr;
110 
111   PetscFunctionBegin;
112   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
113   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
114   PetscFunctionReturn(0);
115 }
116 
117 #undef __FUNCT__
118 #define __FUNCT__ "DMDestroy"
119 /*@
120     DMDestroy - Destroys a vector packer or DMDA.
121 
122     Collective on DM
123 
124     Input Parameter:
125 .   dm - the DM object to destroy
126 
127     Level: developer
128 
129 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
130 
131 @*/
132 PetscErrorCode  DMDestroy(DM *dm)
133 {
134   PetscInt       i, cnt = 0;
135   PetscErrorCode ierr;
136 
137   PetscFunctionBegin;
138   if (!*dm) PetscFunctionReturn(0);
139   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
140 
141   /* count all the circular references of DM and its contained Vecs */
142   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
143     if ((*dm)->localin[i])  {cnt++;}
144     if ((*dm)->globalin[i]) {cnt++;}
145   }
146   if ((*dm)->x) {
147     PetscObject obj;
148     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
149     if (obj == (PetscObject)*dm) cnt++;
150   }
151 
152   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
153   /*
154      Need this test because the dm references the vectors that
155      reference the dm, so destroying the dm calls destroy on the
156      vectors that cause another destroy on the dm
157   */
158   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
159   ((PetscObject) (*dm))->refct = 0;
160   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
161     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
162     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
163   }
164 
165   if ((*dm)->ctx && (*dm)->ctxdestroy) {
166     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
167   }
168   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
169   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
170   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
171   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
172   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
173   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
174   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
175   ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr);
176   /* if memory was published with AMS then destroy it */
177   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
178 
179   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
180   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
181   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
182   PetscFunctionReturn(0);
183 }
184 
185 #undef __FUNCT__
186 #define __FUNCT__ "DMSetUp"
187 /*@
188     DMSetUp - sets up the data structures inside a DM object
189 
190     Collective on DM
191 
192     Input Parameter:
193 .   dm - the DM object to setup
194 
195     Level: developer
196 
197 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
198 
199 @*/
200 PetscErrorCode  DMSetUp(DM dm)
201 {
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
206   if (dm->setupcalled) PetscFunctionReturn(0);
207   if (dm->ops->setup) {
208     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
209   }
210   dm->setupcalled = PETSC_TRUE;
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "DMSetFromOptions"
216 /*@
217     DMSetFromOptions - sets parameters in a DM from the options database
218 
219     Collective on DM
220 
221     Input Parameter:
222 .   dm - the DM object to set options for
223 
224     Options Database:
225 +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
226 .   -dm_vec_type <type>  type of vector to create inside DM
227 .   -dm_mat_type <type>  type of matrix to create inside DM
228 -   -dm_coloring_type <global or ghosted>
229 
230     Level: developer
231 
232 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
233 
234 @*/
235 PetscErrorCode  DMSetFromOptions(DM dm)
236 {
237   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg;
238   PetscErrorCode ierr;
239   char           typeName[256] = MATAIJ;
240 
241   PetscFunctionBegin;
242   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
243   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
244     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
245     ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr);
246     ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr);
247     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);
248     ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
249     if (flg) {
250       ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
251     }
252     ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr);
253     if (flg) {
254       ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
255       ierr = PetscStrallocpy(typeName,&dm->mattype);CHKERRQ(ierr);
256     }
257     ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr);
258     if (dm->ops->setfromoptions) {
259       ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
260     }
261     /* process any options handlers added with PetscObjectAddOptionsHandler() */
262     ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
263   ierr = PetscOptionsEnd();CHKERRQ(ierr);
264   if (flg1) {
265     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
266   }
267   if (flg2) {
268     PetscViewer viewer;
269 
270     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
271     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
272     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
273     ierr = DMView(dm, viewer);CHKERRQ(ierr);
274     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
275   }
276   if (flg3) {
277     PetscViewer viewer;
278 
279     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
280     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
281     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
282     ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr);
283     ierr = DMView(dm, viewer);CHKERRQ(ierr);
284     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
285   }
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "DMView"
291 /*@C
292     DMView - Views a vector packer or DMDA.
293 
294     Collective on DM
295 
296     Input Parameter:
297 +   dm - the DM object to view
298 -   v - the viewer
299 
300     Level: developer
301 
302 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
303 
304 @*/
305 PetscErrorCode  DMView(DM dm,PetscViewer v)
306 {
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
311  if (!v) {
312     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
313   }
314   if (dm->ops->view) {
315     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
316   }
317   PetscFunctionReturn(0);
318 }
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "DMCreateGlobalVector"
322 /*@
323     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
324 
325     Collective on DM
326 
327     Input Parameter:
328 .   dm - the DM object
329 
330     Output Parameter:
331 .   vec - the global vector
332 
333     Level: beginner
334 
335 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
336 
337 @*/
338 PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
339 {
340   PetscErrorCode ierr;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
344   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
345   PetscFunctionReturn(0);
346 }
347 
348 #undef __FUNCT__
349 #define __FUNCT__ "DMCreateLocalVector"
350 /*@
351     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
352 
353     Not Collective
354 
355     Input Parameter:
356 .   dm - the DM object
357 
358     Output Parameter:
359 .   vec - the local vector
360 
361     Level: beginner
362 
363 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
364 
365 @*/
366 PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
367 {
368   PetscErrorCode ierr;
369 
370   PetscFunctionBegin;
371   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
372   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
373   PetscFunctionReturn(0);
374 }
375 
376 #undef __FUNCT__
377 #define __FUNCT__ "DMGetLocalToGlobalMapping"
378 /*@
379    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
380 
381    Collective on DM
382 
383    Input Parameter:
384 .  dm - the DM that provides the mapping
385 
386    Output Parameter:
387 .  ltog - the mapping
388 
389    Level: intermediate
390 
391    Notes:
392    This mapping can then be used by VecSetLocalToGlobalMapping() or
393    MatSetLocalToGlobalMapping().
394 
395 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
396 @*/
397 PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
398 {
399   PetscErrorCode ierr;
400 
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
403   PetscValidPointer(ltog,2);
404   if (!dm->ltogmap) {
405     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
406     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
407   }
408   *ltog = dm->ltogmap;
409   PetscFunctionReturn(0);
410 }
411 
412 #undef __FUNCT__
413 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
414 /*@
415    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
416 
417    Collective on DM
418 
419    Input Parameter:
420 .  da - the distributed array that provides the mapping
421 
422    Output Parameter:
423 .  ltog - the block mapping
424 
425    Level: intermediate
426 
427    Notes:
428    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
429    MatSetLocalToGlobalMappingBlock().
430 
431 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
432 @*/
433 PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
434 {
435   PetscErrorCode ierr;
436 
437   PetscFunctionBegin;
438   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
439   PetscValidPointer(ltog,2);
440   if (!dm->ltogmapb) {
441     PetscInt bs;
442     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
443     if (bs > 1) {
444       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
445       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
446     } else {
447       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
448       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
449     }
450   }
451   *ltog = dm->ltogmapb;
452   PetscFunctionReturn(0);
453 }
454 
455 #undef __FUNCT__
456 #define __FUNCT__ "DMGetBlockSize"
457 /*@
458    DMGetBlockSize - Gets the inherent block size associated with a DM
459 
460    Not Collective
461 
462    Input Parameter:
463 .  dm - the DM with block structure
464 
465    Output Parameter:
466 .  bs - the block size, 1 implies no exploitable block structure
467 
468    Level: intermediate
469 
470 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
471 @*/
472 PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
473 {
474   PetscFunctionBegin;
475   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
476   PetscValidPointer(bs,2);
477   if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet");
478   *bs = dm->bs;
479   PetscFunctionReturn(0);
480 }
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "DMCreateInterpolation"
484 /*@
485     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
486 
487     Collective on DM
488 
489     Input Parameter:
490 +   dm1 - the DM object
491 -   dm2 - the second, finer DM object
492 
493     Output Parameter:
494 +  mat - the interpolation
495 -  vec - the scaling (optional)
496 
497     Level: developer
498 
499     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
500         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
501 
502         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
503         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
504 
505 
506 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
507 
508 @*/
509 PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
510 {
511   PetscErrorCode ierr;
512 
513   PetscFunctionBegin;
514   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
515   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
516   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "DMCreateInjection"
522 /*@
523     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
524 
525     Collective on DM
526 
527     Input Parameter:
528 +   dm1 - the DM object
529 -   dm2 - the second, finer DM object
530 
531     Output Parameter:
532 .   ctx - the injection
533 
534     Level: developer
535 
536    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
537         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
538 
539 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
540 
541 @*/
542 PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
543 {
544   PetscErrorCode ierr;
545 
546   PetscFunctionBegin;
547   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
548   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
549   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
550   PetscFunctionReturn(0);
551 }
552 
553 #undef __FUNCT__
554 #define __FUNCT__ "DMCreateColoring"
555 /*@C
556     DMCreateColoring - Gets coloring for a DMDA or DMComposite
557 
558     Collective on DM
559 
560     Input Parameter:
561 +   dm - the DM object
562 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
563 -   matype - either MATAIJ or MATBAIJ
564 
565     Output Parameter:
566 .   coloring - the coloring
567 
568     Level: developer
569 
570 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
571 
572 @*/
573 PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
574 {
575   PetscErrorCode ierr;
576 
577   PetscFunctionBegin;
578   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
579   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
580   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
581   PetscFunctionReturn(0);
582 }
583 
584 #undef __FUNCT__
585 #define __FUNCT__ "DMCreateMatrix"
586 /*@C
587     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
588 
589     Collective on DM
590 
591     Input Parameter:
592 +   dm - the DM object
593 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
594             any type which inherits from one of these (such as MATAIJ)
595 
596     Output Parameter:
597 .   mat - the empty Jacobian
598 
599     Level: beginner
600 
601     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
602        do not need to do it yourself.
603 
604        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
605        the nonzero pattern call DMDASetMatPreallocateOnly()
606 
607        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
608        internally by PETSc.
609 
610        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
611        the indices for the global numbering for DMDAs which is complicated.
612 
613 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
614 
615 @*/
616 PetscErrorCode  DMCreateMatrix(DM dm,const MatType mtype,Mat *mat)
617 {
618   PetscErrorCode ierr;
619 
620   PetscFunctionBegin;
621   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
622 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
623   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
624 #endif
625   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
626   PetscValidPointer(mat,3);
627   if (dm->mattype) {
628     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
629   } else {
630     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
631   }
632   PetscFunctionReturn(0);
633 }
634 
635 #undef __FUNCT__
636 #define __FUNCT__ "DMSetMatrixPreallocateOnly"
637 /*@
638   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
639     preallocated but the nonzero structure and zero values will not be set.
640 
641   Logically Collective on DMDA
642 
643   Input Parameter:
644 + dm - the DM
645 - only - PETSC_TRUE if only want preallocation
646 
647   Level: developer
648 .seealso DMCreateMatrix()
649 @*/
650 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
651 {
652   PetscFunctionBegin;
653   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
654   dm->prealloc_only = only;
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNCT__
659 #define __FUNCT__ "DMGetWorkArray"
660 /*@C
661   DMGetWorkArray - Gets a work array guaranteed to be at least the input size
662 
663   Not Collective
664 
665   Input Parameters:
666 + dm - the DM object
667 - size - The minium size
668 
669   Output Parameter:
670 . array - the work array
671 
672   Level: developer
673 
674 .seealso DMDestroy(), DMCreate()
675 @*/
676 PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array)
677 {
678   PetscErrorCode ierr;
679 
680   PetscFunctionBegin;
681   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
682   PetscValidPointer(array,3);
683   if (size > dm->workSize) {
684     dm->workSize = size;
685     ierr = PetscFree(dm->workArray);CHKERRQ(ierr);
686     ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr);
687   }
688   *array = dm->workArray;
689   PetscFunctionReturn(0);
690 }
691 
692 #undef __FUNCT__
693 #define __FUNCT__ "DMRefine"
694 /*@
695     DMRefine - Refines a DM object
696 
697     Collective on DM
698 
699     Input Parameter:
700 +   dm - the DM object
701 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
702 
703     Output Parameter:
704 .   dmf - the refined DM
705 
706     Level: developer
707 
708 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
709 
710 @*/
711 PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
712 {
713   PetscErrorCode ierr;
714 
715   PetscFunctionBegin;
716   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
717   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
718   (*dmf)->ops->initialguess = dm->ops->initialguess;
719   (*dmf)->ops->function     = dm->ops->function;
720   (*dmf)->ops->functionj    = dm->ops->functionj;
721   if (dm->ops->jacobian != DMComputeJacobianDefault) {
722     (*dmf)->ops->jacobian     = dm->ops->jacobian;
723   }
724   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
725   (*dmf)->ctx     = dm->ctx;
726   (*dmf)->levelup = dm->levelup + 1;
727   PetscFunctionReturn(0);
728 }
729 
730 #undef __FUNCT__
731 #define __FUNCT__ "DMGetRefineLevel"
732 /*@
733     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
734 
735     Not Collective
736 
737     Input Parameter:
738 .   dm - the DM object
739 
740     Output Parameter:
741 .   level - number of refinements
742 
743     Level: developer
744 
745 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
746 
747 @*/
748 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
749 {
750   PetscFunctionBegin;
751   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
752   *level = dm->levelup;
753   PetscFunctionReturn(0);
754 }
755 
756 #undef __FUNCT__
757 #define __FUNCT__ "DMGlobalToLocalBegin"
758 /*@
759     DMGlobalToLocalBegin - Begins updating local vectors from global vector
760 
761     Neighbor-wise Collective on DM
762 
763     Input Parameters:
764 +   dm - the DM object
765 .   g - the global vector
766 .   mode - INSERT_VALUES or ADD_VALUES
767 -   l - the local vector
768 
769 
770     Level: beginner
771 
772 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
773 
774 @*/
775 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
776 {
777   PetscErrorCode ierr;
778 
779   PetscFunctionBegin;
780   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
781   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
782   PetscFunctionReturn(0);
783 }
784 
785 #undef __FUNCT__
786 #define __FUNCT__ "DMGlobalToLocalEnd"
787 /*@
788     DMGlobalToLocalEnd - Ends updating local vectors from global vector
789 
790     Neighbor-wise Collective on DM
791 
792     Input Parameters:
793 +   dm - the DM object
794 .   g - the global vector
795 .   mode - INSERT_VALUES or ADD_VALUES
796 -   l - the local vector
797 
798 
799     Level: beginner
800 
801 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
802 
803 @*/
804 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
805 {
806   PetscErrorCode ierr;
807 
808   PetscFunctionBegin;
809   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
810   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
811   PetscFunctionReturn(0);
812 }
813 
814 #undef __FUNCT__
815 #define __FUNCT__ "DMLocalToGlobalBegin"
816 /*@
817     DMLocalToGlobalBegin - updates global vectors from local vectors
818 
819     Neighbor-wise Collective on DM
820 
821     Input Parameters:
822 +   dm - the DM object
823 .   l - the local vector
824 .   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
825            base point.
826 - - the global vector
827 
828     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
829            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
830            global array to the final global array with VecAXPY().
831 
832     Level: beginner
833 
834 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
835 
836 @*/
837 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
838 {
839   PetscErrorCode ierr;
840 
841   PetscFunctionBegin;
842   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
843   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
844   PetscFunctionReturn(0);
845 }
846 
847 #undef __FUNCT__
848 #define __FUNCT__ "DMLocalToGlobalEnd"
849 /*@
850     DMLocalToGlobalEnd - updates global vectors from local vectors
851 
852     Neighbor-wise Collective on DM
853 
854     Input Parameters:
855 +   dm - the DM object
856 .   l - the local vector
857 .   mode - INSERT_VALUES or ADD_VALUES
858 -   g - the global vector
859 
860 
861     Level: beginner
862 
863 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
864 
865 @*/
866 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
867 {
868   PetscErrorCode ierr;
869 
870   PetscFunctionBegin;
871   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
872   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "DMComputeJacobianDefault"
878 /*@
879     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
880 
881     Collective on DM
882 
883     Input Parameter:
884 +   dm - the DM object
885 .   x - location to compute Jacobian at; may be ignored for linear problems
886 .   A - matrix that defines the operator for the linear solve
887 -   B - the matrix used to construct the preconditioner
888 
889     Level: developer
890 
891 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
892          DMSetFunction()
893 
894 @*/
895 PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
896 {
897   PetscErrorCode ierr;
898 
899   PetscFunctionBegin;
900   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
901   *stflag = SAME_NONZERO_PATTERN;
902   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
903   if (A != B) {
904     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
905     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
906   }
907   PetscFunctionReturn(0);
908 }
909 
910 #undef __FUNCT__
911 #define __FUNCT__ "DMCoarsen"
912 /*@
913     DMCoarsen - Coarsens a DM object
914 
915     Collective on DM
916 
917     Input Parameter:
918 +   dm - the DM object
919 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
920 
921     Output Parameter:
922 .   dmc - the coarsened DM
923 
924     Level: developer
925 
926 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
927 
928 @*/
929 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
930 {
931   PetscErrorCode ierr;
932 
933   PetscFunctionBegin;
934   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
935   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
936   (*dmc)->ops->initialguess = dm->ops->initialguess;
937   (*dmc)->ops->function     = dm->ops->function;
938   (*dmc)->ops->functionj    = dm->ops->functionj;
939   if (dm->ops->jacobian != DMComputeJacobianDefault) {
940     (*dmc)->ops->jacobian     = dm->ops->jacobian;
941   }
942   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
943   (*dmc)->ctx       = dm->ctx;
944   (*dmc)->leveldown = dm->leveldown + 1;
945   PetscFunctionReturn(0);
946 }
947 
948 #undef __FUNCT__
949 #define __FUNCT__ "DMRefineHierarchy"
950 /*@C
951     DMRefineHierarchy - Refines a DM object, all levels at once
952 
953     Collective on DM
954 
955     Input Parameter:
956 +   dm - the DM object
957 -   nlevels - the number of levels of refinement
958 
959     Output Parameter:
960 .   dmf - the refined DM hierarchy
961 
962     Level: developer
963 
964 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
965 
966 @*/
967 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
968 {
969   PetscErrorCode ierr;
970 
971   PetscFunctionBegin;
972   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
973   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
974   if (nlevels == 0) PetscFunctionReturn(0);
975   if (dm->ops->refinehierarchy) {
976     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
977   } else if (dm->ops->refine) {
978     PetscInt i;
979 
980     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
981     for (i=1; i<nlevels; i++) {
982       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
983     }
984   } else {
985     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
986   }
987   PetscFunctionReturn(0);
988 }
989 
990 #undef __FUNCT__
991 #define __FUNCT__ "DMCoarsenHierarchy"
992 /*@C
993     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
994 
995     Collective on DM
996 
997     Input Parameter:
998 +   dm - the DM object
999 -   nlevels - the number of levels of coarsening
1000 
1001     Output Parameter:
1002 .   dmc - the coarsened DM hierarchy
1003 
1004     Level: developer
1005 
1006 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1007 
1008 @*/
1009 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
1010 {
1011   PetscErrorCode ierr;
1012 
1013   PetscFunctionBegin;
1014   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1015   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1016   if (nlevels == 0) PetscFunctionReturn(0);
1017   PetscValidPointer(dmc,3);
1018   if (dm->ops->coarsenhierarchy) {
1019     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
1020   } else if (dm->ops->coarsen) {
1021     PetscInt i;
1022 
1023     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
1024     for (i=1; i<nlevels; i++) {
1025       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
1026     }
1027   } else {
1028     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
1029   }
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 #undef __FUNCT__
1034 #define __FUNCT__ "DMCreateAggregates"
1035 /*@
1036    DMCreateAggregates - Gets the aggregates that map between
1037    grids associated with two DMs.
1038 
1039    Collective on DM
1040 
1041    Input Parameters:
1042 +  dmc - the coarse grid DM
1043 -  dmf - the fine grid DM
1044 
1045    Output Parameters:
1046 .  rest - the restriction matrix (transpose of the projection matrix)
1047 
1048    Level: intermediate
1049 
1050 .keywords: interpolation, restriction, multigrid
1051 
1052 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
1053 @*/
1054 PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
1055 {
1056   PetscErrorCode ierr;
1057 
1058   PetscFunctionBegin;
1059   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1060   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
1061   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 #undef __FUNCT__
1066 #define __FUNCT__ "DMSetApplicationContextDestroy"
1067 /*@C
1068     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
1069 
1070     Not Collective
1071 
1072     Input Parameters:
1073 +   dm - the DM object
1074 -   destroy - the destroy function
1075 
1076     Level: intermediate
1077 
1078 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1079 
1080 @*/
1081 PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
1082 {
1083   PetscFunctionBegin;
1084   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1085   dm->ctxdestroy = destroy;
1086   PetscFunctionReturn(0);
1087 }
1088 
1089 #undef __FUNCT__
1090 #define __FUNCT__ "DMSetApplicationContext"
1091 /*@
1092     DMSetApplicationContext - Set a user context into a DM object
1093 
1094     Not Collective
1095 
1096     Input Parameters:
1097 +   dm - the DM object
1098 -   ctx - the user context
1099 
1100     Level: intermediate
1101 
1102 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1103 
1104 @*/
1105 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
1106 {
1107   PetscFunctionBegin;
1108   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1109   dm->ctx = ctx;
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 #undef __FUNCT__
1114 #define __FUNCT__ "DMGetApplicationContext"
1115 /*@
1116     DMGetApplicationContext - Gets a user context from a DM object
1117 
1118     Not Collective
1119 
1120     Input Parameter:
1121 .   dm - the DM object
1122 
1123     Output Parameter:
1124 .   ctx - the user context
1125 
1126     Level: intermediate
1127 
1128 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1129 
1130 @*/
1131 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
1132 {
1133   PetscFunctionBegin;
1134   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1135   *(void**)ctx = dm->ctx;
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 #undef __FUNCT__
1140 #define __FUNCT__ "DMSetInitialGuess"
1141 /*@C
1142     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
1143 
1144     Logically Collective on DM
1145 
1146     Input Parameter:
1147 +   dm - the DM object to destroy
1148 -   f - the function to compute the initial guess
1149 
1150     Level: intermediate
1151 
1152 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1153 
1154 @*/
1155 PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
1156 {
1157   PetscFunctionBegin;
1158   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1159   dm->ops->initialguess = f;
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 #undef __FUNCT__
1164 #define __FUNCT__ "DMSetFunction"
1165 /*@C
1166     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
1167 
1168     Logically Collective on DM
1169 
1170     Input Parameter:
1171 +   dm - the DM object
1172 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
1173 
1174     Level: intermediate
1175 
1176     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
1177            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
1178 
1179 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1180          DMSetJacobian()
1181 
1182 @*/
1183 PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1184 {
1185   PetscFunctionBegin;
1186   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1187   dm->ops->function = f;
1188   if (f) {
1189     dm->ops->functionj = f;
1190   }
1191   PetscFunctionReturn(0);
1192 }
1193 
1194 #undef __FUNCT__
1195 #define __FUNCT__ "DMSetJacobian"
1196 /*@C
1197     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
1198 
1199     Logically Collective on DM
1200 
1201     Input Parameter:
1202 +   dm - the DM object to destroy
1203 -   f - the function to compute the matrix entries
1204 
1205     Level: intermediate
1206 
1207 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1208          DMSetFunction()
1209 
1210 @*/
1211 PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
1212 {
1213   PetscFunctionBegin;
1214   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1215   dm->ops->jacobian = f;
1216   PetscFunctionReturn(0);
1217 }
1218 
1219 #undef __FUNCT__
1220 #define __FUNCT__ "DMComputeInitialGuess"
1221 /*@
1222     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
1223 
1224     Collective on DM
1225 
1226     Input Parameter:
1227 +   dm - the DM object to destroy
1228 -   x - the vector to hold the initial guess values
1229 
1230     Level: developer
1231 
1232 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
1233 
1234 @*/
1235 PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
1236 {
1237   PetscErrorCode ierr;
1238   PetscFunctionBegin;
1239   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1240   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
1241   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
1242   PetscFunctionReturn(0);
1243 }
1244 
1245 #undef __FUNCT__
1246 #define __FUNCT__ "DMHasInitialGuess"
1247 /*@
1248     DMHasInitialGuess - does the DM object have an initial guess function
1249 
1250     Not Collective
1251 
1252     Input Parameter:
1253 .   dm - the DM object to destroy
1254 
1255     Output Parameter:
1256 .   flg - PETSC_TRUE if function exists
1257 
1258     Level: developer
1259 
1260 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1261 
1262 @*/
1263 PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
1264 {
1265   PetscFunctionBegin;
1266   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1267   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
1268   PetscFunctionReturn(0);
1269 }
1270 
1271 #undef __FUNCT__
1272 #define __FUNCT__ "DMHasFunction"
1273 /*@
1274     DMHasFunction - does the DM object have a function
1275 
1276     Not Collective
1277 
1278     Input Parameter:
1279 .   dm - the DM object to destroy
1280 
1281     Output Parameter:
1282 .   flg - PETSC_TRUE if function exists
1283 
1284     Level: developer
1285 
1286 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1287 
1288 @*/
1289 PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
1290 {
1291   PetscFunctionBegin;
1292   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1293   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
1294   PetscFunctionReturn(0);
1295 }
1296 
1297 #undef __FUNCT__
1298 #define __FUNCT__ "DMHasJacobian"
1299 /*@
1300     DMHasJacobian - does the DM object have a matrix function
1301 
1302     Not Collective
1303 
1304     Input Parameter:
1305 .   dm - the DM object to destroy
1306 
1307     Output Parameter:
1308 .   flg - PETSC_TRUE if function exists
1309 
1310     Level: developer
1311 
1312 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1313 
1314 @*/
1315 PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
1316 {
1317   PetscFunctionBegin;
1318   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1319   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
1320   PetscFunctionReturn(0);
1321 }
1322 
1323 #undef __FUNCT__
1324 #define __FUNCT__ "DMComputeFunction"
1325 /*@
1326     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
1327 
1328     Collective on DM
1329 
1330     Input Parameter:
1331 +   dm - the DM object to destroy
1332 .   x - the location where the function is evaluationed, may be ignored for linear problems
1333 -   b - the vector to hold the right hand side entries
1334 
1335     Level: developer
1336 
1337 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1338          DMSetJacobian()
1339 
1340 @*/
1341 PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
1342 {
1343   PetscErrorCode ierr;
1344   PetscFunctionBegin;
1345   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1346   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1347   PetscStackPush("DM user function");
1348   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1349   PetscStackPop;
1350   PetscFunctionReturn(0);
1351 }
1352 
1353 
1354 #undef __FUNCT__
1355 #define __FUNCT__ "DMComputeJacobian"
1356 /*@
1357     DMComputeJacobian - compute the matrix entries for the solver
1358 
1359     Collective on DM
1360 
1361     Input Parameter:
1362 +   dm - the DM object
1363 .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
1364 .   A - matrix that defines the operator for the linear solve
1365 -   B - the matrix used to construct the preconditioner
1366 
1367     Level: developer
1368 
1369 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1370          DMSetFunction()
1371 
1372 @*/
1373 PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1374 {
1375   PetscErrorCode ierr;
1376 
1377   PetscFunctionBegin;
1378   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1379   if (!dm->ops->jacobian) {
1380     ISColoring     coloring;
1381     MatFDColoring  fd;
1382 
1383     ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
1384     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1385     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1386     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1387     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1388     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
1389 
1390     dm->fd = fd;
1391     dm->ops->jacobian = DMComputeJacobianDefault;
1392 
1393     /* don't know why this is needed */
1394     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1395   }
1396   if (!x) x = dm->x;
1397   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1398 
1399   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1400   if (x) {
1401     if (!dm->x) {
1402       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1403     }
1404     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1405   }
1406   if (A != B) {
1407     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1408     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1409   }
1410   PetscFunctionReturn(0);
1411 }
1412 
1413 
1414 PetscFList DMList                       = PETSC_NULL;
1415 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1416 
1417 #undef __FUNCT__
1418 #define __FUNCT__ "DMSetType"
1419 /*@C
1420   DMSetType - Builds a DM, for a particular DM implementation.
1421 
1422   Collective on DM
1423 
1424   Input Parameters:
1425 + dm     - The DM object
1426 - method - The name of the DM type
1427 
1428   Options Database Key:
1429 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1430 
1431   Notes:
1432   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1433 
1434   Level: intermediate
1435 
1436 .keywords: DM, set, type
1437 .seealso: DMGetType(), DMCreate()
1438 @*/
1439 PetscErrorCode  DMSetType(DM dm, const DMType method)
1440 {
1441   PetscErrorCode (*r)(DM);
1442   PetscBool      match;
1443   PetscErrorCode ierr;
1444 
1445   PetscFunctionBegin;
1446   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1447   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1448   if (match) PetscFunctionReturn(0);
1449 
1450   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1451   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1452   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1453 
1454   if (dm->ops->destroy) {
1455     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1456   }
1457   ierr = (*r)(dm);CHKERRQ(ierr);
1458   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1459   PetscFunctionReturn(0);
1460 }
1461 
1462 #undef __FUNCT__
1463 #define __FUNCT__ "DMGetType"
1464 /*@C
1465   DMGetType - Gets the DM type name (as a string) from the DM.
1466 
1467   Not Collective
1468 
1469   Input Parameter:
1470 . dm  - The DM
1471 
1472   Output Parameter:
1473 . type - The DM type name
1474 
1475   Level: intermediate
1476 
1477 .keywords: DM, get, type, name
1478 .seealso: DMSetType(), DMCreate()
1479 @*/
1480 PetscErrorCode  DMGetType(DM dm, const DMType *type)
1481 {
1482   PetscErrorCode ierr;
1483 
1484   PetscFunctionBegin;
1485   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1486   PetscValidCharPointer(type,2);
1487   if (!DMRegisterAllCalled) {
1488     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1489   }
1490   *type = ((PetscObject)dm)->type_name;
1491   PetscFunctionReturn(0);
1492 }
1493 
1494 #undef __FUNCT__
1495 #define __FUNCT__ "DMConvert"
1496 /*@C
1497   DMConvert - Converts a DM to another DM, either of the same or different type.
1498 
1499   Collective on DM
1500 
1501   Input Parameters:
1502 + dm - the DM
1503 - newtype - new DM type (use "same" for the same type)
1504 
1505   Output Parameter:
1506 . M - pointer to new DM
1507 
1508   Notes:
1509   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
1510   the MPI communicator of the generated DM is always the same as the communicator
1511   of the input DM.
1512 
1513   Level: intermediate
1514 
1515 .seealso: DMCreate()
1516 @*/
1517 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
1518 {
1519   DM             B;
1520   char           convname[256];
1521   PetscBool      sametype, issame;
1522   PetscErrorCode ierr;
1523 
1524   PetscFunctionBegin;
1525   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1526   PetscValidType(dm,1);
1527   PetscValidPointer(M,3);
1528   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
1529   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
1530   {
1531     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
1532 
1533     /*
1534        Order of precedence:
1535        1) See if a specialized converter is known to the current DM.
1536        2) See if a specialized converter is known to the desired DM class.
1537        3) See if a good general converter is registered for the desired class
1538        4) See if a good general converter is known for the current matrix.
1539        5) Use a really basic converter.
1540     */
1541 
1542     /* 1) See if a specialized converter is known to the current DM and the desired class */
1543     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1544     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1545     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1546     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1547     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1548     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1549     if (conv) goto foundconv;
1550 
1551     /* 2)  See if a specialized converter is known to the desired DM class. */
1552     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
1553     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
1554     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1555     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1556     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1557     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1558     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1559     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1560     if (conv) {
1561       ierr = DMDestroy(&B);CHKERRQ(ierr);
1562       goto foundconv;
1563     }
1564 
1565 #if 0
1566     /* 3) See if a good general converter is registered for the desired class */
1567     conv = B->ops->convertfrom;
1568     ierr = DMDestroy(&B);CHKERRQ(ierr);
1569     if (conv) goto foundconv;
1570 
1571     /* 4) See if a good general converter is known for the current matrix */
1572     if (dm->ops->convert) {
1573       conv = dm->ops->convert;
1574     }
1575     if (conv) goto foundconv;
1576 #endif
1577 
1578     /* 5) Use a really basic converter. */
1579     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
1580 
1581     foundconv:
1582     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1583     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
1584     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1585   }
1586   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
1587   PetscFunctionReturn(0);
1588 }
1589 
1590 /*--------------------------------------------------------------------------------------------------------------------*/
1591 
1592 #undef __FUNCT__
1593 #define __FUNCT__ "DMRegister"
1594 /*@C
1595   DMRegister - See DMRegisterDynamic()
1596 
1597   Level: advanced
1598 @*/
1599 PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1600 {
1601   char fullname[PETSC_MAX_PATH_LEN];
1602   PetscErrorCode ierr;
1603 
1604   PetscFunctionBegin;
1605   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1606   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1607   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1608   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1609   PetscFunctionReturn(0);
1610 }
1611 
1612 
1613 /*--------------------------------------------------------------------------------------------------------------------*/
1614 #undef __FUNCT__
1615 #define __FUNCT__ "DMRegisterDestroy"
1616 /*@C
1617    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1618 
1619    Not Collective
1620 
1621    Level: advanced
1622 
1623 .keywords: DM, register, destroy
1624 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1625 @*/
1626 PetscErrorCode  DMRegisterDestroy(void)
1627 {
1628   PetscErrorCode ierr;
1629 
1630   PetscFunctionBegin;
1631   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1632   DMRegisterAllCalled = PETSC_FALSE;
1633   PetscFunctionReturn(0);
1634 }
1635 
1636 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1637 #include <mex.h>
1638 
1639 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
1640 
1641 #undef __FUNCT__
1642 #define __FUNCT__ "DMComputeFunction_Matlab"
1643 /*
1644    DMComputeFunction_Matlab - Calls the function that has been set with
1645                          DMSetFunctionMatlab().
1646 
1647    For linear problems x is null
1648 
1649 .seealso: DMSetFunction(), DMGetFunction()
1650 */
1651 PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
1652 {
1653   PetscErrorCode    ierr;
1654   DMMatlabContext   *sctx;
1655   int               nlhs = 1,nrhs = 4;
1656   mxArray	    *plhs[1],*prhs[4];
1657   long long int     lx = 0,ly = 0,ls = 0;
1658 
1659   PetscFunctionBegin;
1660   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1661   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1662   PetscCheckSameComm(dm,1,y,3);
1663 
1664   /* call Matlab function in ctx with arguments x and y */
1665   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1666   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1667   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1668   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
1669   prhs[0] =  mxCreateDoubleScalar((double)ls);
1670   prhs[1] =  mxCreateDoubleScalar((double)lx);
1671   prhs[2] =  mxCreateDoubleScalar((double)ly);
1672   prhs[3] =  mxCreateString(sctx->funcname);
1673   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
1674   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
1675   mxDestroyArray(prhs[0]);
1676   mxDestroyArray(prhs[1]);
1677   mxDestroyArray(prhs[2]);
1678   mxDestroyArray(prhs[3]);
1679   mxDestroyArray(plhs[0]);
1680   PetscFunctionReturn(0);
1681 }
1682 
1683 
1684 #undef __FUNCT__
1685 #define __FUNCT__ "DMSetFunctionMatlab"
1686 /*
1687    DMSetFunctionMatlab - Sets the function evaluation routine
1688 
1689 */
1690 PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
1691 {
1692   PetscErrorCode    ierr;
1693   DMMatlabContext   *sctx;
1694 
1695   PetscFunctionBegin;
1696   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1697   /* currently sctx is memory bleed */
1698   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1699   if (!sctx) {
1700     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1701   }
1702   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
1703   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1704   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
1705   PetscFunctionReturn(0);
1706 }
1707 
1708 #undef __FUNCT__
1709 #define __FUNCT__ "DMComputeJacobian_Matlab"
1710 /*
1711    DMComputeJacobian_Matlab - Calls the function that has been set with
1712                          DMSetJacobianMatlab().
1713 
1714    For linear problems x is null
1715 
1716 .seealso: DMSetFunction(), DMGetFunction()
1717 */
1718 PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
1719 {
1720   PetscErrorCode    ierr;
1721   DMMatlabContext   *sctx;
1722   int               nlhs = 2,nrhs = 5;
1723   mxArray	    *plhs[2],*prhs[5];
1724   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
1725 
1726   PetscFunctionBegin;
1727   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1728   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
1729 
1730   /* call MATLAB function in ctx with arguments x, A, and B */
1731   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1732   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1733   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1734   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
1735   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
1736   prhs[0] =  mxCreateDoubleScalar((double)ls);
1737   prhs[1] =  mxCreateDoubleScalar((double)lx);
1738   prhs[2] =  mxCreateDoubleScalar((double)lA);
1739   prhs[3] =  mxCreateDoubleScalar((double)lB);
1740   prhs[4] =  mxCreateString(sctx->jacname);
1741   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
1742   *str    =  (MatStructure) mxGetScalar(plhs[0]);
1743   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
1744   mxDestroyArray(prhs[0]);
1745   mxDestroyArray(prhs[1]);
1746   mxDestroyArray(prhs[2]);
1747   mxDestroyArray(prhs[3]);
1748   mxDestroyArray(prhs[4]);
1749   mxDestroyArray(plhs[0]);
1750   mxDestroyArray(plhs[1]);
1751   PetscFunctionReturn(0);
1752 }
1753 
1754 
1755 #undef __FUNCT__
1756 #define __FUNCT__ "DMSetJacobianMatlab"
1757 /*
1758    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
1759 
1760 */
1761 PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
1762 {
1763   PetscErrorCode    ierr;
1764   DMMatlabContext   *sctx;
1765 
1766   PetscFunctionBegin;
1767   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1768   /* currently sctx is memory bleed */
1769   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1770   if (!sctx) {
1771     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1772   }
1773   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
1774   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1775   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
1776   PetscFunctionReturn(0);
1777 }
1778 #endif
1779 
1780 #undef __FUNCT__
1781 #define __FUNCT__ "DMLoad"
1782 /*@C
1783   DMLoad - Loads a DM that has been stored in binary or HDF5 format
1784   with DMView().
1785 
1786   Collective on PetscViewer
1787 
1788   Input Parameters:
1789 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
1790            some related function before a call to DMLoad().
1791 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
1792            HDF5 file viewer, obtained from PetscViewerHDF5Open()
1793 
1794    Level: intermediate
1795 
1796   Notes:
1797   Defaults to the DM DA.
1798 
1799   Notes for advanced users:
1800   Most users should not need to know the details of the binary storage
1801   format, since DMLoad() and DMView() completely hide these details.
1802   But for anyone who's interested, the standard binary matrix storage
1803   format is
1804 .vb
1805      has not yet been determined
1806 .ve
1807 
1808    In addition, PETSc automatically does the byte swapping for
1809 machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
1810 linux, Windows and the paragon; thus if you write your own binary
1811 read/write routines you have to swap the bytes; see PetscBinaryRead()
1812 and PetscBinaryWrite() to see how this may be done.
1813 
1814   Concepts: vector^loading from file
1815 
1816 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
1817 @*/
1818 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
1819 {
1820   PetscErrorCode ierr;
1821 
1822   PetscFunctionBegin;
1823   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
1824   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1825 
1826   if (!((PetscObject)newdm)->type_name) {
1827     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
1828   }
1829   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1830   PetscFunctionReturn(0);
1831 }
1832 
1833