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