Actual source code: dmshell.c
petsc-3.13.6 2020-09-29
1: #include <petscdmshell.h>
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: void *ctx;
13: } DM_Shell;
15: /*@
16: DMGlobalToLocalBeginDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to begin a global to local scatter
17: Collective
19: Input Arguments:
20: + dm - shell DM
21: . g - global vector
22: . mode - InsertMode
23: - l - local vector
25: Level: advanced
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.
29: .seealso: DMGlobalToLocalEndDefaultShell()
30: @*/
31: PetscErrorCode DMGlobalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
32: {
34: DM_Shell *shell = (DM_Shell*)dm->data;
37: if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
38: VecScatterBegin(shell->gtol,g,l,mode,SCATTER_FORWARD);
39: return(0);
40: }
42: /*@
43: DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to end a global to local scatter
44: Collective
46: Input Arguments:
47: + dm - shell DM
48: . g - global vector
49: . mode - InsertMode
50: - l - local vector
52: Level: advanced
54: .seealso: DMGlobalToLocalBeginDefaultShell()
55: @*/
56: PetscErrorCode DMGlobalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
57: {
59: DM_Shell *shell = (DM_Shell*)dm->data;
62: if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
63: VecScatterEnd(shell->gtol,g,l,mode,SCATTER_FORWARD);
64: return(0);
65: }
67: /*@
68: DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to begin a local to global scatter
69: Collective
71: Input Arguments:
72: + dm - shell DM
73: . l - local vector
74: . mode - InsertMode
75: - g - global vector
77: Level: advanced
79: 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.
81: .seealso: DMLocalToGlobalEndDefaultShell()
82: @*/
83: PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm,Vec l,InsertMode mode,Vec g)
84: {
86: DM_Shell *shell = (DM_Shell*)dm->data;
89: if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
90: VecScatterBegin(shell->ltog,l,g,mode,SCATTER_FORWARD);
91: return(0);
92: }
94: /*@
95: DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to end a local to global scatter
96: Collective
98: Input Arguments:
99: + dm - shell DM
100: . l - local vector
101: . mode - InsertMode
102: - g - global vector
104: Level: advanced
106: .seealso: DMLocalToGlobalBeginDefaultShell()
107: @*/
108: PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm,Vec l,InsertMode mode,Vec g)
109: {
111: DM_Shell *shell = (DM_Shell*)dm->data;
114: if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
115: VecScatterEnd(shell->ltog,l,g,mode,SCATTER_FORWARD);
116: return(0);
117: }
119: /*@
120: DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal VecScatter context set by the user to begin a local to local scatter
121: Collective
123: Input Arguments:
124: + dm - shell DM
125: . g - the original local vector
126: - mode - InsertMode
128: Output Parameter:
129: . l - the local vector with correct ghost values
131: Level: advanced
133: 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.
135: .seealso: DMLocalToLocalEndDefaultShell()
136: @*/
137: PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
138: {
140: DM_Shell *shell = (DM_Shell*)dm->data;
143: if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()");
144: VecScatterBegin(shell->ltol,g,l,mode,SCATTER_FORWARD);
145: return(0);
146: }
148: /*@
149: DMLocalToLocalEndDefaultShell - Uses the LocalToLocal VecScatter context set by the user to end a local to local scatter
150: Collective
152: Input Arguments:
153: + dm - shell DM
154: . g - the original local vector
155: - mode - InsertMode
157: Output Parameter:
158: . l - the local vector with correct ghost values
160: Level: advanced
162: .seealso: DMLocalToLocalBeginDefaultShell()
163: @*/
164: PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l)
165: {
167: DM_Shell *shell = (DM_Shell*)dm->data;
170: if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
171: VecScatterEnd(shell->ltol,g,l,mode,SCATTER_FORWARD);
172: return(0);
173: }
175: static PetscErrorCode DMCreateMatrix_Shell(DM dm,Mat *J)
176: {
178: DM_Shell *shell = (DM_Shell*)dm->data;
179: Mat A;
184: if (!shell->A) {
185: if (shell->Xglobal) {
186: PetscInt m,M;
187: PetscInfo(dm,"Naively creating matrix using global vector distribution without preallocation\n");
188: VecGetSize(shell->Xglobal,&M);
189: VecGetLocalSize(shell->Xglobal,&m);
190: MatCreate(PetscObjectComm((PetscObject)dm),&shell->A);
191: MatSetSizes(shell->A,m,m,M,M);
192: MatSetType(shell->A,dm->mattype);
193: MatSetUp(shell->A);
194: } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector");
195: }
196: A = shell->A;
197: /* the check below is tacky and incomplete */
198: if (dm->mattype) {
199: PetscBool flg,aij,seqaij,mpiaij;
200: PetscObjectTypeCompare((PetscObject)A,dm->mattype,&flg);
201: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij);
202: PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&mpiaij);
203: PetscStrcmp(dm->mattype,MATAIJ,&aij);
204: if (!flg) {
205: 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);
206: }
207: }
208: /* Need to create a copy in order to attach the DM to the matrix */
209: MatDuplicate(A,MAT_SHARE_NONZERO_PATTERN,J);
210: MatSetDM(*J,dm);
211: return(0);
212: }
214: PetscErrorCode DMCreateGlobalVector_Shell(DM dm,Vec *gvec)
215: {
217: DM_Shell *shell = (DM_Shell*)dm->data;
218: Vec X;
223: *gvec = 0;
224: X = shell->Xglobal;
225: if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()");
226: /* Need to create a copy in order to attach the DM to the vector */
227: VecDuplicate(X,gvec);
228: VecZeroEntries(*gvec);
229: VecSetDM(*gvec,dm);
230: return(0);
231: }
233: PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec)
234: {
236: DM_Shell *shell = (DM_Shell*)dm->data;
237: Vec X;
242: *gvec = 0;
243: X = shell->Xlocal;
244: if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()");
245: /* Need to create a copy in order to attach the DM to the vector */
246: VecDuplicate(X,gvec);
247: VecZeroEntries(*gvec);
248: VecSetDM(*gvec,dm);
249: return(0);
250: }
252: /*@
253: DMShellSetContext - set some data to be usable by this DM
255: Collective
257: Input Arguments:
258: + dm - shell DM
259: - ctx - the context
261: Level: advanced
263: .seealso: DMCreateMatrix(), DMShellGetContext()
264: @*/
265: PetscErrorCode DMShellSetContext(DM dm,void *ctx)
266: {
267: DM_Shell *shell = (DM_Shell*)dm->data;
269: PetscBool isshell;
273: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
274: if (!isshell) return(0);
275: shell->ctx = ctx;
276: return(0);
277: }
279: /*@
280: DMShellGetContext - Returns the user-provided context associated to the DM
282: Collective
284: Input Argument:
285: . dm - shell DM
287: Output Argument:
288: . ctx - the context
290: Level: advanced
292: .seealso: DMCreateMatrix(), DMShellSetContext()
293: @*/
294: PetscErrorCode DMShellGetContext(DM dm,void **ctx)
295: {
296: DM_Shell *shell = (DM_Shell*)dm->data;
298: PetscBool isshell;
302: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
303: if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
304: *ctx = shell->ctx;
305: return(0);
306: }
308: /*@
309: DMShellSetMatrix - sets a template matrix associated with the DMShell
311: Collective
313: Input Arguments:
314: + dm - shell DM
315: - J - template matrix
317: Level: advanced
319: Developer Notes:
320: 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.
322: .seealso: DMCreateMatrix(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
323: @*/
324: PetscErrorCode DMShellSetMatrix(DM dm,Mat J)
325: {
326: DM_Shell *shell = (DM_Shell*)dm->data;
328: PetscBool isshell;
329: DM mdm;
334: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
335: if (!isshell) return(0);
336: if (J == shell->A) return(0);
337: MatGetDM(J,&mdm);
338: PetscObjectReference((PetscObject)J);
339: MatDestroy(&shell->A);
340: if (mdm == dm) {
341: MatDuplicate(J,MAT_SHARE_NONZERO_PATTERN,&shell->A);
342: MatSetDM(shell->A,NULL);
343: } else shell->A = J;
344: return(0);
345: }
347: /*@C
348: DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM
350: Logically Collective on dm
352: Input Arguments:
353: + dm - the shell DM
354: - func - the function to create a matrix
356: Level: advanced
358: .seealso: DMCreateMatrix(), DMShellSetMatrix(), DMShellSetContext(), DMShellGetContext()
359: @*/
360: PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*))
361: {
364: dm->ops->creatematrix = func;
365: return(0);
366: }
368: /*@
369: DMShellSetGlobalVector - sets a template global vector associated with the DMShell
371: Logically Collective on dm
373: Input Arguments:
374: + dm - shell DM
375: - X - template vector
377: Level: advanced
379: .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector()
380: @*/
381: PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X)
382: {
383: DM_Shell *shell = (DM_Shell*)dm->data;
385: PetscBool isshell;
386: DM vdm;
391: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
392: if (!isshell) return(0);
393: VecGetDM(X,&vdm);
394: /*
395: if the vector proposed as the new base global vector for the DM is a DM vector associated
396: with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
397: we get a circular dependency that prevents the DM from being destroy when it should be.
398: This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
399: DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
400: to set its input vector (which is associated with the DM) as the base global vector.
401: Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
402: for pointing out the problem.
403: */
404: if (vdm == dm) return(0);
405: PetscObjectReference((PetscObject)X);
406: VecDestroy(&shell->Xglobal);
407: shell->Xglobal = X;
408: return(0);
409: }
411: /*@C
412: DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM
414: Logically Collective
416: Input Arguments:
417: + dm - the shell DM
418: - func - the creation routine
420: Level: advanced
422: .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
423: @*/
424: PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
425: {
428: dm->ops->createglobalvector = func;
429: return(0);
430: }
432: /*@
433: DMShellSetLocalVector - sets a template local vector associated with the DMShell
435: Logically Collective on dm
437: Input Arguments:
438: + dm - shell DM
439: - X - template vector
441: Level: advanced
443: .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector()
444: @*/
445: PetscErrorCode DMShellSetLocalVector(DM dm,Vec X)
446: {
447: DM_Shell *shell = (DM_Shell*)dm->data;
449: PetscBool isshell;
450: DM vdm;
455: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
456: if (!isshell) return(0);
457: VecGetDM(X,&vdm);
458: /*
459: if the vector proposed as the new base global vector for the DM is a DM vector associated
460: with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
461: we get a circular dependency that prevents the DM from being destroy when it should be.
462: This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
463: DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
464: to set its input vector (which is associated with the DM) as the base global vector.
465: Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
466: for pointing out the problem.
467: */
468: if (vdm == dm) return(0);
469: PetscObjectReference((PetscObject)X);
470: VecDestroy(&shell->Xlocal);
471: shell->Xlocal = X;
472: return(0);
473: }
475: /*@C
476: DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM
478: Logically Collective
480: Input Arguments:
481: + dm - the shell DM
482: - func - the creation routine
484: Level: advanced
486: .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
487: @*/
488: PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
489: {
492: dm->ops->createlocalvector = func;
493: return(0);
494: }
496: /*@C
497: DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter
499: Logically Collective on dm
501: Input Arguments
502: + dm - the shell DM
503: . begin - the routine that begins the global to local scatter
504: - end - the routine that ends the global to local scatter
506: Notes:
507: If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then
508: DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers
510: Level: advanced
512: .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
513: @*/
514: PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
517: dm->ops->globaltolocalbegin = begin;
518: dm->ops->globaltolocalend = end;
519: return(0);
520: }
522: /*@C
523: DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter
525: Logically Collective on dm
527: Input Arguments
528: + dm - the shell DM
529: . begin - the routine that begins the local to global scatter
530: - end - the routine that ends the local to global scatter
532: Notes:
533: If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then
534: DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers
536: Level: advanced
538: .seealso: DMShellSetGlobalToLocal()
539: @*/
540: PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
543: dm->ops->localtoglobalbegin = begin;
544: dm->ops->localtoglobalend = end;
545: return(0);
546: }
548: /*@C
549: DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter
551: Logically Collective on dm
553: Input Arguments
554: + dm - the shell DM
555: . begin - the routine that begins the local to local scatter
556: - end - the routine that ends the local to local scatter
558: Notes:
559: If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then
560: DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers
562: Level: advanced
564: .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
565: @*/
566: PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
569: dm->ops->localtolocalbegin = begin;
570: dm->ops->localtolocalend = end;
571: return(0);
572: }
574: /*@
575: DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication
577: Logically Collective on dm
579: Input Arguments
580: + dm - the shell DM
581: - gtol - the global to local VecScatter context
583: Level: advanced
585: .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
586: @*/
587: PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
588: {
589: DM_Shell *shell = (DM_Shell*)dm->data;
595: PetscObjectReference((PetscObject)gtol);
596: VecScatterDestroy(&shell->gtol);
597: shell->gtol = gtol;
598: return(0);
599: }
601: /*@
602: DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication
604: Logically Collective on dm
606: Input Arguments
607: + dm - the shell DM
608: - ltog - the local to global VecScatter context
610: Level: advanced
612: .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell()
613: @*/
614: PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
615: {
616: DM_Shell *shell = (DM_Shell*)dm->data;
622: PetscObjectReference((PetscObject)ltog);
623: VecScatterDestroy(&shell->ltog);
624: shell->ltog = ltog;
625: return(0);
626: }
628: /*@
629: DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication
631: Logically Collective on dm
633: Input Arguments
634: + dm - the shell DM
635: - ltol - the local to local VecScatter context
637: Level: advanced
639: .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
640: @*/
641: PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
642: {
643: DM_Shell *shell = (DM_Shell*)dm->data;
649: PetscObjectReference((PetscObject)ltol);
650: VecScatterDestroy(&shell->ltol);
651: shell->ltol = ltol;
652: return(0);
653: }
655: /*@C
656: DMShellSetCoarsen - Set the routine used to coarsen the shell DM
658: Logically Collective on dm
660: Input Arguments
661: + dm - the shell DM
662: - coarsen - the routine that coarsens the DM
664: Level: advanced
666: .seealso: DMShellSetRefine(), DMCoarsen(), DMShellGetCoarsen(), DMShellSetContext(), DMShellGetContext()
667: @*/
668: PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*))
669: {
671: PetscBool isshell;
675: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
676: if (!isshell) return(0);
677: dm->ops->coarsen = coarsen;
678: return(0);
679: }
681: /*@C
682: DMShellGetCoarsen - Get the routine used to coarsen the shell DM
684: Logically Collective on dm
686: Input Argument:
687: . dm - the shell DM
689: Output Argument:
690: . coarsen - the routine that coarsens the DM
692: Level: advanced
694: .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
695: @*/
696: PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM,MPI_Comm,DM*))
697: {
699: PetscBool isshell;
703: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
704: if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
705: *coarsen = dm->ops->coarsen;
706: return(0);
707: }
709: /*@C
710: DMShellSetRefine - Set the routine used to refine the shell DM
712: Logically Collective on dm
714: Input Arguments
715: + dm - the shell DM
716: - refine - the routine that refines the DM
718: Level: advanced
720: .seealso: DMShellSetCoarsen(), DMRefine(), DMShellGetRefine(), DMShellSetContext(), DMShellGetContext()
721: @*/
722: PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*))
723: {
725: PetscBool isshell;
729: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
730: if (!isshell) return(0);
731: dm->ops->refine = refine;
732: return(0);
733: }
735: /*@C
736: DMShellGetRefine - Get the routine used to refine the shell DM
738: Logically Collective on dm
740: Input Argument:
741: . dm - the shell DM
743: Output Argument:
744: . refine - the routine that refines the DM
746: Level: advanced
748: .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
749: @*/
750: PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM,MPI_Comm,DM*))
751: {
753: PetscBool isshell;
757: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
758: if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
759: *refine = dm->ops->refine;
760: return(0);
761: }
763: /*@C
764: DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator
766: Logically Collective on dm
768: Input Arguments
769: + dm - the shell DM
770: - interp - the routine to create the interpolation
772: Level: advanced
774: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
775: @*/
776: PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*))
777: {
779: PetscBool isshell;
783: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
784: if (!isshell) return(0);
785: dm->ops->createinterpolation = interp;
786: return(0);
787: }
789: /*@C
790: DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator
792: Logically Collective on dm
794: Input Argument:
795: + dm - the shell DM
797: Output Argument:
798: - interp - the routine to create the interpolation
800: Level: advanced
802: .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
803: @*/
804: PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM,DM,Mat*,Vec*))
805: {
807: PetscBool isshell;
811: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
812: if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
813: *interp = dm->ops->createinterpolation;
814: return(0);
815: }
817: /*@C
818: DMShellSetCreateRestriction - Set the routine used to create the restriction operator
820: Logically Collective on dm
822: Input Arguments
823: + dm - the shell DM
824: - striction- the routine to create the restriction
826: Level: advanced
828: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
829: @*/
830: PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*))
831: {
833: PetscBool isshell;
837: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
838: if (!isshell) return(0);
839: dm->ops->createrestriction = restriction;
840: return(0);
841: }
843: /*@C
844: DMShellGetCreateRestriction - Get the routine used to create the restriction operator
846: Logically Collective on dm
848: Input Argument:
849: + dm - the shell DM
851: Output Argument:
852: - restriction - the routine to create the restriction
854: Level: advanced
856: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext()
857: @*/
858: PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM,DM,Mat*))
859: {
861: PetscBool isshell;
865: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
866: if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
867: *restriction = dm->ops->createrestriction;
868: return(0);
869: }
871: /*@C
872: DMShellSetCreateInjection - Set the routine used to create the injection operator
874: Logically Collective on dm
876: Input Arguments
877: + dm - the shell DM
878: - inject - the routine to create the injection
880: Level: advanced
882: .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellGetCreateInjection(), DMShellSetContext(), DMShellGetContext()
883: @*/
884: PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*))
885: {
887: PetscBool isshell;
891: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
892: if (!isshell) return(0);
893: dm->ops->createinjection = inject;
894: return(0);
895: }
897: /*@C
898: DMShellGetCreateInjection - Get the routine used to create the injection operator
900: Logically Collective on dm
902: Input Argument:
903: + dm - the shell DM
905: Output Argument:
906: - inject - the routine to create the injection
908: Level: advanced
910: .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext()
911: @*/
912: PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM,DM,Mat*))
913: {
915: PetscBool isshell;
919: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
920: if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
921: *inject = dm->ops->createinjection;
922: return(0);
923: }
925: /*@C
926: DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM
928: Logically Collective on dm
930: Input Arguments
931: + dm - the shell DM
932: - decomp - the routine to create the decomposition
934: Level: advanced
936: .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext()
937: @*/
938: PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**))
939: {
941: PetscBool isshell;
945: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
946: if (!isshell) return(0);
947: dm->ops->createfielddecomposition = decomp;
948: return(0);
949: }
951: /*@C
952: DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM
954: Logically Collective on dm
956: Input Arguments
957: + dm - the shell DM
958: - decomp - the routine to create the decomposition
960: Level: advanced
962: .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext()
963: @*/
964: PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**))
965: {
967: PetscBool isshell;
971: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
972: if (!isshell) return(0);
973: dm->ops->createdomaindecomposition = decomp;
974: return(0);
975: }
977: /*@C
978: DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM
980: Logically Collective on dm
982: Input Arguments
983: + dm - the shell DM
984: - scatter - the routine to create the scatters
986: Level: advanced
988: .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext()
989: @*/
990: PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**))
991: {
993: PetscBool isshell;
997: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
998: if (!isshell) return(0);
999: dm->ops->createddscatters = scatter;
1000: return(0);
1001: }
1003: /*@C
1004: DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM
1006: Logically Collective on dm
1008: Input Arguments
1009: + dm - the shell DM
1010: - subdm - the routine to create the decomposition
1012: Level: advanced
1014: .seealso: DMCreateSubDM(), DMShellGetCreateSubDM(), DMShellSetContext(), DMShellGetContext()
1015: @*/
1016: PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
1017: {
1019: PetscBool isshell;
1023: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
1024: if (!isshell) return(0);
1025: dm->ops->createsubdm = subdm;
1026: return(0);
1027: }
1029: /*@C
1030: DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM
1032: Logically Collective on dm
1034: Input Argument:
1035: + dm - the shell DM
1037: Output Argument:
1038: - subdm - the routine to create the decomposition
1040: Level: advanced
1042: .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext()
1043: @*/
1044: PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
1045: {
1047: PetscBool isshell;
1051: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
1052: if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
1053: *subdm = dm->ops->createsubdm;
1054: return(0);
1055: }
1057: static PetscErrorCode DMDestroy_Shell(DM dm)
1058: {
1060: DM_Shell *shell = (DM_Shell*)dm->data;
1063: MatDestroy(&shell->A);
1064: VecDestroy(&shell->Xglobal);
1065: VecDestroy(&shell->Xlocal);
1066: VecScatterDestroy(&shell->gtol);
1067: VecScatterDestroy(&shell->ltog);
1068: VecScatterDestroy(&shell->ltol);
1069: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
1070: PetscFree(shell);
1071: return(0);
1072: }
1074: static PetscErrorCode DMView_Shell(DM dm,PetscViewer v)
1075: {
1077: DM_Shell *shell = (DM_Shell*)dm->data;
1080: VecView(shell->Xglobal,v);
1081: return(0);
1082: }
1084: static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v)
1085: {
1087: DM_Shell *shell = (DM_Shell*)dm->data;
1090: VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);
1091: VecLoad(shell->Xglobal,v);
1092: return(0);
1093: }
1095: PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1096: {
1100: if (subdm) {DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);}
1101: DMCreateSectionSubDM(dm, numFields, fields, is, subdm);
1102: return(0);
1103: }
1105: PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
1106: {
1108: DM_Shell *shell;
1111: PetscNewLog(dm,&shell);
1112: dm->data = shell;
1114: dm->ops->destroy = DMDestroy_Shell;
1115: dm->ops->createglobalvector = DMCreateGlobalVector_Shell;
1116: dm->ops->createlocalvector = DMCreateLocalVector_Shell;
1117: dm->ops->creatematrix = DMCreateMatrix_Shell;
1118: dm->ops->view = DMView_Shell;
1119: dm->ops->load = DMLoad_Shell;
1120: dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell;
1121: dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell;
1122: dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell;
1123: dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell;
1124: dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell;
1125: dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell;
1126: dm->ops->createsubdm = DMCreateSubDM_Shell;
1127: return(0);
1128: }
1130: /*@
1131: DMShellCreate - Creates a shell DM object, used to manage user-defined problem data
1133: Collective
1135: Input Parameter:
1136: . comm - the processors that will share the global vector
1138: Output Parameters:
1139: . shell - the shell DM
1141: Level: advanced
1143: .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext()
1144: @*/
1145: PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm)
1146: {
1151: DMCreate(comm,dm);
1152: DMSetType(*dm,DMSHELL);
1153: DMSetUp(*dm);
1154: return(0);
1155: }