Actual source code: dmshell.c
petsc-3.10.5 2019-03-28
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: if (((PetscObject)A)->refct < 2) { /* We have an exclusive reference so we can give it out */
209: PetscObjectReference((PetscObject)A);
210: MatZeroEntries(A);
211: *J = A;
212: } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */
213: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,J);
214: MatZeroEntries(*J);
215: }
216: return(0);
217: }
219: PetscErrorCode DMCreateGlobalVector_Shell(DM dm,Vec *gvec)
220: {
222: DM_Shell *shell = (DM_Shell*)dm->data;
223: Vec X;
228: *gvec = 0;
229: X = shell->Xglobal;
230: if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()");
231: if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */
232: PetscObjectReference((PetscObject)X);
233: VecZeroEntries(X);
234: *gvec = X;
235: } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */
236: VecDuplicate(X,gvec);
237: VecZeroEntries(*gvec);
238: }
239: VecSetDM(*gvec,dm);
240: return(0);
241: }
243: PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec)
244: {
246: DM_Shell *shell = (DM_Shell*)dm->data;
247: Vec X;
252: *gvec = 0;
253: X = shell->Xlocal;
254: if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()");
255: if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */
256: PetscObjectReference((PetscObject)X);
257: VecZeroEntries(X);
258: *gvec = X;
259: } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */
260: VecDuplicate(X,gvec);
261: VecZeroEntries(*gvec);
262: }
263: VecSetDM(*gvec,dm);
264: return(0);
265: }
267: /*@
268: DMShellSetContext - set some data to be usable by this DM
270: Collective
272: Input Arguments:
273: + dm - shell DM
274: - ctx - the context
276: Level: advanced
278: .seealso: DMCreateMatrix(), DMShellGetContext()
279: @*/
280: PetscErrorCode DMShellSetContext(DM dm,void *ctx)
281: {
282: DM_Shell *shell = (DM_Shell*)dm->data;
284: PetscBool isshell;
288: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
289: if (!isshell) return(0);
290: shell->ctx = ctx;
291: return(0);
292: }
294: /*@
295: DMShellGetContext - set some data to be usable by this DM
297: Collective
299: Input Argument:
300: . dm - shell DM
302: Output Argument:
303: . ctx - the context
305: Level: advanced
307: .seealso: DMCreateMatrix(), DMShellSetContext()
308: @*/
309: PetscErrorCode DMShellGetContext(DM dm,void **ctx)
310: {
311: DM_Shell *shell = (DM_Shell*)dm->data;
313: PetscBool isshell;
317: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
318: if (!isshell) return(0);
319: *ctx = shell->ctx;
320: return(0);
321: }
323: /*@
324: DMShellSetMatrix - sets a template matrix associated with the DMShell
326: Collective
328: Input Arguments:
329: + dm - shell DM
330: - J - template matrix
332: Level: advanced
334: .seealso: DMCreateMatrix(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
335: @*/
336: PetscErrorCode DMShellSetMatrix(DM dm,Mat J)
337: {
338: DM_Shell *shell = (DM_Shell*)dm->data;
340: PetscBool isshell;
345: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
346: if (!isshell) return(0);
347: PetscObjectReference((PetscObject)J);
348: MatDestroy(&shell->A);
349: shell->A = J;
350: return(0);
351: }
353: /*@C
354: DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM
356: Logically Collective on DM
358: Input Arguments:
359: + dm - the shell DM
360: - func - the function to create a matrix
362: Level: advanced
364: .seealso: DMCreateMatrix(), DMShellSetMatrix(), DMShellSetContext(), DMShellGetContext()
365: @*/
366: PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*))
367: {
371: dm->ops->creatematrix = func;
372: return(0);
373: }
375: /*@
376: DMShellSetGlobalVector - sets a template global vector associated with the DMShell
378: Logically Collective on DM
380: Input Arguments:
381: + dm - shell DM
382: - X - template vector
384: Level: advanced
386: .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector()
387: @*/
388: PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X)
389: {
390: DM_Shell *shell = (DM_Shell*)dm->data;
392: PetscBool isshell;
393: DM vdm;
398: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
399: if (!isshell) return(0);
400: VecGetDM(X,&vdm);
401: /*
402: if the vector proposed as the new base global vector for the DM is a DM vector associated
403: with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
404: we get a circular dependency that prevents the DM from being destroy when it should be.
405: This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
406: DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
407: to set its input vector (which is associated with the DM) as the base global vector.
408: Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
409: for pointing out the problem.
410: */
411: if (vdm == dm) return(0);
412: PetscObjectReference((PetscObject)X);
413: VecDestroy(&shell->Xglobal);
414: shell->Xglobal = X;
415: return(0);
416: }
418: /*@C
419: DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM
421: Logically Collective
423: Input Arguments:
424: + dm - the shell DM
425: - func - the creation routine
427: Level: advanced
429: .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
430: @*/
431: PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
432: {
436: dm->ops->createglobalvector = func;
437: return(0);
438: }
440: /*@
441: DMShellSetLocalVector - sets a template local vector associated with the DMShell
443: Logically Collective on DM
445: Input Arguments:
446: + dm - shell DM
447: - X - template vector
449: Level: advanced
451: .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector()
452: @*/
453: PetscErrorCode DMShellSetLocalVector(DM dm,Vec X)
454: {
455: DM_Shell *shell = (DM_Shell*)dm->data;
457: PetscBool isshell;
458: DM vdm;
463: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
464: if (!isshell) return(0);
465: VecGetDM(X,&vdm);
466: /*
467: if the vector proposed as the new base global vector for the DM is a DM vector associated
468: with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
469: we get a circular dependency that prevents the DM from being destroy when it should be.
470: This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
471: DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
472: to set its input vector (which is associated with the DM) as the base global vector.
473: Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
474: for pointing out the problem.
475: */
476: if (vdm == dm) return(0);
477: PetscObjectReference((PetscObject)X);
478: VecDestroy(&shell->Xlocal);
479: shell->Xlocal = X;
480: return(0);
481: }
483: /*@C
484: DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM
486: Logically Collective
488: Input Arguments:
489: + dm - the shell DM
490: - func - the creation routine
492: Level: advanced
494: .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
495: @*/
496: PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
497: {
501: dm->ops->createlocalvector = func;
502: return(0);
503: }
505: /*@C
506: DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter
508: Logically Collective on DM
510: Input Arguments
511: + dm - the shell DM
512: . begin - the routine that begins the global to local scatter
513: - end - the routine that ends the global to local scatter
515: Notes:
516: If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then
517: DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers
519: Level: advanced
521: .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
522: @*/
523: PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
525: dm->ops->globaltolocalbegin = begin;
526: dm->ops->globaltolocalend = end;
527: return(0);
528: }
530: /*@C
531: DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter
533: Logically Collective on DM
535: Input Arguments
536: + dm - the shell DM
537: . begin - the routine that begins the local to global scatter
538: - end - the routine that ends the local to global scatter
540: Notes:
541: If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then
542: DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers
544: Level: advanced
546: .seealso: DMShellSetGlobalToLocal()
547: @*/
548: PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
550: dm->ops->localtoglobalbegin = begin;
551: dm->ops->localtoglobalend = end;
552: return(0);
553: }
555: /*@C
556: DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter
558: Logically Collective on DM
560: Input Arguments
561: + dm - the shell DM
562: . begin - the routine that begins the local to local scatter
563: - end - the routine that ends the local to local scatter
565: Notes:
566: If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then
567: DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers
569: Level: advanced
571: .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
572: @*/
573: PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
575: dm->ops->localtolocalbegin = begin;
576: dm->ops->localtolocalend = end;
577: return(0);
578: }
580: /*@
581: DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication
583: Logically Collective on DM
585: Input Arguments
586: + dm - the shell DM
587: - gtol - the global to local VecScatter context
589: Level: advanced
591: .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
592: @*/
593: PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
594: {
595: DM_Shell *shell = (DM_Shell*)dm->data;
599: PetscObjectReference((PetscObject)gtol);
600: /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */
601: VecScatterDestroy(&shell->gtol);
602: shell->gtol = gtol;
603: return(0);
604: }
606: /*@
607: DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication
609: Logically Collective on DM
611: Input Arguments
612: + dm - the shell DM
613: - ltog - the local to global VecScatter context
615: Level: advanced
617: .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell()
618: @*/
619: PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
620: {
621: DM_Shell *shell = (DM_Shell*)dm->data;
625: PetscObjectReference((PetscObject)ltog);
626: /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */
627: VecScatterDestroy(&shell->ltog);
628: shell->ltog = ltog;
629: return(0);
630: }
632: /*@
633: DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication
635: Logically Collective on DM
637: Input Arguments
638: + dm - the shell DM
639: - ltol - the local to local VecScatter context
641: Level: advanced
643: .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
644: @*/
645: PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
646: {
647: DM_Shell *shell = (DM_Shell*)dm->data;
651: PetscObjectReference((PetscObject)ltol);
652: /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */
653: VecScatterDestroy(&shell->ltol);
654: shell->ltol = ltol;
655: return(0);
656: }
658: /*@C
659: DMShellSetCoarsen - Set the routine used to coarsen the shell DM
661: Logically Collective on DM
663: Input Arguments
664: + dm - the shell DM
665: - coarsen - the routine that coarsens the DM
667: Level: advanced
669: .seealso: DMShellSetRefine(), DMCoarsen(), DMShellSetContext(), DMShellGetContext()
670: @*/
671: PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*))
672: {
674: PetscBool isshell;
678: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
679: if (!isshell) return(0);
680: dm->ops->coarsen = coarsen;
681: return(0);
682: }
684: /*@C
685: DMShellSetRefine - Set the routine used to refine the shell DM
687: Logically Collective on DM
689: Input Arguments
690: + dm - the shell DM
691: - refine - the routine that refines the DM
693: Level: advanced
695: .seealso: DMShellSetCoarsen(), DMRefine(), DMShellSetContext(), DMShellGetContext()
696: @*/
697: PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*))
698: {
700: PetscBool isshell;
704: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
705: if (!isshell) return(0);
706: dm->ops->refine = refine;
707: return(0);
708: }
710: /*@C
711: DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator
713: Logically Collective on DM
715: Input Arguments
716: + dm - the shell DM
717: - interp - the routine to create the interpolation
719: Level: advanced
721: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
722: @*/
723: PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*))
724: {
726: PetscBool isshell;
730: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
731: if (!isshell) return(0);
732: dm->ops->createinterpolation = interp;
733: return(0);
734: }
736: /*@C
737: DMShellSetCreateRestriction - Set the routine used to create the restriction operator
739: Logically Collective on DM
741: Input Arguments
742: + dm - the shell DM
743: - striction- the routine to create the restriction
745: Level: advanced
747: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext()
748: @*/
749: PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*))
750: {
752: PetscBool isshell;
756: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
757: if (!isshell) return(0);
758: dm->ops->createrestriction = restriction;
759: return(0);
760: }
762: /*@C
763: DMShellSetCreateInjection - Set the routine used to create the injection operator
765: Logically Collective on DM
767: Input Arguments
768: + dm - the shell DM
769: - inject - the routine to create the injection
771: Level: advanced
773: .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext()
774: @*/
775: PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*))
776: {
778: PetscBool isshell;
782: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
783: if (!isshell) return(0);
784: dm->ops->getinjection = inject;
785: return(0);
786: }
788: static PetscErrorCode DMHasCreateInjection_Shell(DM dm, PetscBool *flg)
789: {
791: PetscBool isshell;
796: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
797: if (!isshell) {
798: *flg = PETSC_FALSE;
799: } else {
800: *flg = dm->ops->getinjection ? PETSC_TRUE : PETSC_FALSE;
801: }
802: return(0);
803: }
804: /*@C
805: DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM
807: Logically Collective on DM
809: Input Arguments
810: + dm - the shell DM
811: - decomp - the routine to create the decomposition
813: Level: advanced
815: .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext()
816: @*/
817: PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**))
818: {
820: PetscBool isshell;
824: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
825: if (!isshell) return(0);
826: dm->ops->createfielddecomposition = decomp;
827: return(0);
828: }
830: /*@C
831: DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM
833: Logically Collective on DM
835: Input Arguments
836: + dm - the shell DM
837: - decomp - the routine to create the decomposition
839: Level: advanced
841: .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext()
842: @*/
843: PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**))
844: {
846: PetscBool isshell;
850: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
851: if (!isshell) return(0);
852: dm->ops->createdomaindecomposition = decomp;
853: return(0);
854: }
856: /*@C
857: DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM
859: Logically Collective on DM
861: Input Arguments
862: + dm - the shell DM
863: - scatter - the routine to create the scatters
865: Level: advanced
867: .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext()
868: @*/
869: PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**))
870: {
872: PetscBool isshell;
876: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
877: if (!isshell) return(0);
878: dm->ops->createddscatters = scatter;
879: return(0);
880: }
882: /*@C
883: DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM
885: Logically Collective on DM
887: Input Arguments
888: + dm - the shell DM
889: - subdm - the routine to create the decomposition
891: Level: advanced
893: .seealso: DMCreateSubDM(), DMShellSetContext(), DMShellGetContext()
894: @*/
895: PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
896: {
898: PetscBool isshell;
902: PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
903: if (!isshell) return(0);
904: dm->ops->createsubdm = subdm;
905: return(0);
906: }
908: static PetscErrorCode DMDestroy_Shell(DM dm)
909: {
911: DM_Shell *shell = (DM_Shell*)dm->data;
914: MatDestroy(&shell->A);
915: VecDestroy(&shell->Xglobal);
916: VecDestroy(&shell->Xlocal);
917: VecScatterDestroy(&shell->gtol);
918: VecScatterDestroy(&shell->ltog);
919: VecScatterDestroy(&shell->ltol);
920: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
921: PetscFree(shell);
922: return(0);
923: }
925: static PetscErrorCode DMView_Shell(DM dm,PetscViewer v)
926: {
928: DM_Shell *shell = (DM_Shell*)dm->data;
931: VecView(shell->Xglobal,v);
932: return(0);
933: }
935: static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v)
936: {
938: DM_Shell *shell = (DM_Shell*)dm->data;
941: VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);
942: VecLoad(shell->Xglobal,v);
943: return(0);
944: }
946: PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
947: {
951: if (subdm) {DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);}
952: DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);
953: return(0);
954: }
956: PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
957: {
959: DM_Shell *shell;
962: PetscNewLog(dm,&shell);
963: dm->data = shell;
965: PetscObjectChangeTypeName((PetscObject)dm,DMSHELL);
967: dm->ops->destroy = DMDestroy_Shell;
968: dm->ops->createglobalvector = DMCreateGlobalVector_Shell;
969: dm->ops->createlocalvector = DMCreateLocalVector_Shell;
970: dm->ops->creatematrix = DMCreateMatrix_Shell;
971: dm->ops->view = DMView_Shell;
972: dm->ops->load = DMLoad_Shell;
973: dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell;
974: dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell;
975: dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell;
976: dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell;
977: dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell;
978: dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell;
979: dm->ops->createsubdm = DMCreateSubDM_Shell;
980: dm->ops->hascreateinjection = DMHasCreateInjection_Shell;
981: return(0);
982: }
984: /*@
985: DMShellCreate - Creates a shell DM object, used to manage user-defined problem data
987: Collective on MPI_Comm
989: Input Parameter:
990: . comm - the processors that will share the global vector
992: Output Parameters:
993: . shell - the shell DM
995: Level: advanced
997: .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext()
998: @*/
999: PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm)
1000: {
1005: DMCreate(comm,dm);
1006: DMSetType(*dm,DMSHELL);
1007: DMSetUp(*dm);
1008: return(0);
1009: }