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