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