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