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