xref: /petsc/src/dm/interface/dm.c (revision ff218e97a57ed641f3ebc93f697e38ef0f3aa217)
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 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix()
441 
442 @*/
443 PetscErrorCode  DMGetInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
444 {
445   PetscErrorCode ierr;
446 
447   PetscFunctionBegin;
448   ierr = (*dm1->ops->getinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "DMGetInjection"
454 /*@
455     DMGetInjection - Gets injection matrix between two DMDA or DMComposite objects
456 
457     Collective on DM
458 
459     Input Parameter:
460 +   dm1 - the DM object
461 -   dm2 - the second, finer DM object
462 
463     Output Parameter:
464 .   ctx - the injection
465 
466     Level: developer
467 
468 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMGetInterpolation()
469 
470 @*/
471 PetscErrorCode  DMGetInjection(DM dm1,DM dm2,VecScatter *ctx)
472 {
473   PetscErrorCode ierr;
474 
475   PetscFunctionBegin;
476   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
477   PetscFunctionReturn(0);
478 }
479 
480 #undef __FUNCT__
481 #define __FUNCT__ "DMGetColoring"
482 /*@C
483     DMGetColoring - Gets coloring for a DMDA or DMComposite
484 
485     Collective on DM
486 
487     Input Parameter:
488 +   dm - the DM object
489 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
490 -   matype - either MATAIJ or MATBAIJ
491 
492     Output Parameter:
493 .   coloring - the coloring
494 
495     Level: developer
496 
497 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
498 
499 @*/
500 PetscErrorCode  DMGetColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
501 {
502   PetscErrorCode ierr;
503 
504   PetscFunctionBegin;
505   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
506   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "DMGetMatrix"
512 /*@C
513     DMGetMatrix - Gets empty Jacobian for a DMDA or DMComposite
514 
515     Collective on DM
516 
517     Input Parameter:
518 +   dm - the DM object
519 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
520             any type which inherits from one of these (such as MATAIJ)
521 
522     Output Parameter:
523 .   mat - the empty Jacobian
524 
525     Level: developer
526 
527     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
528        do not need to do it yourself.
529 
530        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
531        the nonzero pattern call DMDASetMatPreallocateOnly()
532 
533        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
534        internally by PETSc.
535 
536        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
537        the indices for the global numbering for DMDAs which is complicated.
538 
539 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
540 
541 @*/
542 PetscErrorCode  DMGetMatrix(DM dm,const MatType mtype,Mat *mat)
543 {
544   PetscErrorCode ierr;
545   char           ttype[256];
546   PetscBool      flg;
547 
548   PetscFunctionBegin;
549 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
550   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
551 #endif
552   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
553   PetscValidPointer(mat,3);
554   ierr = PetscStrncpy(ttype,mtype,sizeof(ttype));CHKERRQ(ierr);
555   ttype[sizeof(ttype)-1] = 0;
556   ierr = PetscOptionsBegin(((PetscObject)dm)->comm,((PetscObject)dm)->prefix,"DM options","Mat");CHKERRQ(ierr);
557   ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,ttype,ttype,sizeof(ttype),&flg);CHKERRQ(ierr);
558   ierr = PetscOptionsEnd();
559   if (flg || mtype) {
560     ierr = (*dm->ops->getmatrix)(dm,ttype,mat);CHKERRQ(ierr);
561   } else {                      /* Let the implementation decide */
562     ierr = (*dm->ops->getmatrix)(dm,PETSC_NULL,mat);CHKERRQ(ierr);
563   }
564   PetscFunctionReturn(0);
565 }
566 
567 #undef __FUNCT__
568 #define __FUNCT__ "DMSetMatrixPreallocateOnly"
569 /*@
570   DMSetMatrixPreallocateOnly - When DMGetMatrix() is called the matrix will be properly
571     preallocated but the nonzero structure and zero values will not be set.
572 
573   Logically Collective on DMDA
574 
575   Input Parameter:
576 + dm - the DM
577 - only - PETSC_TRUE if only want preallocation
578 
579   Level: developer
580 .seealso DMGetMatrix()
581 @*/
582 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
583 {
584   PetscFunctionBegin;
585   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
586   dm->prealloc_only = only;
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "DMRefine"
592 /*@
593     DMRefine - Refines a DM object
594 
595     Collective on DM
596 
597     Input Parameter:
598 +   dm - the DM object
599 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
600 
601     Output Parameter:
602 .   dmf - the refined DM
603 
604     Level: developer
605 
606 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
607 
608 @*/
609 PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
610 {
611   PetscErrorCode ierr;
612 
613   PetscFunctionBegin;
614   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
615   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
616   (*dmf)->ops->initialguess = dm->ops->initialguess;
617   (*dmf)->ops->function     = dm->ops->function;
618   (*dmf)->ops->functionj    = dm->ops->functionj;
619   if (dm->ops->jacobian != DMComputeJacobianDefault) {
620     (*dmf)->ops->jacobian     = dm->ops->jacobian;
621   }
622   (*dmf)->ctx     = dm->ctx;
623   (*dmf)->levelup = dm->levelup + 1;
624   PetscFunctionReturn(0);
625 }
626 
627 #undef __FUNCT__
628 #define __FUNCT__ "DMGetRefineLevel"
629 /*@
630     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
631 
632     Not Collective
633 
634     Input Parameter:
635 .   dm - the DM object
636 
637     Output Parameter:
638 .   level - number of refinements
639 
640     Level: developer
641 
642 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
643 
644 @*/
645 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
646 {
647   PetscFunctionBegin;
648   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
649   *level = dm->levelup;
650   PetscFunctionReturn(0);
651 }
652 
653 #undef __FUNCT__
654 #define __FUNCT__ "DMGlobalToLocalBegin"
655 /*@
656     DMGlobalToLocalBegin - Begins updating local vectors from global vector
657 
658     Neighbor-wise Collective on DM
659 
660     Input Parameters:
661 +   dm - the DM object
662 .   g - the global vector
663 .   mode - INSERT_VALUES or ADD_VALUES
664 -   l - the local vector
665 
666 
667     Level: beginner
668 
669 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
670 
671 @*/
672 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
673 {
674   PetscErrorCode ierr;
675 
676   PetscFunctionBegin;
677   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr);
678   PetscFunctionReturn(0);
679 }
680 
681 #undef __FUNCT__
682 #define __FUNCT__ "DMGlobalToLocalEnd"
683 /*@
684     DMGlobalToLocalEnd - Ends updating local vectors from global vector
685 
686     Neighbor-wise Collective on DM
687 
688     Input Parameters:
689 +   dm - the DM object
690 .   g - the global vector
691 .   mode - INSERT_VALUES or ADD_VALUES
692 -   l - the local vector
693 
694 
695     Level: beginner
696 
697 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
698 
699 @*/
700 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
701 {
702   PetscErrorCode ierr;
703 
704   PetscFunctionBegin;
705   ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr);
706   PetscFunctionReturn(0);
707 }
708 
709 #undef __FUNCT__
710 #define __FUNCT__ "DMLocalToGlobalBegin"
711 /*@
712     DMLocalToGlobalBegin - updates global vectors from local vectors
713 
714     Neighbor-wise Collective on DM
715 
716     Input Parameters:
717 +   dm - the DM object
718 .   l - the local vector
719 .   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
720            base point.
721 - - the global vector
722 
723     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
724            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
725            global array to the final global array with VecAXPY().
726 
727     Level: beginner
728 
729 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
730 
731 @*/
732 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
733 {
734   PetscErrorCode ierr;
735 
736   PetscFunctionBegin;
737   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode,g);CHKERRQ(ierr);
738   PetscFunctionReturn(0);
739 }
740 
741 #undef __FUNCT__
742 #define __FUNCT__ "DMLocalToGlobalEnd"
743 /*@
744     DMLocalToGlobalEnd - updates global vectors from local vectors
745 
746     Neighbor-wise Collective on DM
747 
748     Input Parameters:
749 +   dm - the DM object
750 .   l - the local vector
751 .   mode - INSERT_VALUES or ADD_VALUES
752 -   g - the global vector
753 
754 
755     Level: beginner
756 
757 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
758 
759 @*/
760 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
761 {
762   PetscErrorCode ierr;
763 
764   PetscFunctionBegin;
765   ierr = (*dm->ops->localtoglobalend)(dm,l,mode,g);CHKERRQ(ierr);
766   PetscFunctionReturn(0);
767 }
768 
769 #undef __FUNCT__
770 #define __FUNCT__ "DMComputeJacobianDefault"
771 /*@
772     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
773 
774     Collective on DM
775 
776     Input Parameter:
777 +   dm - the DM object
778 .   x - location to compute Jacobian at; may be ignored for linear problems
779 .   A - matrix that defines the operator for the linear solve
780 -   B - the matrix used to construct the preconditioner
781 
782     Level: developer
783 
784 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
785          DMSetFunction()
786 
787 @*/
788 PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
789 {
790   PetscErrorCode ierr;
791   PetscFunctionBegin;
792   *stflag = SAME_NONZERO_PATTERN;
793   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
794   if (A != B) {
795     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
796     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
797   }
798   PetscFunctionReturn(0);
799 }
800 
801 #undef __FUNCT__
802 #define __FUNCT__ "DMCoarsen"
803 /*@
804     DMCoarsen - Coarsens a DM object
805 
806     Collective on DM
807 
808     Input Parameter:
809 +   dm - the DM object
810 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
811 
812     Output Parameter:
813 .   dmc - the coarsened DM
814 
815     Level: developer
816 
817 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
818 
819 @*/
820 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
821 {
822   PetscErrorCode ierr;
823 
824   PetscFunctionBegin;
825   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
826   (*dmc)->ops->initialguess = dm->ops->initialguess;
827   (*dmc)->ops->function     = dm->ops->function;
828   (*dmc)->ops->functionj    = dm->ops->functionj;
829   if (dm->ops->jacobian != DMComputeJacobianDefault) {
830     (*dmc)->ops->jacobian     = dm->ops->jacobian;
831   }
832   (*dmc)->ctx       = dm->ctx;
833   (*dmc)->leveldown = dm->leveldown + 1;
834   PetscFunctionReturn(0);
835 }
836 
837 #undef __FUNCT__
838 #define __FUNCT__ "DMRefineHierarchy"
839 /*@C
840     DMRefineHierarchy - Refines a DM object, all levels at once
841 
842     Collective on DM
843 
844     Input Parameter:
845 +   dm - the DM object
846 -   nlevels - the number of levels of refinement
847 
848     Output Parameter:
849 .   dmf - the refined DM hierarchy
850 
851     Level: developer
852 
853 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
854 
855 @*/
856 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
857 {
858   PetscErrorCode ierr;
859 
860   PetscFunctionBegin;
861   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
862   if (nlevels == 0) PetscFunctionReturn(0);
863   if (dm->ops->refinehierarchy) {
864     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
865   } else if (dm->ops->refine) {
866     PetscInt i;
867 
868     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
869     for (i=1; i<nlevels; i++) {
870       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
871     }
872   } else {
873     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
874   }
875   PetscFunctionReturn(0);
876 }
877 
878 #undef __FUNCT__
879 #define __FUNCT__ "DMCoarsenHierarchy"
880 /*@C
881     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
882 
883     Collective on DM
884 
885     Input Parameter:
886 +   dm - the DM object
887 -   nlevels - the number of levels of coarsening
888 
889     Output Parameter:
890 .   dmc - the coarsened DM hierarchy
891 
892     Level: developer
893 
894 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
895 
896 @*/
897 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
898 {
899   PetscErrorCode ierr;
900 
901   PetscFunctionBegin;
902   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
903   if (nlevels == 0) PetscFunctionReturn(0);
904   PetscValidPointer(dmc,3);
905   if (dm->ops->coarsenhierarchy) {
906     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
907   } else if (dm->ops->coarsen) {
908     PetscInt i;
909 
910     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
911     for (i=1; i<nlevels; i++) {
912       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
913     }
914   } else {
915     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
916   }
917   PetscFunctionReturn(0);
918 }
919 
920 #undef __FUNCT__
921 #define __FUNCT__ "DMGetAggregates"
922 /*@
923    DMGetAggregates - Gets the aggregates that map between
924    grids associated with two DMs.
925 
926    Collective on DM
927 
928    Input Parameters:
929 +  dmc - the coarse grid DM
930 -  dmf - the fine grid DM
931 
932    Output Parameters:
933 .  rest - the restriction matrix (transpose of the projection matrix)
934 
935    Level: intermediate
936 
937 .keywords: interpolation, restriction, multigrid
938 
939 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation()
940 @*/
941 PetscErrorCode  DMGetAggregates(DM dmc, DM dmf, Mat *rest)
942 {
943   PetscErrorCode ierr;
944 
945   PetscFunctionBegin;
946   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
947   PetscFunctionReturn(0);
948 }
949 
950 #undef __FUNCT__
951 #define __FUNCT__ "DMSetApplicationContext"
952 /*@
953     DMSetApplicationContext - Set a user context into a DM object
954 
955     Not Collective
956 
957     Input Parameters:
958 +   dm - the DM object
959 -   ctx - the user context
960 
961     Level: intermediate
962 
963 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext()
964 
965 @*/
966 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
967 {
968   PetscFunctionBegin;
969   dm->ctx = ctx;
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "DMGetApplicationContext"
975 /*@
976     DMGetApplicationContext - Gets a user context from a DM object
977 
978     Not Collective
979 
980     Input Parameter:
981 .   dm - the DM object
982 
983     Output Parameter:
984 .   ctx - the user context
985 
986     Level: intermediate
987 
988 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext()
989 
990 @*/
991 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
992 {
993   PetscFunctionBegin;
994   *(void**)ctx = dm->ctx;
995   PetscFunctionReturn(0);
996 }
997 
998 #undef __FUNCT__
999 #define __FUNCT__ "DMSetInitialGuess"
1000 /*@
1001     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
1002 
1003     Logically Collective on DM
1004 
1005     Input Parameter:
1006 +   dm - the DM object to destroy
1007 -   f - the function to compute the initial guess
1008 
1009     Level: intermediate
1010 
1011 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1012 
1013 @*/
1014 PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
1015 {
1016   PetscFunctionBegin;
1017   dm->ops->initialguess = f;
1018   PetscFunctionReturn(0);
1019 }
1020 
1021 #undef __FUNCT__
1022 #define __FUNCT__ "DMSetFunction"
1023 /*@
1024     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
1025 
1026     Logically Collective on DM
1027 
1028     Input Parameter:
1029 +   dm - the DM object
1030 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
1031 
1032     Level: intermediate
1033 
1034     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
1035            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
1036 
1037 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1038          DMSetJacobian()
1039 
1040 @*/
1041 PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1042 {
1043   PetscFunctionBegin;
1044   dm->ops->function = f;
1045   if (f) {
1046     dm->ops->functionj = f;
1047   }
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 #undef __FUNCT__
1052 #define __FUNCT__ "DMSetJacobian"
1053 /*@
1054     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
1055 
1056     Logically Collective on DM
1057 
1058     Input Parameter:
1059 +   dm - the DM object to destroy
1060 -   f - the function to compute the matrix entries
1061 
1062     Level: intermediate
1063 
1064 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1065          DMSetFunction()
1066 
1067 @*/
1068 PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
1069 {
1070   PetscFunctionBegin;
1071   dm->ops->jacobian = f;
1072   PetscFunctionReturn(0);
1073 }
1074 
1075 #undef __FUNCT__
1076 #define __FUNCT__ "DMComputeInitialGuess"
1077 /*@
1078     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
1079 
1080     Collective on DM
1081 
1082     Input Parameter:
1083 +   dm - the DM object to destroy
1084 -   x - the vector to hold the initial guess values
1085 
1086     Level: developer
1087 
1088 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
1089 
1090 @*/
1091 PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
1092 {
1093   PetscErrorCode ierr;
1094   PetscFunctionBegin;
1095   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
1096   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 #undef __FUNCT__
1101 #define __FUNCT__ "DMHasInitialGuess"
1102 /*@
1103     DMHasInitialGuess - does the DM object have an initial guess function
1104 
1105     Not Collective
1106 
1107     Input Parameter:
1108 .   dm - the DM object to destroy
1109 
1110     Output Parameter:
1111 .   flg - PETSC_TRUE if function exists
1112 
1113     Level: developer
1114 
1115 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1116 
1117 @*/
1118 PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
1119 {
1120   PetscFunctionBegin;
1121   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 #undef __FUNCT__
1126 #define __FUNCT__ "DMHasFunction"
1127 /*@
1128     DMHasFunction - does the DM object have a function
1129 
1130     Not Collective
1131 
1132     Input Parameter:
1133 .   dm - the DM object to destroy
1134 
1135     Output Parameter:
1136 .   flg - PETSC_TRUE if function exists
1137 
1138     Level: developer
1139 
1140 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1141 
1142 @*/
1143 PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
1144 {
1145   PetscFunctionBegin;
1146   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 #undef __FUNCT__
1151 #define __FUNCT__ "DMHasJacobian"
1152 /*@
1153     DMHasJacobian - does the DM object have a matrix function
1154 
1155     Not Collective
1156 
1157     Input Parameter:
1158 .   dm - the DM object to destroy
1159 
1160     Output Parameter:
1161 .   flg - PETSC_TRUE if function exists
1162 
1163     Level: developer
1164 
1165 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1166 
1167 @*/
1168 PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
1169 {
1170   PetscFunctionBegin;
1171   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 #undef __FUNCT__
1176 #define __FUNCT__ "DMComputeFunction"
1177 /*@
1178     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
1179 
1180     Collective on DM
1181 
1182     Input Parameter:
1183 +   dm - the DM object to destroy
1184 .   x - the location where the function is evaluationed, may be ignored for linear problems
1185 -   b - the vector to hold the right hand side entries
1186 
1187     Level: developer
1188 
1189 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1190          DMSetJacobian()
1191 
1192 @*/
1193 PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
1194 {
1195   PetscErrorCode ierr;
1196   PetscFunctionBegin;
1197   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1198   PetscStackPush("DM user function");
1199   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1200   PetscStackPop;
1201   PetscFunctionReturn(0);
1202 }
1203 
1204 
1205 #undef __FUNCT__
1206 #define __FUNCT__ "DMComputeJacobian"
1207 /*@
1208     DMComputeJacobian - compute the matrix entries for the solver
1209 
1210     Collective on DM
1211 
1212     Input Parameter:
1213 +   dm - the DM object
1214 .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
1215 .   A - matrix that defines the operator for the linear solve
1216 -   B - the matrix used to construct the preconditioner
1217 
1218     Level: developer
1219 
1220 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1221          DMSetFunction()
1222 
1223 @*/
1224 PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1225 {
1226   PetscErrorCode ierr;
1227 
1228   PetscFunctionBegin;
1229   if (!dm->ops->jacobian) {
1230     ISColoring     coloring;
1231     MatFDColoring  fd;
1232 
1233     ierr = DMGetColoring(dm,IS_COLORING_GLOBAL,MATAIJ,&coloring);CHKERRQ(ierr);
1234     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1235     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1236     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1237 
1238     dm->fd = fd;
1239     dm->ops->jacobian = DMComputeJacobianDefault;
1240 
1241     /* don't know why this is needed */
1242     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1243   }
1244   if (!x) x = dm->x;
1245   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1246 
1247   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1248   if (x) {
1249     if (!dm->x) {
1250       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1251     }
1252     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1253   }
1254   if (A != B) {
1255     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1256     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1257   }
1258   PetscFunctionReturn(0);
1259 }
1260 
1261 
1262 PetscFList DMList                       = PETSC_NULL;
1263 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1264 
1265 #undef __FUNCT__
1266 #define __FUNCT__ "DMSetType"
1267 /*@C
1268   DMSetType - Builds a DM, for a particular DM implementation.
1269 
1270   Collective on DM
1271 
1272   Input Parameters:
1273 + dm     - The DM object
1274 - method - The name of the DM type
1275 
1276   Options Database Key:
1277 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1278 
1279   Notes:
1280   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1281 
1282   Level: intermediate
1283 
1284 .keywords: DM, set, type
1285 .seealso: DMGetType(), DMCreate()
1286 @*/
1287 PetscErrorCode  DMSetType(DM dm, const DMType method)
1288 {
1289   PetscErrorCode (*r)(DM);
1290   PetscBool      match;
1291   PetscErrorCode ierr;
1292 
1293   PetscFunctionBegin;
1294   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1295   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1296   if (match) PetscFunctionReturn(0);
1297 
1298   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1299   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1300   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1301 
1302   if (dm->ops->destroy) {
1303     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1304   }
1305   ierr = (*r)(dm);CHKERRQ(ierr);
1306   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1307   PetscFunctionReturn(0);
1308 }
1309 
1310 #undef __FUNCT__
1311 #define __FUNCT__ "DMGetType"
1312 /*@C
1313   DMGetType - Gets the DM type name (as a string) from the DM.
1314 
1315   Not Collective
1316 
1317   Input Parameter:
1318 . dm  - The DM
1319 
1320   Output Parameter:
1321 . type - The DM type name
1322 
1323   Level: intermediate
1324 
1325 .keywords: DM, get, type, name
1326 .seealso: DMSetType(), DMCreate()
1327 @*/
1328 PetscErrorCode  DMGetType(DM dm, const DMType *type)
1329 {
1330   PetscErrorCode ierr;
1331 
1332   PetscFunctionBegin;
1333   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1334   PetscValidCharPointer(type,2);
1335   if (!DMRegisterAllCalled) {
1336     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1337   }
1338   *type = ((PetscObject)dm)->type_name;
1339   PetscFunctionReturn(0);
1340 }
1341 
1342 #undef __FUNCT__
1343 #define __FUNCT__ "DMConvert"
1344 /*@C
1345   DMConvert - Converts a DM to another DM, either of the same or different type.
1346 
1347   Collective on DM
1348 
1349   Input Parameters:
1350 + dm - the DM
1351 - newtype - new DM type (use "same" for the same type)
1352 
1353   Output Parameter:
1354 . M - pointer to new DM
1355 
1356   Notes:
1357   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
1358   the MPI communicator of the generated DM is always the same as the communicator
1359   of the input DM.
1360 
1361   Level: intermediate
1362 
1363 .seealso: DMCreate()
1364 @*/
1365 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
1366 {
1367   DM             B;
1368   char           convname[256];
1369   PetscBool      sametype, issame;
1370   PetscErrorCode ierr;
1371 
1372   PetscFunctionBegin;
1373   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1374   PetscValidType(dm,1);
1375   PetscValidPointer(M,3);
1376   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
1377   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
1378   {
1379     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
1380 
1381     /*
1382        Order of precedence:
1383        1) See if a specialized converter is known to the current DM.
1384        2) See if a specialized converter is known to the desired DM class.
1385        3) See if a good general converter is registered for the desired class
1386        4) See if a good general converter is known for the current matrix.
1387        5) Use a really basic converter.
1388     */
1389 
1390     /* 1) See if a specialized converter is known to the current DM and the desired class */
1391     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1392     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1393     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1394     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1395     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1396     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1397     if (conv) goto foundconv;
1398 
1399     /* 2)  See if a specialized converter is known to the desired DM class. */
1400     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
1401     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
1402     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1403     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1404     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1405     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1406     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1407     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1408     if (conv) {
1409       ierr = DMDestroy(&B);CHKERRQ(ierr);
1410       goto foundconv;
1411     }
1412 
1413 #if 0
1414     /* 3) See if a good general converter is registered for the desired class */
1415     conv = B->ops->convertfrom;
1416     ierr = DMDestroy(&B);CHKERRQ(ierr);
1417     if (conv) goto foundconv;
1418 
1419     /* 4) See if a good general converter is known for the current matrix */
1420     if (dm->ops->convert) {
1421       conv = dm->ops->convert;
1422     }
1423     if (conv) goto foundconv;
1424 #endif
1425 
1426     /* 5) Use a really basic converter. */
1427     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
1428 
1429     foundconv:
1430     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1431     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
1432     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1433   }
1434   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
1435   PetscFunctionReturn(0);
1436 }
1437 
1438 /*--------------------------------------------------------------------------------------------------------------------*/
1439 
1440 #undef __FUNCT__
1441 #define __FUNCT__ "DMRegister"
1442 /*@C
1443   DMRegister - See DMRegisterDynamic()
1444 
1445   Level: advanced
1446 @*/
1447 PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1448 {
1449   char fullname[PETSC_MAX_PATH_LEN];
1450   PetscErrorCode ierr;
1451 
1452   PetscFunctionBegin;
1453   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1454   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1455   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1456   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1457   PetscFunctionReturn(0);
1458 }
1459 
1460 
1461 /*--------------------------------------------------------------------------------------------------------------------*/
1462 #undef __FUNCT__
1463 #define __FUNCT__ "DMRegisterDestroy"
1464 /*@C
1465    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1466 
1467    Not Collective
1468 
1469    Level: advanced
1470 
1471 .keywords: DM, register, destroy
1472 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1473 @*/
1474 PetscErrorCode  DMRegisterDestroy(void)
1475 {
1476   PetscErrorCode ierr;
1477 
1478   PetscFunctionBegin;
1479   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1480   DMRegisterAllCalled = PETSC_FALSE;
1481   PetscFunctionReturn(0);
1482 }
1483 
1484 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1485 #include <mex.h>
1486 
1487 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
1488 
1489 #undef __FUNCT__
1490 #define __FUNCT__ "DMComputeFunction_Matlab"
1491 /*
1492    DMComputeFunction_Matlab - Calls the function that has been set with
1493                          DMSetFunctionMatlab().
1494 
1495    For linear problems x is null
1496 
1497 .seealso: DMSetFunction(), DMGetFunction()
1498 */
1499 PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
1500 {
1501   PetscErrorCode    ierr;
1502   DMMatlabContext   *sctx;
1503   int               nlhs = 1,nrhs = 4;
1504   mxArray	    *plhs[1],*prhs[4];
1505   long long int     lx = 0,ly = 0,ls = 0;
1506 
1507   PetscFunctionBegin;
1508   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1509   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1510   PetscCheckSameComm(dm,1,y,3);
1511 
1512   /* call Matlab function in ctx with arguments x and y */
1513   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1514   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1515   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1516   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
1517   prhs[0] =  mxCreateDoubleScalar((double)ls);
1518   prhs[1] =  mxCreateDoubleScalar((double)lx);
1519   prhs[2] =  mxCreateDoubleScalar((double)ly);
1520   prhs[3] =  mxCreateString(sctx->funcname);
1521   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
1522   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
1523   mxDestroyArray(prhs[0]);
1524   mxDestroyArray(prhs[1]);
1525   mxDestroyArray(prhs[2]);
1526   mxDestroyArray(prhs[3]);
1527   mxDestroyArray(plhs[0]);
1528   PetscFunctionReturn(0);
1529 }
1530 
1531 
1532 #undef __FUNCT__
1533 #define __FUNCT__ "DMSetFunctionMatlab"
1534 /*
1535    DMSetFunctionMatlab - Sets the function evaluation routine
1536 
1537 */
1538 PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
1539 {
1540   PetscErrorCode    ierr;
1541   DMMatlabContext   *sctx;
1542 
1543   PetscFunctionBegin;
1544   /* currently sctx is memory bleed */
1545   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1546   if (!sctx) {
1547     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1548   }
1549   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
1550   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1551   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
1552   PetscFunctionReturn(0);
1553 }
1554 
1555 #undef __FUNCT__
1556 #define __FUNCT__ "DMComputeJacobian_Matlab"
1557 /*
1558    DMComputeJacobian_Matlab - Calls the function that has been set with
1559                          DMSetJacobianMatlab().
1560 
1561    For linear problems x is null
1562 
1563 .seealso: DMSetFunction(), DMGetFunction()
1564 */
1565 PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
1566 {
1567   PetscErrorCode    ierr;
1568   DMMatlabContext   *sctx;
1569   int               nlhs = 2,nrhs = 5;
1570   mxArray	    *plhs[2],*prhs[5];
1571   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
1572 
1573   PetscFunctionBegin;
1574   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1575   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
1576 
1577   /* call MATLAB function in ctx with arguments x, A, and B */
1578   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1579   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1580   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1581   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
1582   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
1583   prhs[0] =  mxCreateDoubleScalar((double)ls);
1584   prhs[1] =  mxCreateDoubleScalar((double)lx);
1585   prhs[2] =  mxCreateDoubleScalar((double)lA);
1586   prhs[3] =  mxCreateDoubleScalar((double)lB);
1587   prhs[4] =  mxCreateString(sctx->jacname);
1588   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
1589   *str    =  (MatStructure) mxGetScalar(plhs[0]);
1590   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
1591   mxDestroyArray(prhs[0]);
1592   mxDestroyArray(prhs[1]);
1593   mxDestroyArray(prhs[2]);
1594   mxDestroyArray(prhs[3]);
1595   mxDestroyArray(prhs[4]);
1596   mxDestroyArray(plhs[0]);
1597   mxDestroyArray(plhs[1]);
1598   PetscFunctionReturn(0);
1599 }
1600 
1601 
1602 #undef __FUNCT__
1603 #define __FUNCT__ "DMSetJacobianMatlab"
1604 /*
1605    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
1606 
1607 */
1608 PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
1609 {
1610   PetscErrorCode    ierr;
1611   DMMatlabContext   *sctx;
1612 
1613   PetscFunctionBegin;
1614   /* currently sctx is memory bleed */
1615   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1616   if (!sctx) {
1617     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1618   }
1619   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
1620   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1621   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
1622   PetscFunctionReturn(0);
1623 }
1624 #endif
1625 
1626 #undef __FUNCT__
1627 #define __FUNCT__ "DMLoad"
1628 /*@C
1629   DMLoad - Loads a DM that has been stored in binary or HDF5 format
1630   with DMView().
1631 
1632   Collective on PetscViewer
1633 
1634   Input Parameters:
1635 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
1636            some related function before a call to DMLoad().
1637 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
1638            HDF5 file viewer, obtained from PetscViewerHDF5Open()
1639 
1640    Level: intermediate
1641 
1642   Notes:
1643   Defaults to the DM DA.
1644 
1645   Notes for advanced users:
1646   Most users should not need to know the details of the binary storage
1647   format, since DMLoad() and DMView() completely hide these details.
1648   But for anyone who's interested, the standard binary matrix storage
1649   format is
1650 .vb
1651      has not yet been determined
1652 .ve
1653 
1654    In addition, PETSc automatically does the byte swapping for
1655 machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
1656 linux, Windows and the paragon; thus if you write your own binary
1657 read/write routines you have to swap the bytes; see PetscBinaryRead()
1658 and PetscBinaryWrite() to see how this may be done.
1659 
1660   Concepts: vector^loading from file
1661 
1662 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
1663 @*/
1664 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
1665 {
1666   PetscErrorCode ierr;
1667 
1668   PetscFunctionBegin;
1669   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
1670   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1671 
1672   if (!((PetscObject)newdm)->type_name) {
1673     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
1674   }
1675   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1676   PetscFunctionReturn(0);
1677 }
1678 
1679