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