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