xref: /petsc/src/dm/impls/shell/dmshell.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
1 #include <petscdmshell.h>       /*I    "petscdmshell.h"  I*/
2 #include <petscmat.h>
3 #include <petsc/private/dmimpl.h>
4 
5 typedef struct  {
6   Vec        Xglobal;
7   Vec        Xlocal;
8   Mat        A;
9   VecScatter gtol;
10   VecScatter ltog;
11   VecScatter ltol;
12   void       *ctx;
13 } DM_Shell;
14 
15 /*@
16    DMGlobalToLocalBeginDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to begin a global to local scatter
17    Collective
18 
19    Input Arguments:
20 +  dm - shell DM
21 .  g - global vector
22 .  mode - InsertMode
23 -  l - local vector
24 
25    Level: advanced
26 
27    Note:  This is not normally called directly by user code, generally user code calls DMGlobalToLocalBegin() and DMGlobalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function.
28 
29 .seealso: DMGlobalToLocalEndDefaultShell()
30 @*/
31 PetscErrorCode DMGlobalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
32 {
33   PetscErrorCode ierr;
34   DM_Shell       *shell = (DM_Shell*)dm->data;
35 
36   PetscFunctionBegin;
37   if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
38   ierr = VecScatterBegin(shell->gtol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
39   PetscFunctionReturn(0);
40 }
41 
42 /*@
43    DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to end a global to local scatter
44    Collective
45 
46    Input Arguments:
47 +  dm - shell DM
48 .  g - global vector
49 .  mode - InsertMode
50 -  l - local vector
51 
52    Level: advanced
53 
54 .seealso: DMGlobalToLocalBeginDefaultShell()
55 @*/
56 PetscErrorCode DMGlobalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
57 {
58   PetscErrorCode ierr;
59   DM_Shell       *shell = (DM_Shell*)dm->data;
60 
61   PetscFunctionBegin;
62    if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
63   ierr = VecScatterEnd(shell->gtol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
64   PetscFunctionReturn(0);
65 }
66 
67 /*@
68    DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to begin a local to global scatter
69    Collective
70 
71    Input Arguments:
72 +  dm - shell DM
73 .  l - local vector
74 .  mode - InsertMode
75 -  g - global vector
76 
77    Level: advanced
78 
79    Note:  This is not normally called directly by user code, generally user code calls DMLocalToGlobalBegin() and DMLocalToGlobalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function.
80 
81 .seealso: DMLocalToGlobalEndDefaultShell()
82 @*/
83 PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm,Vec l,InsertMode mode,Vec g)
84 {
85   PetscErrorCode ierr;
86   DM_Shell       *shell = (DM_Shell*)dm->data;
87 
88   PetscFunctionBegin;
89   if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
90   ierr = VecScatterBegin(shell->ltog,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr);
91   PetscFunctionReturn(0);
92 }
93 
94 /*@
95    DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to end a local to global scatter
96    Collective
97 
98    Input Arguments:
99 +  dm - shell DM
100 .  l - local vector
101 .  mode - InsertMode
102 -  g - global vector
103 
104    Level: advanced
105 
106 .seealso: DMLocalToGlobalBeginDefaultShell()
107 @*/
108 PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm,Vec l,InsertMode mode,Vec g)
109 {
110   PetscErrorCode ierr;
111   DM_Shell       *shell = (DM_Shell*)dm->data;
112 
113   PetscFunctionBegin;
114    if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
115   ierr = VecScatterEnd(shell->ltog,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118 
119 /*@
120    DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal VecScatter context set by the user to begin a local to local scatter
121    Collective
122 
123    Input Arguments:
124 +  dm - shell DM
125 .  g - the original local vector
126 -  mode - InsertMode
127 
128    Output Parameter:
129 .  l  - the local vector with correct ghost values
130 
131    Level: advanced
132 
133    Note:  This is not normally called directly by user code, generally user code calls DMLocalToLocalBegin() and DMLocalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToLocal() then those routines might have reason to call this function.
134 
135 .seealso: DMLocalToLocalEndDefaultShell()
136 @*/
137 PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
138 {
139   PetscErrorCode ierr;
140   DM_Shell       *shell = (DM_Shell*)dm->data;
141 
142   PetscFunctionBegin;
143   if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()");
144   ierr = VecScatterBegin(shell->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 /*@
149    DMLocalToLocalEndDefaultShell - Uses the LocalToLocal VecScatter context set by the user to end a local to local scatter
150    Collective
151 
152    Input Arguments:
153 +  dm - shell DM
154 .  g - the original local vector
155 -  mode - InsertMode
156 
157    Output Parameter:
158 .  l  - the local vector with correct ghost values
159 
160    Level: advanced
161 
162 .seealso: DMLocalToLocalBeginDefaultShell()
163 @*/
164 PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
165 {
166   PetscErrorCode ierr;
167   DM_Shell       *shell = (DM_Shell*)dm->data;
168 
169   PetscFunctionBegin;
170    if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
171   ierr = VecScatterEnd(shell->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
172   PetscFunctionReturn(0);
173 }
174 
175 static PetscErrorCode DMCreateMatrix_Shell(DM dm,Mat *J)
176 {
177   PetscErrorCode ierr;
178   DM_Shell       *shell = (DM_Shell*)dm->data;
179   Mat            A;
180 
181   PetscFunctionBegin;
182   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
183   PetscValidPointer(J,3);
184   if (!shell->A) {
185     if (shell->Xglobal) {
186       PetscInt m,M;
187       ierr = PetscInfo(dm,"Naively creating matrix using global vector distribution without preallocation\n");CHKERRQ(ierr);
188       ierr = VecGetSize(shell->Xglobal,&M);CHKERRQ(ierr);
189       ierr = VecGetLocalSize(shell->Xglobal,&m);CHKERRQ(ierr);
190       ierr = MatCreate(PetscObjectComm((PetscObject)dm),&shell->A);CHKERRQ(ierr);
191       ierr = MatSetSizes(shell->A,m,m,M,M);CHKERRQ(ierr);
192       ierr = MatSetType(shell->A,dm->mattype);CHKERRQ(ierr);
193       ierr = MatSetUp(shell->A);CHKERRQ(ierr);
194     } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector");
195   }
196   A = shell->A;
197   ierr = MatDuplicate(A,MAT_SHARE_NONZERO_PATTERN,J);CHKERRQ(ierr);
198   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
199   PetscFunctionReturn(0);
200 }
201 
202 PetscErrorCode DMCreateGlobalVector_Shell(DM dm,Vec *gvec)
203 {
204   PetscErrorCode ierr;
205   DM_Shell       *shell = (DM_Shell*)dm->data;
206   Vec            X;
207 
208   PetscFunctionBegin;
209   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
210   PetscValidPointer(gvec,2);
211   *gvec = NULL;
212   X     = shell->Xglobal;
213   if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()");
214   /* Need to create a copy in order to attach the DM to the vector */
215   ierr = VecDuplicate(X,gvec);CHKERRQ(ierr);
216   ierr = VecZeroEntries(*gvec);CHKERRQ(ierr);
217   ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr);
218   PetscFunctionReturn(0);
219 }
220 
221 PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec)
222 {
223   PetscErrorCode ierr;
224   DM_Shell       *shell = (DM_Shell*)dm->data;
225   Vec            X;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
229   PetscValidPointer(gvec,2);
230   *gvec = NULL;
231   X     = shell->Xlocal;
232   if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()");
233   /* Need to create a copy in order to attach the DM to the vector */
234   ierr = VecDuplicate(X,gvec);CHKERRQ(ierr);
235   ierr = VecZeroEntries(*gvec);CHKERRQ(ierr);
236   ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr);
237   PetscFunctionReturn(0);
238 }
239 
240 /*@
241    DMShellSetContext - set some data to be usable by this DM
242 
243    Collective
244 
245    Input Arguments:
246 +  dm - shell DM
247 -  ctx - the context
248 
249    Level: advanced
250 
251 .seealso: DMCreateMatrix(), DMShellGetContext()
252 @*/
253 PetscErrorCode DMShellSetContext(DM dm,void *ctx)
254 {
255   DM_Shell       *shell = (DM_Shell*)dm->data;
256   PetscErrorCode ierr;
257   PetscBool      isshell;
258 
259   PetscFunctionBegin;
260   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
261   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
262   if (!isshell) PetscFunctionReturn(0);
263   shell->ctx = ctx;
264   PetscFunctionReturn(0);
265 }
266 
267 /*@
268    DMShellGetContext - Returns the user-provided context associated to the DM
269 
270    Collective
271 
272    Input Argument:
273 .  dm - shell DM
274 
275    Output Argument:
276 .  ctx - the context
277 
278    Level: advanced
279 
280 .seealso: DMCreateMatrix(), DMShellSetContext()
281 @*/
282 PetscErrorCode DMShellGetContext(DM dm,void **ctx)
283 {
284   DM_Shell       *shell = (DM_Shell*)dm->data;
285   PetscErrorCode ierr;
286   PetscBool      isshell;
287 
288   PetscFunctionBegin;
289   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
290   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
291   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
292   *ctx = shell->ctx;
293   PetscFunctionReturn(0);
294 }
295 
296 /*@
297    DMShellSetMatrix - sets a template matrix associated with the DMShell
298 
299    Collective
300 
301    Input Arguments:
302 +  dm - shell DM
303 -  J - template matrix
304 
305    Level: advanced
306 
307    Developer Notes:
308     To avoid circular references, if J is already associated to the same DM, then MatDuplicate(SHARE_NONZERO_PATTERN) is called, followed by removing the DM reference from the private template.
309 
310 .seealso: DMCreateMatrix(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
311 @*/
312 PetscErrorCode DMShellSetMatrix(DM dm,Mat J)
313 {
314   DM_Shell       *shell = (DM_Shell*)dm->data;
315   PetscErrorCode ierr;
316   PetscBool      isshell;
317   DM             mdm;
318 
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
321   PetscValidHeaderSpecific(J,MAT_CLASSID,2);
322   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
323   if (!isshell) PetscFunctionReturn(0);
324   if (J == shell->A) PetscFunctionReturn(0);
325   ierr = MatGetDM(J,&mdm);CHKERRQ(ierr);
326   ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
327   ierr = MatDestroy(&shell->A);CHKERRQ(ierr);
328   if (mdm == dm) {
329     ierr = MatDuplicate(J,MAT_SHARE_NONZERO_PATTERN,&shell->A);CHKERRQ(ierr);
330     ierr = MatSetDM(shell->A,NULL);CHKERRQ(ierr);
331   } else shell->A = J;
332   PetscFunctionReturn(0);
333 }
334 
335 /*@C
336    DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM
337 
338    Logically Collective on dm
339 
340    Input Arguments:
341 +  dm - the shell DM
342 -  func - the function to create a matrix
343 
344    Level: advanced
345 
346 .seealso: DMCreateMatrix(), DMShellSetMatrix(), DMShellSetContext(), DMShellGetContext()
347 @*/
348 PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*))
349 {
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
352   dm->ops->creatematrix = func;
353   PetscFunctionReturn(0);
354 }
355 
356 /*@
357    DMShellSetGlobalVector - sets a template global vector associated with the DMShell
358 
359    Logically Collective on dm
360 
361    Input Arguments:
362 +  dm - shell DM
363 -  X - template vector
364 
365    Level: advanced
366 
367 .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector()
368 @*/
369 PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X)
370 {
371   DM_Shell       *shell = (DM_Shell*)dm->data;
372   PetscErrorCode ierr;
373   PetscBool      isshell;
374   DM             vdm;
375 
376   PetscFunctionBegin;
377   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
378   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
379   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
380   if (!isshell) PetscFunctionReturn(0);
381   ierr = VecGetDM(X,&vdm);CHKERRQ(ierr);
382   /*
383       if the vector proposed as the new base global vector for the DM is a DM vector associated
384       with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
385       we get a circular dependency that prevents the DM from being destroy when it should be.
386       This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
387       DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
388       to set its input vector (which is associated with the DM) as the base global vector.
389       Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
390       for pointing out the problem.
391    */
392   if (vdm == dm) PetscFunctionReturn(0);
393   ierr           = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
394   ierr           = VecDestroy(&shell->Xglobal);CHKERRQ(ierr);
395   shell->Xglobal = X;
396   PetscFunctionReturn(0);
397 }
398 
399 /*@C
400    DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM
401 
402    Logically Collective
403 
404    Input Arguments:
405 +  dm - the shell DM
406 -  func - the creation routine
407 
408    Level: advanced
409 
410 .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
411 @*/
412 PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
413 {
414   PetscFunctionBegin;
415   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
416   dm->ops->createglobalvector = func;
417   PetscFunctionReturn(0);
418 }
419 
420 /*@
421    DMShellSetLocalVector - sets a template local vector associated with the DMShell
422 
423    Logically Collective on dm
424 
425    Input Arguments:
426 +  dm - shell DM
427 -  X - template vector
428 
429    Level: advanced
430 
431 .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector()
432 @*/
433 PetscErrorCode DMShellSetLocalVector(DM dm,Vec X)
434 {
435   DM_Shell       *shell = (DM_Shell*)dm->data;
436   PetscErrorCode ierr;
437   PetscBool      isshell;
438   DM             vdm;
439 
440   PetscFunctionBegin;
441   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
442   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
443   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
444   if (!isshell) PetscFunctionReturn(0);
445   ierr = VecGetDM(X,&vdm);CHKERRQ(ierr);
446   /*
447       if the vector proposed as the new base global vector for the DM is a DM vector associated
448       with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
449       we get a circular dependency that prevents the DM from being destroy when it should be.
450       This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
451       DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
452       to set its input vector (which is associated with the DM) as the base global vector.
453       Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
454       for pointing out the problem.
455    */
456   if (vdm == dm) PetscFunctionReturn(0);
457   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
458   ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr);
459   shell->Xlocal = X;
460   PetscFunctionReturn(0);
461 }
462 
463 /*@C
464    DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM
465 
466    Logically Collective
467 
468    Input Arguments:
469 +  dm - the shell DM
470 -  func - the creation routine
471 
472    Level: advanced
473 
474 .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
475 @*/
476 PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
477 {
478   PetscFunctionBegin;
479   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
480   dm->ops->createlocalvector = func;
481   PetscFunctionReturn(0);
482 }
483 
484 /*@C
485    DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter
486 
487    Logically Collective on dm
488 
489    Input Arguments
490 +  dm - the shell DM
491 .  begin - the routine that begins the global to local scatter
492 -  end - the routine that ends the global to local scatter
493 
494    Notes:
495     If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then
496    DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers
497 
498    Level: advanced
499 
500 .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
501 @*/
502 PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
503   PetscFunctionBegin;
504   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
505   dm->ops->globaltolocalbegin = begin;
506   dm->ops->globaltolocalend = end;
507   PetscFunctionReturn(0);
508 }
509 
510 /*@C
511    DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter
512 
513    Logically Collective on dm
514 
515    Input Arguments
516 +  dm - the shell DM
517 .  begin - the routine that begins the local to global scatter
518 -  end - the routine that ends the local to global scatter
519 
520    Notes:
521     If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then
522    DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers
523 
524    Level: advanced
525 
526 .seealso: DMShellSetGlobalToLocal()
527 @*/
528 PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
531   dm->ops->localtoglobalbegin = begin;
532   dm->ops->localtoglobalend = end;
533   PetscFunctionReturn(0);
534 }
535 
536 /*@C
537    DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter
538 
539    Logically Collective on dm
540 
541    Input Arguments
542 +  dm - the shell DM
543 .  begin - the routine that begins the local to local scatter
544 -  end - the routine that ends the local to local scatter
545 
546    Notes:
547     If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then
548    DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers
549 
550    Level: advanced
551 
552 .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
553 @*/
554 PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
555   PetscFunctionBegin;
556   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
557   dm->ops->localtolocalbegin = begin;
558   dm->ops->localtolocalend = end;
559   PetscFunctionReturn(0);
560 }
561 
562 /*@
563    DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication
564 
565    Logically Collective on dm
566 
567    Input Arguments
568 +  dm - the shell DM
569 -  gtol - the global to local VecScatter context
570 
571    Level: advanced
572 
573 .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
574 @*/
575 PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
576 {
577   DM_Shell       *shell = (DM_Shell*)dm->data;
578   PetscErrorCode ierr;
579 
580   PetscFunctionBegin;
581   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
582   PetscValidHeaderSpecific(gtol,PETSCSF_CLASSID,2);
583   ierr = PetscObjectReference((PetscObject)gtol);CHKERRQ(ierr);
584   ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr);
585   shell->gtol = gtol;
586   PetscFunctionReturn(0);
587 }
588 
589 /*@
590    DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication
591 
592    Logically Collective on dm
593 
594    Input Arguments
595 +  dm - the shell DM
596 -  ltog - the local to global VecScatter context
597 
598    Level: advanced
599 
600 .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell()
601 @*/
602 PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
603 {
604   DM_Shell       *shell = (DM_Shell*)dm->data;
605   PetscErrorCode ierr;
606 
607   PetscFunctionBegin;
608   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
609   PetscValidHeaderSpecific(ltog,PETSCSF_CLASSID,2);
610   ierr = PetscObjectReference((PetscObject)ltog);CHKERRQ(ierr);
611   ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr);
612   shell->ltog = ltog;
613   PetscFunctionReturn(0);
614 }
615 
616 /*@
617    DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication
618 
619    Logically Collective on dm
620 
621    Input Arguments
622 +  dm - the shell DM
623 -  ltol - the local to local VecScatter context
624 
625    Level: advanced
626 
627 .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
628 @*/
629 PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
630 {
631   DM_Shell       *shell = (DM_Shell*)dm->data;
632   PetscErrorCode ierr;
633 
634   PetscFunctionBegin;
635   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
636   PetscValidHeaderSpecific(ltol,PETSCSF_CLASSID,2);
637   ierr = PetscObjectReference((PetscObject)ltol);CHKERRQ(ierr);
638   ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr);
639   shell->ltol = ltol;
640   PetscFunctionReturn(0);
641 }
642 
643 /*@C
644    DMShellSetCoarsen - Set the routine used to coarsen the shell DM
645 
646    Logically Collective on dm
647 
648    Input Arguments
649 +  dm - the shell DM
650 -  coarsen - the routine that coarsens the DM
651 
652    Level: advanced
653 
654 .seealso: DMShellSetRefine(), DMCoarsen(), DMShellGetCoarsen(), DMShellSetContext(), DMShellGetContext()
655 @*/
656 PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*))
657 {
658   PetscErrorCode ierr;
659   PetscBool      isshell;
660 
661   PetscFunctionBegin;
662   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
663   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
664   if (!isshell) PetscFunctionReturn(0);
665   dm->ops->coarsen = coarsen;
666   PetscFunctionReturn(0);
667 }
668 
669 /*@C
670    DMShellGetCoarsen - Get the routine used to coarsen the shell DM
671 
672    Logically Collective on dm
673 
674    Input Argument:
675 .  dm - the shell DM
676 
677    Output Argument:
678 .  coarsen - the routine that coarsens the DM
679 
680    Level: advanced
681 
682 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
683 @*/
684 PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM,MPI_Comm,DM*))
685 {
686   PetscErrorCode ierr;
687   PetscBool      isshell;
688 
689   PetscFunctionBegin;
690   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
691   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
692   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
693   *coarsen = dm->ops->coarsen;
694   PetscFunctionReturn(0);
695 }
696 
697 /*@C
698    DMShellSetRefine - Set the routine used to refine the shell DM
699 
700    Logically Collective on dm
701 
702    Input Arguments
703 +  dm - the shell DM
704 -  refine - the routine that refines the DM
705 
706    Level: advanced
707 
708 .seealso: DMShellSetCoarsen(), DMRefine(), DMShellGetRefine(), DMShellSetContext(), DMShellGetContext()
709 @*/
710 PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*))
711 {
712   PetscErrorCode ierr;
713   PetscBool      isshell;
714 
715   PetscFunctionBegin;
716   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
717   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
718   if (!isshell) PetscFunctionReturn(0);
719   dm->ops->refine = refine;
720   PetscFunctionReturn(0);
721 }
722 
723 /*@C
724    DMShellGetRefine - Get the routine used to refine the shell DM
725 
726    Logically Collective on dm
727 
728    Input Argument:
729 .  dm - the shell DM
730 
731    Output Argument:
732 .  refine - the routine that refines the DM
733 
734    Level: advanced
735 
736 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
737 @*/
738 PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM,MPI_Comm,DM*))
739 {
740   PetscErrorCode ierr;
741   PetscBool      isshell;
742 
743   PetscFunctionBegin;
744   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
745   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
746   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
747   *refine = dm->ops->refine;
748   PetscFunctionReturn(0);
749 }
750 
751 /*@C
752    DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator
753 
754    Logically Collective on dm
755 
756    Input Arguments
757 +  dm - the shell DM
758 -  interp - the routine to create the interpolation
759 
760    Level: advanced
761 
762 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
763 @*/
764 PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*))
765 {
766   PetscErrorCode ierr;
767   PetscBool      isshell;
768 
769   PetscFunctionBegin;
770   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
771   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
772   if (!isshell) PetscFunctionReturn(0);
773   dm->ops->createinterpolation = interp;
774   PetscFunctionReturn(0);
775 }
776 
777 /*@C
778    DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator
779 
780    Logically Collective on dm
781 
782    Input Argument:
783 +  dm - the shell DM
784 
785    Output Argument:
786 -  interp - the routine to create the interpolation
787 
788    Level: advanced
789 
790 .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
791 @*/
792 PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM,DM,Mat*,Vec*))
793 {
794   PetscErrorCode ierr;
795   PetscBool      isshell;
796 
797   PetscFunctionBegin;
798   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
799   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
800   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
801   *interp = dm->ops->createinterpolation;
802   PetscFunctionReturn(0);
803 }
804 
805 /*@C
806    DMShellSetCreateRestriction - Set the routine used to create the restriction operator
807 
808    Logically Collective on dm
809 
810    Input Arguments
811 +  dm - the shell DM
812 -  striction- the routine to create the restriction
813 
814    Level: advanced
815 
816 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
817 @*/
818 PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*))
819 {
820   PetscErrorCode ierr;
821   PetscBool      isshell;
822 
823   PetscFunctionBegin;
824   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
825   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
826   if (!isshell) PetscFunctionReturn(0);
827   dm->ops->createrestriction = restriction;
828   PetscFunctionReturn(0);
829 }
830 
831 /*@C
832    DMShellGetCreateRestriction - Get the routine used to create the restriction operator
833 
834    Logically Collective on dm
835 
836    Input Argument:
837 +  dm - the shell DM
838 
839    Output Argument:
840 -  restriction - the routine to create the restriction
841 
842    Level: advanced
843 
844 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext()
845 @*/
846 PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM,DM,Mat*))
847 {
848   PetscErrorCode ierr;
849   PetscBool      isshell;
850 
851   PetscFunctionBegin;
852   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
853   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
854   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
855   *restriction = dm->ops->createrestriction;
856   PetscFunctionReturn(0);
857 }
858 
859 /*@C
860    DMShellSetCreateInjection - Set the routine used to create the injection operator
861 
862    Logically Collective on dm
863 
864    Input Arguments
865 +  dm - the shell DM
866 -  inject - the routine to create the injection
867 
868    Level: advanced
869 
870 .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellGetCreateInjection(), DMShellSetContext(), DMShellGetContext()
871 @*/
872 PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*))
873 {
874   PetscErrorCode ierr;
875   PetscBool      isshell;
876 
877   PetscFunctionBegin;
878   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
879   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
880   if (!isshell) PetscFunctionReturn(0);
881   dm->ops->createinjection = inject;
882   PetscFunctionReturn(0);
883 }
884 
885 /*@C
886    DMShellGetCreateInjection - Get the routine used to create the injection operator
887 
888    Logically Collective on dm
889 
890    Input Argument:
891 +  dm - the shell DM
892 
893    Output Argument:
894 -  inject - the routine to create the injection
895 
896    Level: advanced
897 
898 .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext()
899 @*/
900 PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM,DM,Mat*))
901 {
902   PetscErrorCode ierr;
903   PetscBool      isshell;
904 
905   PetscFunctionBegin;
906   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
907   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
908   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
909   *inject = dm->ops->createinjection;
910   PetscFunctionReturn(0);
911 }
912 
913 /*@C
914    DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM
915 
916    Logically Collective on dm
917 
918    Input Arguments
919 +  dm - the shell DM
920 -  decomp - the routine to create the decomposition
921 
922    Level: advanced
923 
924 .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext()
925 @*/
926 PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**))
927 {
928   PetscErrorCode ierr;
929   PetscBool      isshell;
930 
931   PetscFunctionBegin;
932   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
933   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
934   if (!isshell) PetscFunctionReturn(0);
935   dm->ops->createfielddecomposition = decomp;
936   PetscFunctionReturn(0);
937 }
938 
939 /*@C
940    DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM
941 
942    Logically Collective on dm
943 
944    Input Arguments
945 +  dm - the shell DM
946 -  decomp - the routine to create the decomposition
947 
948    Level: advanced
949 
950 .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext()
951 @*/
952 PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**))
953 {
954   PetscErrorCode ierr;
955   PetscBool      isshell;
956 
957   PetscFunctionBegin;
958   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
959   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
960   if (!isshell) PetscFunctionReturn(0);
961   dm->ops->createdomaindecomposition = decomp;
962   PetscFunctionReturn(0);
963 }
964 
965 /*@C
966    DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM
967 
968    Logically Collective on dm
969 
970    Input Arguments
971 +  dm - the shell DM
972 -  scatter - the routine to create the scatters
973 
974    Level: advanced
975 
976 .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext()
977 @*/
978 PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**))
979 {
980   PetscErrorCode ierr;
981   PetscBool      isshell;
982 
983   PetscFunctionBegin;
984   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
985   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
986   if (!isshell) PetscFunctionReturn(0);
987   dm->ops->createddscatters = scatter;
988   PetscFunctionReturn(0);
989 }
990 
991 /*@C
992    DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM
993 
994    Logically Collective on dm
995 
996    Input Arguments
997 +  dm - the shell DM
998 -  subdm - the routine to create the decomposition
999 
1000    Level: advanced
1001 
1002 .seealso: DMCreateSubDM(), DMShellGetCreateSubDM(), DMShellSetContext(), DMShellGetContext()
1003 @*/
1004 PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
1005 {
1006   PetscErrorCode ierr;
1007   PetscBool      isshell;
1008 
1009   PetscFunctionBegin;
1010   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1011   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
1012   if (!isshell) PetscFunctionReturn(0);
1013   dm->ops->createsubdm = subdm;
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 /*@C
1018    DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM
1019 
1020    Logically Collective on dm
1021 
1022    Input Argument:
1023 +  dm - the shell DM
1024 
1025    Output Argument:
1026 -  subdm - the routine to create the decomposition
1027 
1028    Level: advanced
1029 
1030 .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext()
1031 @*/
1032 PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
1033 {
1034   PetscErrorCode ierr;
1035   PetscBool      isshell;
1036 
1037   PetscFunctionBegin;
1038   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1039   ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr);
1040   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
1041   *subdm = dm->ops->createsubdm;
1042   PetscFunctionReturn(0);
1043 }
1044 
1045 static PetscErrorCode DMDestroy_Shell(DM dm)
1046 {
1047   PetscErrorCode ierr;
1048   DM_Shell       *shell = (DM_Shell*)dm->data;
1049 
1050   PetscFunctionBegin;
1051   ierr = MatDestroy(&shell->A);CHKERRQ(ierr);
1052   ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr);
1053   ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr);
1054   ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr);
1055   ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr);
1056   ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr);
1057   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
1058   ierr = PetscFree(shell);CHKERRQ(ierr);
1059   PetscFunctionReturn(0);
1060 }
1061 
1062 static PetscErrorCode DMView_Shell(DM dm,PetscViewer v)
1063 {
1064   PetscErrorCode ierr;
1065   DM_Shell       *shell = (DM_Shell*)dm->data;
1066 
1067   PetscFunctionBegin;
1068   ierr = VecView(shell->Xglobal,v);CHKERRQ(ierr);
1069   PetscFunctionReturn(0);
1070 }
1071 
1072 static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v)
1073 {
1074   PetscErrorCode ierr;
1075   DM_Shell       *shell = (DM_Shell*)dm->data;
1076 
1077   PetscFunctionBegin;
1078   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);CHKERRQ(ierr);
1079   ierr = VecLoad(shell->Xglobal,v);CHKERRQ(ierr);
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1084 {
1085   PetscErrorCode ierr;
1086 
1087   PetscFunctionBegin;
1088   if (subdm) {ierr = DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);CHKERRQ(ierr);}
1089   ierr = DMCreateSectionSubDM(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1090   PetscFunctionReturn(0);
1091 }
1092 
1093 PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
1094 {
1095   PetscErrorCode ierr;
1096   DM_Shell       *shell;
1097 
1098   PetscFunctionBegin;
1099   ierr     = PetscNewLog(dm,&shell);CHKERRQ(ierr);
1100   dm->data = shell;
1101 
1102   dm->ops->destroy            = DMDestroy_Shell;
1103   dm->ops->createglobalvector = DMCreateGlobalVector_Shell;
1104   dm->ops->createlocalvector  = DMCreateLocalVector_Shell;
1105   dm->ops->creatematrix       = DMCreateMatrix_Shell;
1106   dm->ops->view               = DMView_Shell;
1107   dm->ops->load               = DMLoad_Shell;
1108   dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell;
1109   dm->ops->globaltolocalend   = DMGlobalToLocalEndDefaultShell;
1110   dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell;
1111   dm->ops->localtoglobalend   = DMLocalToGlobalEndDefaultShell;
1112   dm->ops->localtolocalbegin  = DMLocalToLocalBeginDefaultShell;
1113   dm->ops->localtolocalend    = DMLocalToLocalEndDefaultShell;
1114   dm->ops->createsubdm        = DMCreateSubDM_Shell;
1115   ierr = DMSetMatType(dm,MATDENSE);CHKERRQ(ierr);
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 /*@
1120     DMShellCreate - Creates a shell DM object, used to manage user-defined problem data
1121 
1122     Collective
1123 
1124     Input Parameter:
1125 .   comm - the processors that will share the global vector
1126 
1127     Output Parameters:
1128 .   shell - the shell DM
1129 
1130     Level: advanced
1131 
1132 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext()
1133 @*/
1134 PetscErrorCode  DMShellCreate(MPI_Comm comm,DM *dm)
1135 {
1136   PetscErrorCode ierr;
1137 
1138   PetscFunctionBegin;
1139   PetscValidPointer(dm,2);
1140   ierr = DMCreate(comm,dm);CHKERRQ(ierr);
1141   ierr = DMSetType(*dm,DMSHELL);CHKERRQ(ierr);
1142   ierr = DMSetUp(*dm);CHKERRQ(ierr);
1143   PetscFunctionReturn(0);
1144 }
1145