Actual source code: dmshell.c
petsc-3.6.4 2016-04-12
1: #include <petscdmshell.h> /*I "petscdmshell.h" I*/
2: #include <petscmat.h>
3: #include <petsc/private/dmimpl.h>
5: typedef struct {
6: Vec Xglobal;
7: Vec Xlocal;
8: Mat A;
9: VecScatter gtol;
10: VecScatter ltog;
11: VecScatter ltol;
12: } DM_Shell;
16: /*@
17: DMGlobalToLocalBeginDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to begin a global to local scatter
18: Collective
20: Input Arguments:
21: + dm - shell DM
22: . g - global vector
23: . mode - InsertMode
24: - l - local vector
26: Level: advanced
28: 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: .seealso: DMGlobalToLocalEndDefaultShell()
31: @*/
32: PetscErrorCode DMGlobalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
33: {
35: DM_Shell *shell = (DM_Shell*)dm->data;
38: if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
39: VecScatterBegin(shell->gtol,g,l,mode,SCATTER_FORWARD);
40: return(0);
41: }
45: /*@
46: DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to end a global to local scatter
47: Collective
49: Input Arguments:
50: + dm - shell DM
51: . g - global vector
52: . mode - InsertMode
53: - l - local vector
55: Level: advanced
57: .seealso: DMGlobalToLocalBeginDefaultShell()
58: @*/
59: PetscErrorCode DMGlobalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
60: {
62: DM_Shell *shell = (DM_Shell*)dm->data;
65: if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
66: VecScatterEnd(shell->gtol,g,l,mode,SCATTER_FORWARD);
67: return(0);
68: }
72: /*@
73: DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to begin a local to global scatter
74: Collective
76: Input Arguments:
77: + dm - shell DM
78: . l - local vector
79: . mode - InsertMode
80: - g - global vector
82: Level: advanced
84: 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: .seealso: DMLocalToGlobalEndDefaultShell()
87: @*/
88: PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm,Vec l,InsertMode mode,Vec g)
89: {
91: DM_Shell *shell = (DM_Shell*)dm->data;
94: if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
95: VecScatterBegin(shell->ltog,l,g,mode,SCATTER_FORWARD);
96: return(0);
97: }
101: /*@
102: DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to end a local to global scatter
103: Collective
105: Input Arguments:
106: + dm - shell DM
107: . l - local vector
108: . mode - InsertMode
109: - g - global vector
111: Level: advanced
113: .seealso: DMLocalToGlobalBeginDefaultShell()
114: @*/
115: PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm,Vec l,InsertMode mode,Vec g)
116: {
118: DM_Shell *shell = (DM_Shell*)dm->data;
121: if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
122: VecScatterEnd(shell->ltog,l,g,mode,SCATTER_FORWARD);
123: return(0);
124: }
128: /*@
129: DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal VecScatter context set by the user to begin a local to local scatter
130: Collective
132: Input Arguments:
133: + dm - shell DM
134: . g - the original local vector
135: - mode - InsertMode
137: Output Parameter:
138: . l - the local vector with correct ghost values
140: Level: advanced
142: 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: .seealso: DMLocalToLocalEndDefaultShell()
145: @*/
146: PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
147: {
149: DM_Shell *shell = (DM_Shell*)dm->data;
152: if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()");
153: VecScatterBegin(shell->ltol,g,l,mode,SCATTER_FORWARD);
154: return(0);
155: }
159: /*@
160: DMLocalToLocalEndDefaultShell - Uses the LocalToLocal VecScatter context set by the user to end a local to local scatter
161: Collective
163: Input Arguments:
164: + dm - shell DM
165: . g - the original local vector
166: - mode - InsertMode
168: Output Parameter:
169: . l - the local vector with correct ghost values
171: Level: advanced
173: .seealso: DMLocalToLocalBeginDefaultShell()
174: @*/
175: PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
176: {
178: DM_Shell *shell = (DM_Shell*)dm->data;
181: if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
182: VecScatterEnd(shell->ltol,g,l,mode,SCATTER_FORWARD);
183: return(0);
184: }
188: static PetscErrorCode DMCreateMatrix_Shell(DM dm,Mat *J)
189: {
191: DM_Shell *shell = (DM_Shell*)dm->data;
192: Mat A;
197: if (!shell->A) {
198: if (shell->Xglobal) {
199: PetscInt m,M;
200: PetscInfo(dm,"Naively creating matrix using global vector distribution without preallocation\n");
201: VecGetSize(shell->Xglobal,&M);
202: VecGetLocalSize(shell->Xglobal,&m);
203: MatCreate(PetscObjectComm((PetscObject)dm),&shell->A);
204: MatSetSizes(shell->A,m,m,M,M);
205: MatSetType(shell->A,dm->mattype);
206: MatSetUp(shell->A);
207: } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector");
208: }
209: A = shell->A;
210: /* the check below is tacky and incomplete */
211: if (dm->mattype) {
212: PetscBool flg,aij,seqaij,mpiaij;
213: PetscObjectTypeCompare((PetscObject)A,dm->mattype,&flg);
214: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij);
215: PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&mpiaij);
216: PetscStrcmp(dm->mattype,MATAIJ,&aij);
217: if (!flg) {
218: 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);
219: }
220: }
221: if (((PetscObject)A)->refct < 2) { /* We have an exclusive reference so we can give it out */
222: PetscObjectReference((PetscObject)A);
223: MatZeroEntries(A);
224: *J = A;
225: } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */
226: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,J);
227: MatZeroEntries(*J);
228: }
229: return(0);
230: }
234: PetscErrorCode DMCreateGlobalVector_Shell(DM dm,Vec *gvec)
235: {
237: DM_Shell *shell = (DM_Shell*)dm->data;
238: Vec X;
243: *gvec = 0;
244: X = shell->Xglobal;
245: if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()");
246: if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */
247: PetscObjectReference((PetscObject)X);
248: VecZeroEntries(X);
249: *gvec = X;
250: } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */
251: VecDuplicate(X,gvec);
252: VecZeroEntries(*gvec);
253: }
254: VecSetDM(*gvec,dm);
255: return(0);
256: }
260: PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec)
261: {
263: DM_Shell *shell = (DM_Shell*)dm->data;
264: Vec X;
269: *gvec = 0;
270: X = shell->Xlocal;
271: if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()");
272: if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */
273: PetscObjectReference((PetscObject)X);
274: VecZeroEntries(X);
275: *gvec = X;
276: } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */
277: VecDuplicate(X,gvec);
278: VecZeroEntries(*gvec);
279: }
280: VecSetDM(*gvec,dm);
281: return(0);
282: }
286: /*@
287: DMShellSetMatrix - sets a template matrix associated with the DMShell
289: Collective
291: Input Arguments:
292: + dm - shell DM
293: - J - template matrix
295: Level: advanced
297: .seealso: DMCreateMatrix(), DMShellSetCreateMatrix()
298: @*/
299: PetscErrorCode DMShellSetMatrix(DM dm,Mat J)
300: {
301: DM_Shell *shell = (DM_Shell*)dm->data;
303: PetscBool isshell;
308: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
309: if (!isshell) return(0);
310: PetscObjectReference((PetscObject)J);
311: MatDestroy(&shell->A);
312: shell->A = J;
313: return(0);
314: }
318: /*@C
319: DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM
321: Logically Collective on DM
323: Input Arguments:
324: + dm - the shell DM
325: - func - the function to create a matrix
327: Level: advanced
329: .seealso: DMCreateMatrix(), DMShellSetMatrix()
330: @*/
331: PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*))
332: {
336: dm->ops->creatematrix = func;
337: return(0);
338: }
342: /*@
343: DMShellSetGlobalVector - sets a template global vector associated with the DMShell
345: Logically Collective on DM
347: Input Arguments:
348: + dm - shell DM
349: - X - template vector
351: Level: advanced
353: .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector()
354: @*/
355: PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X)
356: {
357: DM_Shell *shell = (DM_Shell*)dm->data;
359: PetscBool isshell;
364: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
365: if (!isshell) return(0);
366: PetscObjectReference((PetscObject)X);
367: VecDestroy(&shell->Xglobal);
368: shell->Xglobal = X;
369: return(0);
370: }
374: /*@C
375: DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM
377: Logically Collective
379: Input Arguments:
380: + dm - the shell DM
381: - func - the creation routine
383: Level: advanced
385: .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix()
386: @*/
387: PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
388: {
392: dm->ops->createglobalvector = func;
393: return(0);
394: }
398: /*@
399: DMShellSetLocalVector - sets a template local vector associated with the DMShell
401: Logically Collective on DM
403: Input Arguments:
404: + dm - shell DM
405: - X - template vector
407: Level: advanced
409: .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector()
410: @*/
411: PetscErrorCode DMShellSetLocalVector(DM dm,Vec X)
412: {
413: DM_Shell *shell = (DM_Shell*)dm->data;
415: PetscBool isshell;
420: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
421: if (!isshell) return(0);
422: PetscObjectReference((PetscObject)X);
423: VecDestroy(&shell->Xlocal);
424: shell->Xlocal = X;
425: return(0);
426: }
430: /*@C
431: DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM
433: Logically Collective
435: Input Arguments:
436: + dm - the shell DM
437: - func - the creation routine
439: Level: advanced
441: .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix()
442: @*/
443: PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
444: {
448: dm->ops->createlocalvector = func;
449: return(0);
450: }
454: /*@C
455: DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter
457: Logically Collective on DM
459: Input Arguments
460: + dm - the shell DM
461: . begin - the routine that begins the global to local scatter
462: - end - the routine that ends the global to local scatter
464: Notes: If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then
465: DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers
467: Level: advanced
469: .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
470: @*/
471: PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
473: dm->ops->globaltolocalbegin = begin;
474: dm->ops->globaltolocalend = end;
475: return(0);
476: }
480: /*@C
481: DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter
483: Logically Collective on DM
485: Input Arguments
486: + dm - the shell DM
487: . begin - the routine that begins the local to global scatter
488: - end - the routine that ends the local to global scatter
490: Notes: If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then
491: DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers
493: Level: advanced
495: .seealso: DMShellSetGlobalToLocal()
496: @*/
497: PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
499: dm->ops->localtoglobalbegin = begin;
500: dm->ops->localtoglobalend = end;
501: return(0);
502: }
506: /*@C
507: DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter
509: Logically Collective on DM
511: Input Arguments
512: + dm - the shell DM
513: . begin - the routine that begins the local to local scatter
514: - end - the routine that ends the local to local scatter
516: Notes: If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then
517: DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers
519: Level: advanced
521: .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
522: @*/
523: PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
525: dm->ops->localtolocalbegin = begin;
526: dm->ops->localtolocalend = end;
527: return(0);
528: }
532: /*@
533: DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication
535: Logically Collective on DM
537: Input Arguments
538: + dm - the shell DM
539: - gtol - the global to local VecScatter context
541: Level: advanced
543: .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
544: @*/
545: PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
546: {
547: DM_Shell *shell = (DM_Shell*)dm->data;
551: PetscObjectReference((PetscObject)gtol);
552: /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */
553: VecScatterDestroy(&shell->gtol);
554: shell->gtol = gtol;
555: return(0);
556: }
560: /*@
561: DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication
563: Logically Collective on DM
565: Input Arguments
566: + dm - the shell DM
567: - ltog - the local to global VecScatter context
569: Level: advanced
571: .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell()
572: @*/
573: PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
574: {
575: DM_Shell *shell = (DM_Shell*)dm->data;
579: PetscObjectReference((PetscObject)ltog);
580: /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */
581: VecScatterDestroy(&shell->ltog);
582: shell->ltog = ltog;
583: return(0);
584: }
588: /*@
589: DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication
591: Logically Collective on DM
593: Input Arguments
594: + dm - the shell DM
595: - ltol - the local to local VecScatter context
597: Level: advanced
599: .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
600: @*/
601: PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
602: {
603: DM_Shell *shell = (DM_Shell*)dm->data;
607: PetscObjectReference((PetscObject)ltol);
608: /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */
609: VecScatterDestroy(&shell->ltol);
610: shell->ltol = ltol;
611: return(0);
612: }
616: /*@C
617: DMShellSetCoarsen - Set the routine used to coarsen the shell DM
619: Logically Collective on DM
621: Input Arguments
622: + dm - the shell DM
623: - coarsen - the routine that coarsens the DM
625: Level: advanced
627: .seealso: DMShellSetRefine(), DMCoarsen()
628: @*/
629: PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*))
630: {
632: PetscBool isshell;
636: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
637: if (!isshell) return(0);
638: dm->ops->coarsen = coarsen;
639: return(0);
640: }
644: /*@C
645: DMShellSetRefine - Set the routine used to refine the shell DM
647: Logically Collective on DM
649: Input Arguments
650: + dm - the shell DM
651: - refine - the routine that refines the DM
653: Level: advanced
655: .seealso: DMShellSetCoarsen(), DMRefine()
656: @*/
657: PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*))
658: {
660: PetscBool isshell;
664: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
665: if (!isshell) return(0);
666: dm->ops->refine = refine;
667: return(0);
668: }
672: /*@C
673: DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator
675: Logically Collective on DM
677: Input Arguments
678: + dm - the shell DM
679: - interp - the routine to create the interpolation
681: Level: advanced
683: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation()
684: @*/
685: PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*))
686: {
688: PetscBool isshell;
692: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
693: if (!isshell) return(0);
694: dm->ops->createinterpolation = interp;
695: return(0);
696: }
700: /*@C
701: DMShellSetCreateInjection - Set the routine used to create the injection operator
703: Logically Collective on DM
705: Input Arguments
706: + dm - the shell DM
707: - inject - the routine to create the injection
709: Level: advanced
711: .seealso: DMShellSetCreateInterpolation(), DMCreateInjection()
712: @*/
713: PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*))
714: {
716: PetscBool isshell;
720: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
721: if (!isshell) return(0);
722: dm->ops->getinjection = inject;
723: return(0);
724: }
728: /*@C
729: DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM
731: Logically Collective on DM
733: Input Arguments
734: + dm - the shell DM
735: - decomp - the routine to create the decomposition
737: Level: advanced
739: .seealso: DMCreateFieldDecomposition()
740: @*/
741: PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**))
742: {
744: PetscBool isshell;
748: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
749: if (!isshell) return(0);
750: dm->ops->createfielddecomposition = decomp;
751: return(0);
752: }
756: static PetscErrorCode DMDestroy_Shell(DM dm)
757: {
759: DM_Shell *shell = (DM_Shell*)dm->data;
762: MatDestroy(&shell->A);
763: VecDestroy(&shell->Xglobal);
764: VecDestroy(&shell->Xlocal);
765: VecScatterDestroy(&shell->gtol);
766: VecScatterDestroy(&shell->ltog);
767: VecScatterDestroy(&shell->ltol);
768: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
769: PetscFree(shell);
770: return(0);
771: }
775: static PetscErrorCode DMView_Shell(DM dm,PetscViewer v)
776: {
778: DM_Shell *shell = (DM_Shell*)dm->data;
781: VecView(shell->Xglobal,v);
782: return(0);
783: }
787: static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v)
788: {
790: DM_Shell *shell = (DM_Shell*)dm->data;
793: VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);
794: VecLoad(shell->Xglobal,v);
795: return(0);
796: }
800: PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
801: {
805: if (subdm) {DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);}
806: DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);
807: return(0);
808: }
812: PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
813: {
815: DM_Shell *shell;
818: PetscNewLog(dm,&shell);
819: dm->data = shell;
821: PetscObjectChangeTypeName((PetscObject)dm,DMSHELL);
823: dm->ops->destroy = DMDestroy_Shell;
824: dm->ops->createglobalvector = DMCreateGlobalVector_Shell;
825: dm->ops->createlocalvector = DMCreateLocalVector_Shell;
826: dm->ops->creatematrix = DMCreateMatrix_Shell;
827: dm->ops->view = DMView_Shell;
828: dm->ops->load = DMLoad_Shell;
829: dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell;
830: dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell;
831: dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell;
832: dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell;
833: dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell;
834: dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell;
835: dm->ops->createsubdm = DMCreateSubDM_Shell;
836: return(0);
837: }
841: /*@
842: DMShellCreate - Creates a shell DM object, used to manage user-defined problem data
844: Collective on MPI_Comm
846: Input Parameter:
847: . comm - the processors that will share the global vector
849: Output Parameters:
850: . shell - the shell DM
852: Level: advanced
854: .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector()
855: @*/
856: PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm)
857: {
862: DMCreate(comm,dm);
863: DMSetType(*dm,DMSHELL);
864: DMSetUp(*dm);
865: return(0);
866: }