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