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