xref: /petsc/src/dm/interface/dm.c (revision bc0cc02b50436274feb0da0b98fa7f6ad991e737)
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 
693 #undef __FUNCT__
694 #define __FUNCT__ "DMRefine"
695 /*@
696     DMRefine - Refines a DM object
697 
698     Collective on DM
699 
700     Input Parameter:
701 +   dm - the DM object
702 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
703 
704     Output Parameter:
705 .   dmf - the refined DM
706 
707     Level: developer
708 
709 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
710 
711 @*/
712 PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
713 {
714   PetscErrorCode ierr;
715 
716   PetscFunctionBegin;
717   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
718   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
719   (*dmf)->ops->initialguess = dm->ops->initialguess;
720   (*dmf)->ops->function     = dm->ops->function;
721   (*dmf)->ops->functionj    = dm->ops->functionj;
722   if (dm->ops->jacobian != DMComputeJacobianDefault) {
723     (*dmf)->ops->jacobian     = dm->ops->jacobian;
724   }
725   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
726   (*dmf)->ctx     = dm->ctx;
727   (*dmf)->levelup = dm->levelup + 1;
728   PetscFunctionReturn(0);
729 }
730 
731 #undef __FUNCT__
732 #define __FUNCT__ "DMGetRefineLevel"
733 /*@
734     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
735 
736     Not Collective
737 
738     Input Parameter:
739 .   dm - the DM object
740 
741     Output Parameter:
742 .   level - number of refinements
743 
744     Level: developer
745 
746 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
747 
748 @*/
749 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
750 {
751   PetscFunctionBegin;
752   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
753   *level = dm->levelup;
754   PetscFunctionReturn(0);
755 }
756 
757 #undef __FUNCT__
758 #define __FUNCT__ "DMGlobalToLocalBegin"
759 /*@
760     DMGlobalToLocalBegin - Begins updating local vectors from global vector
761 
762     Neighbor-wise Collective on DM
763 
764     Input Parameters:
765 +   dm - the DM object
766 .   g - the global vector
767 .   mode - INSERT_VALUES or ADD_VALUES
768 -   l - the local vector
769 
770 
771     Level: beginner
772 
773 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
774 
775 @*/
776 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
777 {
778   PetscErrorCode ierr;
779 
780   PetscFunctionBegin;
781   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
782   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
783   PetscFunctionReturn(0);
784 }
785 
786 #undef __FUNCT__
787 #define __FUNCT__ "DMGlobalToLocalEnd"
788 /*@
789     DMGlobalToLocalEnd - Ends updating local vectors from global vector
790 
791     Neighbor-wise Collective on DM
792 
793     Input Parameters:
794 +   dm - the DM object
795 .   g - the global vector
796 .   mode - INSERT_VALUES or ADD_VALUES
797 -   l - the local vector
798 
799 
800     Level: beginner
801 
802 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
803 
804 @*/
805 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
806 {
807   PetscErrorCode ierr;
808 
809   PetscFunctionBegin;
810   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
811   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
812   PetscFunctionReturn(0);
813 }
814 
815 #undef __FUNCT__
816 #define __FUNCT__ "DMLocalToGlobalBegin"
817 /*@
818     DMLocalToGlobalBegin - updates global vectors from local vectors
819 
820     Neighbor-wise Collective on DM
821 
822     Input Parameters:
823 +   dm - the DM object
824 .   l - the local vector
825 .   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
826            base point.
827 - - the global vector
828 
829     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
830            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
831            global array to the final global array with VecAXPY().
832 
833     Level: beginner
834 
835 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
836 
837 @*/
838 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
839 {
840   PetscErrorCode ierr;
841 
842   PetscFunctionBegin;
843   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
844   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
845   PetscFunctionReturn(0);
846 }
847 
848 #undef __FUNCT__
849 #define __FUNCT__ "DMLocalToGlobalEnd"
850 /*@
851     DMLocalToGlobalEnd - updates global vectors from local vectors
852 
853     Neighbor-wise Collective on DM
854 
855     Input Parameters:
856 +   dm - the DM object
857 .   l - the local vector
858 .   mode - INSERT_VALUES or ADD_VALUES
859 -   g - the global vector
860 
861 
862     Level: beginner
863 
864 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
865 
866 @*/
867 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
868 {
869   PetscErrorCode ierr;
870 
871   PetscFunctionBegin;
872   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
873   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
874   PetscFunctionReturn(0);
875 }
876 
877 #undef __FUNCT__
878 #define __FUNCT__ "DMComputeJacobianDefault"
879 /*@
880     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
881 
882     Collective on DM
883 
884     Input Parameter:
885 +   dm - the DM object
886 .   x - location to compute Jacobian at; may be ignored for linear problems
887 .   A - matrix that defines the operator for the linear solve
888 -   B - the matrix used to construct the preconditioner
889 
890     Level: developer
891 
892 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
893          DMSetFunction()
894 
895 @*/
896 PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
897 {
898   PetscErrorCode ierr;
899 
900   PetscFunctionBegin;
901   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
902   *stflag = SAME_NONZERO_PATTERN;
903   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
904   if (A != B) {
905     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
906     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
907   }
908   PetscFunctionReturn(0);
909 }
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "DMCoarsen"
913 /*@
914     DMCoarsen - Coarsens a DM object
915 
916     Collective on DM
917 
918     Input Parameter:
919 +   dm - the DM object
920 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
921 
922     Output Parameter:
923 .   dmc - the coarsened DM
924 
925     Level: developer
926 
927 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
928 
929 @*/
930 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
931 {
932   PetscErrorCode ierr;
933 
934   PetscFunctionBegin;
935   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
936   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
937   (*dmc)->ops->initialguess = dm->ops->initialguess;
938   (*dmc)->ops->function     = dm->ops->function;
939   (*dmc)->ops->functionj    = dm->ops->functionj;
940   if (dm->ops->jacobian != DMComputeJacobianDefault) {
941     (*dmc)->ops->jacobian     = dm->ops->jacobian;
942   }
943   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
944   (*dmc)->ctx       = dm->ctx;
945   (*dmc)->leveldown = dm->leveldown + 1;
946   PetscFunctionReturn(0);
947 }
948 
949 #undef __FUNCT__
950 #define __FUNCT__ "DMGetCoarsenLevel"
951 /*@
952     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
953 
954     Not Collective
955 
956     Input Parameter:
957 .   dm - the DM object
958 
959     Output Parameter:
960 .   level - number of coarsenings
961 
962     Level: developer
963 
964 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
965 
966 @*/
967 PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
968 {
969   PetscFunctionBegin;
970   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
971   *level = dm->leveldown;
972   PetscFunctionReturn(0);
973 }
974 
975 
976 
977 #undef __FUNCT__
978 #define __FUNCT__ "DMRefineHierarchy"
979 /*@C
980     DMRefineHierarchy - Refines a DM object, all levels at once
981 
982     Collective on DM
983 
984     Input Parameter:
985 +   dm - the DM object
986 -   nlevels - the number of levels of refinement
987 
988     Output Parameter:
989 .   dmf - the refined DM hierarchy
990 
991     Level: developer
992 
993 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
994 
995 @*/
996 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
997 {
998   PetscErrorCode ierr;
999 
1000   PetscFunctionBegin;
1001   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1002   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1003   if (nlevels == 0) PetscFunctionReturn(0);
1004   if (dm->ops->refinehierarchy) {
1005     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
1006   } else if (dm->ops->refine) {
1007     PetscInt i;
1008 
1009     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
1010     for (i=1; i<nlevels; i++) {
1011       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
1012     }
1013   } else {
1014     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
1015   }
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 #undef __FUNCT__
1020 #define __FUNCT__ "DMCoarsenHierarchy"
1021 /*@C
1022     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
1023 
1024     Collective on DM
1025 
1026     Input Parameter:
1027 +   dm - the DM object
1028 -   nlevels - the number of levels of coarsening
1029 
1030     Output Parameter:
1031 .   dmc - the coarsened DM hierarchy
1032 
1033     Level: developer
1034 
1035 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1036 
1037 @*/
1038 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
1039 {
1040   PetscErrorCode ierr;
1041 
1042   PetscFunctionBegin;
1043   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1044   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1045   if (nlevels == 0) PetscFunctionReturn(0);
1046   PetscValidPointer(dmc,3);
1047   if (dm->ops->coarsenhierarchy) {
1048     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
1049   } else if (dm->ops->coarsen) {
1050     PetscInt i;
1051 
1052     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
1053     for (i=1; i<nlevels; i++) {
1054       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
1055     }
1056   } else {
1057     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
1058   }
1059   PetscFunctionReturn(0);
1060 }
1061 
1062 #undef __FUNCT__
1063 #define __FUNCT__ "DMCreateAggregates"
1064 /*@
1065    DMCreateAggregates - Gets the aggregates that map between
1066    grids associated with two DMs.
1067 
1068    Collective on DM
1069 
1070    Input Parameters:
1071 +  dmc - the coarse grid DM
1072 -  dmf - the fine grid DM
1073 
1074    Output Parameters:
1075 .  rest - the restriction matrix (transpose of the projection matrix)
1076 
1077    Level: intermediate
1078 
1079 .keywords: interpolation, restriction, multigrid
1080 
1081 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
1082 @*/
1083 PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
1084 {
1085   PetscErrorCode ierr;
1086 
1087   PetscFunctionBegin;
1088   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1089   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
1090   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
1091   PetscFunctionReturn(0);
1092 }
1093 
1094 #undef __FUNCT__
1095 #define __FUNCT__ "DMSetApplicationContextDestroy"
1096 /*@C
1097     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
1098 
1099     Not Collective
1100 
1101     Input Parameters:
1102 +   dm - the DM object
1103 -   destroy - the destroy function
1104 
1105     Level: intermediate
1106 
1107 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1108 
1109 @*/
1110 PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
1111 {
1112   PetscFunctionBegin;
1113   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1114   dm->ctxdestroy = destroy;
1115   PetscFunctionReturn(0);
1116 }
1117 
1118 #undef __FUNCT__
1119 #define __FUNCT__ "DMSetApplicationContext"
1120 /*@
1121     DMSetApplicationContext - Set a user context into a DM object
1122 
1123     Not Collective
1124 
1125     Input Parameters:
1126 +   dm - the DM object
1127 -   ctx - the user context
1128 
1129     Level: intermediate
1130 
1131 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1132 
1133 @*/
1134 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
1135 {
1136   PetscFunctionBegin;
1137   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1138   dm->ctx = ctx;
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 #undef __FUNCT__
1143 #define __FUNCT__ "DMGetApplicationContext"
1144 /*@
1145     DMGetApplicationContext - Gets a user context from a DM object
1146 
1147     Not Collective
1148 
1149     Input Parameter:
1150 .   dm - the DM object
1151 
1152     Output Parameter:
1153 .   ctx - the user context
1154 
1155     Level: intermediate
1156 
1157 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1158 
1159 @*/
1160 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
1161 {
1162   PetscFunctionBegin;
1163   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1164   *(void**)ctx = dm->ctx;
1165   PetscFunctionReturn(0);
1166 }
1167 
1168 #undef __FUNCT__
1169 #define __FUNCT__ "DMSetInitialGuess"
1170 /*@C
1171     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
1172 
1173     Logically Collective on DM
1174 
1175     Input Parameter:
1176 +   dm - the DM object to destroy
1177 -   f - the function to compute the initial guess
1178 
1179     Level: intermediate
1180 
1181 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1182 
1183 @*/
1184 PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
1185 {
1186   PetscFunctionBegin;
1187   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1188   dm->ops->initialguess = f;
1189   PetscFunctionReturn(0);
1190 }
1191 
1192 #undef __FUNCT__
1193 #define __FUNCT__ "DMSetFunction"
1194 /*@C
1195     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
1196 
1197     Logically Collective on DM
1198 
1199     Input Parameter:
1200 +   dm - the DM object
1201 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
1202 
1203     Level: intermediate
1204 
1205     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
1206            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
1207 
1208 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1209          DMSetJacobian()
1210 
1211 @*/
1212 PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1213 {
1214   PetscFunctionBegin;
1215   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1216   dm->ops->function = f;
1217   if (f) {
1218     dm->ops->functionj = f;
1219   }
1220   PetscFunctionReturn(0);
1221 }
1222 
1223 #undef __FUNCT__
1224 #define __FUNCT__ "DMSetJacobian"
1225 /*@C
1226     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
1227 
1228     Logically Collective on DM
1229 
1230     Input Parameter:
1231 +   dm - the DM object to destroy
1232 -   f - the function to compute the matrix entries
1233 
1234     Level: intermediate
1235 
1236 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1237          DMSetFunction()
1238 
1239 @*/
1240 PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
1241 {
1242   PetscFunctionBegin;
1243   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1244   dm->ops->jacobian = f;
1245   PetscFunctionReturn(0);
1246 }
1247 
1248 #undef __FUNCT__
1249 #define __FUNCT__ "DMComputeInitialGuess"
1250 /*@
1251     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
1252 
1253     Collective on DM
1254 
1255     Input Parameter:
1256 +   dm - the DM object to destroy
1257 -   x - the vector to hold the initial guess values
1258 
1259     Level: developer
1260 
1261 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
1262 
1263 @*/
1264 PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
1265 {
1266   PetscErrorCode ierr;
1267   PetscFunctionBegin;
1268   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1269   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
1270   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
1271   PetscFunctionReturn(0);
1272 }
1273 
1274 #undef __FUNCT__
1275 #define __FUNCT__ "DMHasInitialGuess"
1276 /*@
1277     DMHasInitialGuess - does the DM object have an initial guess function
1278 
1279     Not Collective
1280 
1281     Input Parameter:
1282 .   dm - the DM object to destroy
1283 
1284     Output Parameter:
1285 .   flg - PETSC_TRUE if function exists
1286 
1287     Level: developer
1288 
1289 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1290 
1291 @*/
1292 PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
1293 {
1294   PetscFunctionBegin;
1295   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1296   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
1297   PetscFunctionReturn(0);
1298 }
1299 
1300 #undef __FUNCT__
1301 #define __FUNCT__ "DMHasFunction"
1302 /*@
1303     DMHasFunction - does the DM object have a function
1304 
1305     Not Collective
1306 
1307     Input Parameter:
1308 .   dm - the DM object to destroy
1309 
1310     Output Parameter:
1311 .   flg - PETSC_TRUE if function exists
1312 
1313     Level: developer
1314 
1315 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1316 
1317 @*/
1318 PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
1319 {
1320   PetscFunctionBegin;
1321   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1322   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
1323   PetscFunctionReturn(0);
1324 }
1325 
1326 #undef __FUNCT__
1327 #define __FUNCT__ "DMHasJacobian"
1328 /*@
1329     DMHasJacobian - does the DM object have a matrix function
1330 
1331     Not Collective
1332 
1333     Input Parameter:
1334 .   dm - the DM object to destroy
1335 
1336     Output Parameter:
1337 .   flg - PETSC_TRUE if function exists
1338 
1339     Level: developer
1340 
1341 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1342 
1343 @*/
1344 PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
1345 {
1346   PetscFunctionBegin;
1347   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1348   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
1349   PetscFunctionReturn(0);
1350 }
1351 
1352 #undef __FUNCT__
1353 #define __FUNCT__ "DMComputeFunction"
1354 /*@
1355     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
1356 
1357     Collective on DM
1358 
1359     Input Parameter:
1360 +   dm - the DM object to destroy
1361 .   x - the location where the function is evaluationed, may be ignored for linear problems
1362 -   b - the vector to hold the right hand side entries
1363 
1364     Level: developer
1365 
1366 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1367          DMSetJacobian()
1368 
1369 @*/
1370 PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
1371 {
1372   PetscErrorCode ierr;
1373   PetscFunctionBegin;
1374   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1375   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1376   PetscStackPush("DM user function");
1377   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1378   PetscStackPop;
1379   PetscFunctionReturn(0);
1380 }
1381 
1382 
1383 #undef __FUNCT__
1384 #define __FUNCT__ "DMComputeJacobian"
1385 /*@
1386     DMComputeJacobian - compute the matrix entries for the solver
1387 
1388     Collective on DM
1389 
1390     Input Parameter:
1391 +   dm - the DM object
1392 .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
1393 .   A - matrix that defines the operator for the linear solve
1394 -   B - the matrix used to construct the preconditioner
1395 
1396     Level: developer
1397 
1398 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1399          DMSetFunction()
1400 
1401 @*/
1402 PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1403 {
1404   PetscErrorCode ierr;
1405 
1406   PetscFunctionBegin;
1407   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1408   if (!dm->ops->jacobian) {
1409     ISColoring     coloring;
1410     MatFDColoring  fd;
1411 
1412     ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
1413     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1414     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1415     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1416     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1417     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
1418 
1419     dm->fd = fd;
1420     dm->ops->jacobian = DMComputeJacobianDefault;
1421 
1422     /* don't know why this is needed */
1423     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1424   }
1425   if (!x) x = dm->x;
1426   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1427 
1428   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1429   if (x) {
1430     if (!dm->x) {
1431       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1432     }
1433     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1434   }
1435   if (A != B) {
1436     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1437     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1438   }
1439   PetscFunctionReturn(0);
1440 }
1441 
1442 
1443 PetscFList DMList                       = PETSC_NULL;
1444 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1445 
1446 #undef __FUNCT__
1447 #define __FUNCT__ "DMSetType"
1448 /*@C
1449   DMSetType - Builds a DM, for a particular DM implementation.
1450 
1451   Collective on DM
1452 
1453   Input Parameters:
1454 + dm     - The DM object
1455 - method - The name of the DM type
1456 
1457   Options Database Key:
1458 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1459 
1460   Notes:
1461   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1462 
1463   Level: intermediate
1464 
1465 .keywords: DM, set, type
1466 .seealso: DMGetType(), DMCreate()
1467 @*/
1468 PetscErrorCode  DMSetType(DM dm, const DMType method)
1469 {
1470   PetscErrorCode (*r)(DM);
1471   PetscBool      match;
1472   PetscErrorCode ierr;
1473 
1474   PetscFunctionBegin;
1475   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1476   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1477   if (match) PetscFunctionReturn(0);
1478 
1479   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1480   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1481   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1482 
1483   if (dm->ops->destroy) {
1484     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1485   }
1486   ierr = (*r)(dm);CHKERRQ(ierr);
1487   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1488   PetscFunctionReturn(0);
1489 }
1490 
1491 #undef __FUNCT__
1492 #define __FUNCT__ "DMGetType"
1493 /*@C
1494   DMGetType - Gets the DM type name (as a string) from the DM.
1495 
1496   Not Collective
1497 
1498   Input Parameter:
1499 . dm  - The DM
1500 
1501   Output Parameter:
1502 . type - The DM type name
1503 
1504   Level: intermediate
1505 
1506 .keywords: DM, get, type, name
1507 .seealso: DMSetType(), DMCreate()
1508 @*/
1509 PetscErrorCode  DMGetType(DM dm, const DMType *type)
1510 {
1511   PetscErrorCode ierr;
1512 
1513   PetscFunctionBegin;
1514   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1515   PetscValidCharPointer(type,2);
1516   if (!DMRegisterAllCalled) {
1517     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1518   }
1519   *type = ((PetscObject)dm)->type_name;
1520   PetscFunctionReturn(0);
1521 }
1522 
1523 #undef __FUNCT__
1524 #define __FUNCT__ "DMConvert"
1525 /*@C
1526   DMConvert - Converts a DM to another DM, either of the same or different type.
1527 
1528   Collective on DM
1529 
1530   Input Parameters:
1531 + dm - the DM
1532 - newtype - new DM type (use "same" for the same type)
1533 
1534   Output Parameter:
1535 . M - pointer to new DM
1536 
1537   Notes:
1538   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
1539   the MPI communicator of the generated DM is always the same as the communicator
1540   of the input DM.
1541 
1542   Level: intermediate
1543 
1544 .seealso: DMCreate()
1545 @*/
1546 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
1547 {
1548   DM             B;
1549   char           convname[256];
1550   PetscBool      sametype, issame;
1551   PetscErrorCode ierr;
1552 
1553   PetscFunctionBegin;
1554   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1555   PetscValidType(dm,1);
1556   PetscValidPointer(M,3);
1557   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
1558   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
1559   {
1560     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
1561 
1562     /*
1563        Order of precedence:
1564        1) See if a specialized converter is known to the current DM.
1565        2) See if a specialized converter is known to the desired DM class.
1566        3) See if a good general converter is registered for the desired class
1567        4) See if a good general converter is known for the current matrix.
1568        5) Use a really basic converter.
1569     */
1570 
1571     /* 1) See if a specialized converter is known to the current DM and the desired class */
1572     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1573     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1574     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1575     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1576     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1577     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1578     if (conv) goto foundconv;
1579 
1580     /* 2)  See if a specialized converter is known to the desired DM class. */
1581     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
1582     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
1583     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1584     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1585     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1586     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1587     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1588     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1589     if (conv) {
1590       ierr = DMDestroy(&B);CHKERRQ(ierr);
1591       goto foundconv;
1592     }
1593 
1594 #if 0
1595     /* 3) See if a good general converter is registered for the desired class */
1596     conv = B->ops->convertfrom;
1597     ierr = DMDestroy(&B);CHKERRQ(ierr);
1598     if (conv) goto foundconv;
1599 
1600     /* 4) See if a good general converter is known for the current matrix */
1601     if (dm->ops->convert) {
1602       conv = dm->ops->convert;
1603     }
1604     if (conv) goto foundconv;
1605 #endif
1606 
1607     /* 5) Use a really basic converter. */
1608     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
1609 
1610     foundconv:
1611     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1612     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
1613     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1614   }
1615   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
1616   PetscFunctionReturn(0);
1617 }
1618 
1619 /*--------------------------------------------------------------------------------------------------------------------*/
1620 
1621 #undef __FUNCT__
1622 #define __FUNCT__ "DMRegister"
1623 /*@C
1624   DMRegister - See DMRegisterDynamic()
1625 
1626   Level: advanced
1627 @*/
1628 PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1629 {
1630   char fullname[PETSC_MAX_PATH_LEN];
1631   PetscErrorCode ierr;
1632 
1633   PetscFunctionBegin;
1634   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1635   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1636   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1637   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1638   PetscFunctionReturn(0);
1639 }
1640 
1641 
1642 /*--------------------------------------------------------------------------------------------------------------------*/
1643 #undef __FUNCT__
1644 #define __FUNCT__ "DMRegisterDestroy"
1645 /*@C
1646    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1647 
1648    Not Collective
1649 
1650    Level: advanced
1651 
1652 .keywords: DM, register, destroy
1653 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1654 @*/
1655 PetscErrorCode  DMRegisterDestroy(void)
1656 {
1657   PetscErrorCode ierr;
1658 
1659   PetscFunctionBegin;
1660   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1661   DMRegisterAllCalled = PETSC_FALSE;
1662   PetscFunctionReturn(0);
1663 }
1664 
1665 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1666 #include <mex.h>
1667 
1668 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
1669 
1670 #undef __FUNCT__
1671 #define __FUNCT__ "DMComputeFunction_Matlab"
1672 /*
1673    DMComputeFunction_Matlab - Calls the function that has been set with
1674                          DMSetFunctionMatlab().
1675 
1676    For linear problems x is null
1677 
1678 .seealso: DMSetFunction(), DMGetFunction()
1679 */
1680 PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
1681 {
1682   PetscErrorCode    ierr;
1683   DMMatlabContext   *sctx;
1684   int               nlhs = 1,nrhs = 4;
1685   mxArray	    *plhs[1],*prhs[4];
1686   long long int     lx = 0,ly = 0,ls = 0;
1687 
1688   PetscFunctionBegin;
1689   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1690   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1691   PetscCheckSameComm(dm,1,y,3);
1692 
1693   /* call Matlab function in ctx with arguments x and y */
1694   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1695   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1696   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1697   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
1698   prhs[0] =  mxCreateDoubleScalar((double)ls);
1699   prhs[1] =  mxCreateDoubleScalar((double)lx);
1700   prhs[2] =  mxCreateDoubleScalar((double)ly);
1701   prhs[3] =  mxCreateString(sctx->funcname);
1702   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
1703   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
1704   mxDestroyArray(prhs[0]);
1705   mxDestroyArray(prhs[1]);
1706   mxDestroyArray(prhs[2]);
1707   mxDestroyArray(prhs[3]);
1708   mxDestroyArray(plhs[0]);
1709   PetscFunctionReturn(0);
1710 }
1711 
1712 
1713 #undef __FUNCT__
1714 #define __FUNCT__ "DMSetFunctionMatlab"
1715 /*
1716    DMSetFunctionMatlab - Sets the function evaluation routine
1717 
1718 */
1719 PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
1720 {
1721   PetscErrorCode    ierr;
1722   DMMatlabContext   *sctx;
1723 
1724   PetscFunctionBegin;
1725   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1726   /* currently sctx is memory bleed */
1727   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1728   if (!sctx) {
1729     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1730   }
1731   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
1732   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1733   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
1734   PetscFunctionReturn(0);
1735 }
1736 
1737 #undef __FUNCT__
1738 #define __FUNCT__ "DMComputeJacobian_Matlab"
1739 /*
1740    DMComputeJacobian_Matlab - Calls the function that has been set with
1741                          DMSetJacobianMatlab().
1742 
1743    For linear problems x is null
1744 
1745 .seealso: DMSetFunction(), DMGetFunction()
1746 */
1747 PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
1748 {
1749   PetscErrorCode    ierr;
1750   DMMatlabContext   *sctx;
1751   int               nlhs = 2,nrhs = 5;
1752   mxArray	    *plhs[2],*prhs[5];
1753   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
1754 
1755   PetscFunctionBegin;
1756   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1757   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
1758 
1759   /* call MATLAB function in ctx with arguments x, A, and B */
1760   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1761   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1762   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1763   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
1764   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
1765   prhs[0] =  mxCreateDoubleScalar((double)ls);
1766   prhs[1] =  mxCreateDoubleScalar((double)lx);
1767   prhs[2] =  mxCreateDoubleScalar((double)lA);
1768   prhs[3] =  mxCreateDoubleScalar((double)lB);
1769   prhs[4] =  mxCreateString(sctx->jacname);
1770   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
1771   *str    =  (MatStructure) mxGetScalar(plhs[0]);
1772   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
1773   mxDestroyArray(prhs[0]);
1774   mxDestroyArray(prhs[1]);
1775   mxDestroyArray(prhs[2]);
1776   mxDestroyArray(prhs[3]);
1777   mxDestroyArray(prhs[4]);
1778   mxDestroyArray(plhs[0]);
1779   mxDestroyArray(plhs[1]);
1780   PetscFunctionReturn(0);
1781 }
1782 
1783 
1784 #undef __FUNCT__
1785 #define __FUNCT__ "DMSetJacobianMatlab"
1786 /*
1787    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
1788 
1789 */
1790 PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
1791 {
1792   PetscErrorCode    ierr;
1793   DMMatlabContext   *sctx;
1794 
1795   PetscFunctionBegin;
1796   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1797   /* currently sctx is memory bleed */
1798   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1799   if (!sctx) {
1800     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1801   }
1802   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
1803   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1804   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
1805   PetscFunctionReturn(0);
1806 }
1807 #endif
1808 
1809 #undef __FUNCT__
1810 #define __FUNCT__ "DMLoad"
1811 /*@C
1812   DMLoad - Loads a DM that has been stored in binary or HDF5 format
1813   with DMView().
1814 
1815   Collective on PetscViewer
1816 
1817   Input Parameters:
1818 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
1819            some related function before a call to DMLoad().
1820 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
1821            HDF5 file viewer, obtained from PetscViewerHDF5Open()
1822 
1823    Level: intermediate
1824 
1825   Notes:
1826   Defaults to the DM DA.
1827 
1828   Notes for advanced users:
1829   Most users should not need to know the details of the binary storage
1830   format, since DMLoad() and DMView() completely hide these details.
1831   But for anyone who's interested, the standard binary matrix storage
1832   format is
1833 .vb
1834      has not yet been determined
1835 .ve
1836 
1837    In addition, PETSc automatically does the byte swapping for
1838 machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
1839 linux, Windows and the paragon; thus if you write your own binary
1840 read/write routines you have to swap the bytes; see PetscBinaryRead()
1841 and PetscBinaryWrite() to see how this may be done.
1842 
1843   Concepts: vector^loading from file
1844 
1845 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
1846 @*/
1847 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
1848 {
1849   PetscErrorCode ierr;
1850 
1851   PetscFunctionBegin;
1852   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
1853   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1854 
1855   if (!((PetscObject)newdm)->type_name) {
1856     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
1857   }
1858   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1859   PetscFunctionReturn(0);
1860 }
1861 
1862