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