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