Actual source code: dmshell.c
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: PetscErrorCode (*destroyctx)(void *);
14: } 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 Parameters:
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: {
34: DM_Shell *shell = (DM_Shell *)dm->data;
36: PetscFunctionBegin;
37: PetscCheck(shell->gtol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
38: PetscCall(VecScatterBegin(shell->gtol, g, l, mode, SCATTER_FORWARD));
39: PetscFunctionReturn(PETSC_SUCCESS);
40: }
42: /*@
43: DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to end a global to local scatter
44: Collective
46: Input Parameters:
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: {
58: DM_Shell *shell = (DM_Shell *)dm->data;
60: PetscFunctionBegin;
61: PetscCheck(shell->gtol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
62: PetscCall(VecScatterEnd(shell->gtol, g, l, mode, SCATTER_FORWARD));
63: PetscFunctionReturn(PETSC_SUCCESS);
64: }
66: /*@
67: DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to begin a local to global scatter
68: Collective
70: Input Parameters:
71: + dm - shell DM
72: . l - local vector
73: . mode - InsertMode
74: - g - global vector
76: Level: advanced
78: 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.
80: .seealso: `DMLocalToGlobalEndDefaultShell()`
81: @*/
82: PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm, Vec l, InsertMode mode, Vec g)
83: {
84: DM_Shell *shell = (DM_Shell *)dm->data;
86: PetscFunctionBegin;
87: PetscCheck(shell->ltog, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
88: PetscCall(VecScatterBegin(shell->ltog, l, g, mode, SCATTER_FORWARD));
89: PetscFunctionReturn(PETSC_SUCCESS);
90: }
92: /*@
93: DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to end a local to global scatter
94: Collective
96: Input Parameters:
97: + dm - shell DM
98: . l - local vector
99: . mode - InsertMode
100: - g - global vector
102: Level: advanced
104: .seealso: `DMLocalToGlobalBeginDefaultShell()`
105: @*/
106: PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm, Vec l, InsertMode mode, Vec g)
107: {
108: DM_Shell *shell = (DM_Shell *)dm->data;
110: PetscFunctionBegin;
111: PetscCheck(shell->ltog, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
112: PetscCall(VecScatterEnd(shell->ltog, l, g, mode, SCATTER_FORWARD));
113: PetscFunctionReturn(PETSC_SUCCESS);
114: }
116: /*@
117: DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal VecScatter context set by the user to begin a local to local scatter
118: Collective
120: Input Parameters:
121: + dm - shell DM
122: . g - the original local vector
123: - mode - InsertMode
125: Output Parameter:
126: . l - the local vector with correct ghost values
128: Level: advanced
130: 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.
132: .seealso: `DMLocalToLocalEndDefaultShell()`
133: @*/
134: PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm, Vec g, InsertMode mode, Vec l)
135: {
136: DM_Shell *shell = (DM_Shell *)dm->data;
138: PetscFunctionBegin;
139: PetscCheck(shell->ltol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()");
140: PetscCall(VecScatterBegin(shell->ltol, g, l, mode, SCATTER_FORWARD));
141: PetscFunctionReturn(PETSC_SUCCESS);
142: }
144: /*@
145: DMLocalToLocalEndDefaultShell - Uses the LocalToLocal VecScatter context set by the user to end a local to local scatter
146: Collective
148: Input Parameters:
149: + dm - shell DM
150: . g - the original local vector
151: - mode - InsertMode
153: Output Parameter:
154: . l - the local vector with correct ghost values
156: Level: advanced
158: .seealso: `DMLocalToLocalBeginDefaultShell()`
159: @*/
160: PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm, Vec g, InsertMode mode, Vec l)
161: {
162: DM_Shell *shell = (DM_Shell *)dm->data;
164: PetscFunctionBegin;
165: PetscCheck(shell->ltol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
166: PetscCall(VecScatterEnd(shell->ltol, g, l, mode, SCATTER_FORWARD));
167: PetscFunctionReturn(PETSC_SUCCESS);
168: }
170: static PetscErrorCode DMCreateMatrix_Shell(DM dm, Mat *J)
171: {
172: DM_Shell *shell = (DM_Shell *)dm->data;
173: Mat A;
175: PetscFunctionBegin;
178: if (!shell->A) {
179: if (shell->Xglobal) {
180: PetscInt m, M;
181: PetscCall(PetscInfo(dm, "Naively creating matrix using global vector distribution without preallocation\n"));
182: PetscCall(VecGetSize(shell->Xglobal, &M));
183: PetscCall(VecGetLocalSize(shell->Xglobal, &m));
184: PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), &shell->A));
185: PetscCall(MatSetSizes(shell->A, m, m, M, M));
186: PetscCall(MatSetType(shell->A, dm->mattype));
187: PetscCall(MatSetUp(shell->A));
188: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector");
189: }
190: A = shell->A;
191: PetscCall(MatDuplicate(A, MAT_SHARE_NONZERO_PATTERN, J));
192: PetscCall(MatSetDM(*J, dm));
193: PetscFunctionReturn(PETSC_SUCCESS);
194: }
196: PetscErrorCode DMCreateGlobalVector_Shell(DM dm, Vec *gvec)
197: {
198: DM_Shell *shell = (DM_Shell *)dm->data;
199: Vec X;
201: PetscFunctionBegin;
204: *gvec = NULL;
205: X = shell->Xglobal;
206: PetscCheck(X, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()");
207: /* Need to create a copy in order to attach the DM to the vector */
208: PetscCall(VecDuplicate(X, gvec));
209: PetscCall(VecZeroEntries(*gvec));
210: PetscCall(VecSetDM(*gvec, dm));
211: PetscFunctionReturn(PETSC_SUCCESS);
212: }
214: PetscErrorCode DMCreateLocalVector_Shell(DM dm, Vec *gvec)
215: {
216: DM_Shell *shell = (DM_Shell *)dm->data;
217: Vec X;
219: PetscFunctionBegin;
222: *gvec = NULL;
223: X = shell->Xlocal;
224: PetscCheck(X, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()");
225: /* Need to create a copy in order to attach the DM to the vector */
226: PetscCall(VecDuplicate(X, gvec));
227: PetscCall(VecZeroEntries(*gvec));
228: PetscCall(VecSetDM(*gvec, dm));
229: PetscFunctionReturn(PETSC_SUCCESS);
230: }
232: /*@
233: DMShellSetDestroyContext - set a function that destroys the context provided with `DMShellSetContext()`
235: Collective
237: Input Parameters:
238: . ctx - the context
240: Level: advanced
242: .seealso: `DMShellSetContext()`, `DMShellGetContext()`
243: @*/
244: PetscErrorCode DMShellSetDestroyContext(DM dm, PetscErrorCode (*destroyctx)(void *))
245: {
246: DM_Shell *shell = (DM_Shell *)dm->data;
247: PetscBool isshell;
249: PetscFunctionBegin;
251: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
252: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
253: shell->destroyctx = destroyctx;
254: PetscFunctionReturn(PETSC_SUCCESS);
255: }
257: /*@
258: DMShellSetContext - set some data to be usable by this DM
260: Collective
262: Input Parameters:
263: + dm - shell DM
264: - ctx - the context
266: Level: advanced
268: .seealso: `DMCreateMatrix()`, `DMShellGetContext()`
269: @*/
270: PetscErrorCode DMShellSetContext(DM dm, void *ctx)
271: {
272: DM_Shell *shell = (DM_Shell *)dm->data;
273: PetscBool isshell;
275: PetscFunctionBegin;
277: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
278: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
279: shell->ctx = ctx;
280: PetscFunctionReturn(PETSC_SUCCESS);
281: }
283: /*@
284: DMShellGetContext - Returns the user-provided context associated to the DM
286: Collective
288: Input Parameter:
289: . dm - shell DM
291: Output Parameter:
292: . ctx - the context
294: Level: advanced
296: .seealso: `DMCreateMatrix()`, `DMShellSetContext()`
297: @*/
298: PetscErrorCode DMShellGetContext(DM dm, void *ctx)
299: {
300: DM_Shell *shell = (DM_Shell *)dm->data;
301: PetscBool isshell;
303: PetscFunctionBegin;
305: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
306: PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
307: *(void **)ctx = shell->ctx;
308: PetscFunctionReturn(PETSC_SUCCESS);
309: }
311: /*@
312: DMShellSetMatrix - sets a template matrix associated with the DMShell
314: Collective
316: Input Parameters:
317: + dm - shell DM
318: - J - template matrix
320: Level: advanced
322: Developer Notes:
323: 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.
325: .seealso: `DMCreateMatrix()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
326: @*/
327: PetscErrorCode DMShellSetMatrix(DM dm, Mat J)
328: {
329: DM_Shell *shell = (DM_Shell *)dm->data;
330: PetscBool isshell;
331: DM mdm;
333: PetscFunctionBegin;
336: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
337: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
338: if (J == shell->A) PetscFunctionReturn(PETSC_SUCCESS);
339: PetscCall(MatGetDM(J, &mdm));
340: PetscCall(PetscObjectReference((PetscObject)J));
341: PetscCall(MatDestroy(&shell->A));
342: if (mdm == dm) {
343: PetscCall(MatDuplicate(J, MAT_SHARE_NONZERO_PATTERN, &shell->A));
344: PetscCall(MatSetDM(shell->A, NULL));
345: } else shell->A = J;
346: PetscFunctionReturn(PETSC_SUCCESS);
347: }
349: /*@C
350: DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM
352: Logically Collective on dm
354: Input Parameters:
355: + dm - the shell DM
356: - func - the function to create a matrix
358: Level: advanced
360: .seealso: `DMCreateMatrix()`, `DMShellSetMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
361: @*/
362: PetscErrorCode DMShellSetCreateMatrix(DM dm, PetscErrorCode (*func)(DM, Mat *))
363: {
364: PetscFunctionBegin;
366: dm->ops->creatematrix = func;
367: PetscFunctionReturn(PETSC_SUCCESS);
368: }
370: /*@
371: DMShellSetGlobalVector - sets a template global vector associated with the DMShell
373: Logically Collective on dm
375: Input Parameters:
376: + dm - shell DM
377: - X - template vector
379: Level: advanced
381: .seealso: `DMCreateGlobalVector()`, `DMShellSetMatrix()`, `DMShellSetCreateGlobalVector()`
382: @*/
383: PetscErrorCode DMShellSetGlobalVector(DM dm, Vec X)
384: {
385: DM_Shell *shell = (DM_Shell *)dm->data;
386: PetscBool isshell;
387: DM vdm;
389: PetscFunctionBegin;
392: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
393: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
394: PetscCall(VecGetDM(X, &vdm));
395: /*
396: if the vector proposed as the new base global vector for the DM is a DM vector associated
397: with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
398: we get a circular dependency that prevents the DM from being destroy when it should be.
399: This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
400: DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
401: to set its input vector (which is associated with the DM) as the base global vector.
402: Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
403: for pointing out the problem.
404: */
405: if (vdm == dm) PetscFunctionReturn(PETSC_SUCCESS);
406: PetscCall(PetscObjectReference((PetscObject)X));
407: PetscCall(VecDestroy(&shell->Xglobal));
408: shell->Xglobal = X;
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: /*@
413: DMShellGetGlobalVector - Returns the template global vector associated with the DMShell, or NULL if it was not set
415: Not collective
417: Input Parameters:
418: + dm - shell DM
419: - X - template vector
421: Level: advanced
423: .seealso: `DMShellSetGlobalVector()`, `DMShellSetCreateGlobalVector()`, `DMCreateGlobalVector()`
424: @*/
425: PetscErrorCode DMShellGetGlobalVector(DM dm, Vec *X)
426: {
427: DM_Shell *shell = (DM_Shell *)dm->data;
428: PetscBool isshell;
430: PetscFunctionBegin;
433: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
434: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
435: *X = shell->Xglobal;
436: PetscFunctionReturn(PETSC_SUCCESS);
437: }
439: /*@C
440: DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM
442: Logically Collective
444: Input Parameters:
445: + dm - the shell DM
446: - func - the creation routine
448: Level: advanced
450: .seealso: `DMShellSetGlobalVector()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
451: @*/
452: PetscErrorCode DMShellSetCreateGlobalVector(DM dm, PetscErrorCode (*func)(DM, Vec *))
453: {
454: PetscFunctionBegin;
456: dm->ops->createglobalvector = func;
457: PetscFunctionReturn(PETSC_SUCCESS);
458: }
460: /*@
461: DMShellSetLocalVector - sets a template local vector associated with the DMShell
463: Logically Collective on dm
465: Input Parameters:
466: + dm - shell DM
467: - X - template vector
469: Level: advanced
471: .seealso: `DMCreateLocalVector()`, `DMShellSetMatrix()`, `DMShellSetCreateLocalVector()`
472: @*/
473: PetscErrorCode DMShellSetLocalVector(DM dm, Vec X)
474: {
475: DM_Shell *shell = (DM_Shell *)dm->data;
476: PetscBool isshell;
477: DM vdm;
479: PetscFunctionBegin;
482: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
483: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
484: PetscCall(VecGetDM(X, &vdm));
485: /*
486: if the vector proposed as the new base global vector for the DM is a DM vector associated
487: with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
488: we get a circular dependency that prevents the DM from being destroy when it should be.
489: This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
490: DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
491: to set its input vector (which is associated with the DM) as the base global vector.
492: Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
493: for pointing out the problem.
494: */
495: if (vdm == dm) PetscFunctionReturn(PETSC_SUCCESS);
496: PetscCall(PetscObjectReference((PetscObject)X));
497: PetscCall(VecDestroy(&shell->Xlocal));
498: shell->Xlocal = X;
499: PetscFunctionReturn(PETSC_SUCCESS);
500: }
502: /*@C
503: DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM
505: Logically Collective
507: Input Parameters:
508: + dm - the shell DM
509: - func - the creation routine
511: Level: advanced
513: .seealso: `DMShellSetLocalVector()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
514: @*/
515: PetscErrorCode DMShellSetCreateLocalVector(DM dm, PetscErrorCode (*func)(DM, Vec *))
516: {
517: PetscFunctionBegin;
519: dm->ops->createlocalvector = func;
520: PetscFunctionReturn(PETSC_SUCCESS);
521: }
523: /*@C
524: DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter
526: Logically Collective on dm
528: Input Parameters
529: + dm - the shell DM
530: . begin - the routine that begins the global to local scatter
531: - end - the routine that ends the global to local scatter
533: Notes:
534: If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then
535: DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers
537: Level: advanced
539: .seealso: `DMShellSetLocalToGlobal()`, `DMGlobalToLocalBeginDefaultShell()`, `DMGlobalToLocalEndDefaultShell()`
540: @*/
541: PetscErrorCode DMShellSetGlobalToLocal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec))
542: {
543: PetscFunctionBegin;
545: dm->ops->globaltolocalbegin = begin;
546: dm->ops->globaltolocalend = end;
547: PetscFunctionReturn(PETSC_SUCCESS);
548: }
550: /*@C
551: DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter
553: Logically Collective on dm
555: Input Parameters
556: + dm - the shell DM
557: . begin - the routine that begins the local to global scatter
558: - end - the routine that ends the local to global scatter
560: Notes:
561: If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then
562: DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers
564: Level: advanced
566: .seealso: `DMShellSetGlobalToLocal()`
567: @*/
568: PetscErrorCode DMShellSetLocalToGlobal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec))
569: {
570: PetscFunctionBegin;
572: dm->ops->localtoglobalbegin = begin;
573: dm->ops->localtoglobalend = end;
574: PetscFunctionReturn(PETSC_SUCCESS);
575: }
577: /*@C
578: DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter
580: Logically Collective on dm
582: Input Parameters
583: + dm - the shell DM
584: . begin - the routine that begins the local to local scatter
585: - end - the routine that ends the local to local scatter
587: Notes:
588: If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then
589: DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers
591: Level: advanced
593: .seealso: `DMShellSetGlobalToLocal()`, `DMLocalToLocalBeginDefaultShell()`, `DMLocalToLocalEndDefaultShell()`
594: @*/
595: PetscErrorCode DMShellSetLocalToLocal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec))
596: {
597: PetscFunctionBegin;
599: dm->ops->localtolocalbegin = begin;
600: dm->ops->localtolocalend = end;
601: PetscFunctionReturn(PETSC_SUCCESS);
602: }
604: /*@
605: DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication
607: Logically Collective on dm
609: Input Parameters
610: + dm - the shell DM
611: - gtol - the global to local VecScatter context
613: Level: advanced
615: .seealso: `DMShellSetGlobalToLocal()`, `DMGlobalToLocalBeginDefaultShell()`, `DMGlobalToLocalEndDefaultShell()`
616: @*/
617: PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
618: {
619: DM_Shell *shell = (DM_Shell *)dm->data;
621: PetscFunctionBegin;
624: PetscCall(PetscObjectReference((PetscObject)gtol));
625: PetscCall(VecScatterDestroy(&shell->gtol));
626: shell->gtol = gtol;
627: PetscFunctionReturn(PETSC_SUCCESS);
628: }
630: /*@
631: DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication
633: Logically Collective on dm
635: Input Parameters
636: + dm - the shell DM
637: - ltog - the local to global VecScatter context
639: Level: advanced
641: .seealso: `DMShellSetLocalToGlobal()`, `DMLocalToGlobalBeginDefaultShell()`, `DMLocalToGlobalEndDefaultShell()`
642: @*/
643: PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
644: {
645: DM_Shell *shell = (DM_Shell *)dm->data;
647: PetscFunctionBegin;
650: PetscCall(PetscObjectReference((PetscObject)ltog));
651: PetscCall(VecScatterDestroy(&shell->ltog));
652: shell->ltog = ltog;
653: PetscFunctionReturn(PETSC_SUCCESS);
654: }
656: /*@
657: DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication
659: Logically Collective on dm
661: Input Parameters
662: + dm - the shell DM
663: - ltol - the local to local VecScatter context
665: Level: advanced
667: .seealso: `DMShellSetLocalToLocal()`, `DMLocalToLocalBeginDefaultShell()`, `DMLocalToLocalEndDefaultShell()`
668: @*/
669: PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
670: {
671: DM_Shell *shell = (DM_Shell *)dm->data;
673: PetscFunctionBegin;
676: PetscCall(PetscObjectReference((PetscObject)ltol));
677: PetscCall(VecScatterDestroy(&shell->ltol));
678: shell->ltol = ltol;
679: PetscFunctionReturn(PETSC_SUCCESS);
680: }
682: /*@C
683: DMShellSetCoarsen - Set the routine used to coarsen the shell DM
685: Logically Collective on dm
687: Input Parameters
688: + dm - the shell DM
689: - coarsen - the routine that coarsens the DM
691: Level: advanced
693: .seealso: `DMShellSetRefine()`, `DMCoarsen()`, `DMShellGetCoarsen()`, `DMShellSetContext()`, `DMShellGetContext()`
694: @*/
695: PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *))
696: {
697: PetscBool isshell;
699: PetscFunctionBegin;
701: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
702: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
703: dm->ops->coarsen = coarsen;
704: PetscFunctionReturn(PETSC_SUCCESS);
705: }
707: /*@C
708: DMShellGetCoarsen - Get the routine used to coarsen the shell DM
710: Logically Collective on dm
712: Input Parameter:
713: . dm - the shell DM
715: Output Parameter:
716: . coarsen - the routine that coarsens the DM
718: Level: advanced
720: .seealso: `DMShellSetCoarsen()`, `DMCoarsen()`, `DMShellSetRefine()`, `DMRefine()`
721: @*/
722: PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM, MPI_Comm, DM *))
723: {
724: PetscBool isshell;
726: PetscFunctionBegin;
728: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
729: PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
730: *coarsen = dm->ops->coarsen;
731: PetscFunctionReturn(PETSC_SUCCESS);
732: }
734: /*@C
735: DMShellSetRefine - Set the routine used to refine the shell DM
737: Logically Collective on dm
739: Input Parameters
740: + dm - the shell DM
741: - refine - the routine that refines the DM
743: Level: advanced
745: .seealso: `DMShellSetCoarsen()`, `DMRefine()`, `DMShellGetRefine()`, `DMShellSetContext()`, `DMShellGetContext()`
746: @*/
747: PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM, MPI_Comm, DM *))
748: {
749: PetscBool isshell;
751: PetscFunctionBegin;
753: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
754: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
755: dm->ops->refine = refine;
756: PetscFunctionReturn(PETSC_SUCCESS);
757: }
759: /*@C
760: DMShellGetRefine - Get the routine used to refine the shell DM
762: Logically Collective on dm
764: Input Parameter:
765: . dm - the shell DM
767: Output Parameter:
768: . refine - the routine that refines the DM
770: Level: advanced
772: .seealso: `DMShellSetCoarsen()`, `DMCoarsen()`, `DMShellSetRefine()`, `DMRefine()`
773: @*/
774: PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM, MPI_Comm, DM *))
775: {
776: PetscBool isshell;
778: PetscFunctionBegin;
780: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
781: PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
782: *refine = dm->ops->refine;
783: PetscFunctionReturn(PETSC_SUCCESS);
784: }
786: /*@C
787: DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator
789: Logically Collective on dm
791: Input Parameters
792: + dm - the shell DM
793: - interp - the routine to create the interpolation
795: Level: advanced
797: .seealso: `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateInterpolation()`, `DMShellSetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()`
798: @*/
799: PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM, DM, Mat *, Vec *))
800: {
801: PetscBool isshell;
803: PetscFunctionBegin;
805: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
806: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
807: dm->ops->createinterpolation = interp;
808: PetscFunctionReturn(PETSC_SUCCESS);
809: }
811: /*@C
812: DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator
814: Logically Collective on dm
816: Input Parameter:
817: . dm - the shell DM
819: Output Parameter:
820: . interp - the routine to create the interpolation
822: Level: advanced
824: .seealso: `DMShellGetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()`
825: @*/
826: PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM, DM, Mat *, Vec *))
827: {
828: PetscBool isshell;
830: PetscFunctionBegin;
832: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
833: PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
834: *interp = dm->ops->createinterpolation;
835: PetscFunctionReturn(PETSC_SUCCESS);
836: }
838: /*@C
839: DMShellSetCreateRestriction - Set the routine used to create the restriction operator
841: Logically Collective on dm
843: Input Parameters
844: + dm - the shell DM
845: - striction- the routine to create the restriction
847: Level: advanced
849: .seealso: `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()`
850: @*/
851: PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM, DM, Mat *))
852: {
853: PetscBool isshell;
855: PetscFunctionBegin;
857: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
858: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
859: dm->ops->createrestriction = restriction;
860: PetscFunctionReturn(PETSC_SUCCESS);
861: }
863: /*@C
864: DMShellGetCreateRestriction - Get the routine used to create the restriction operator
866: Logically Collective on dm
868: Input Parameter:
869: . dm - the shell DM
871: Output Parameter:
872: . restriction - the routine to create the restriction
874: Level: advanced
876: .seealso: `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellSetContext()`, `DMShellGetContext()`
877: @*/
878: PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM, DM, Mat *))
879: {
880: PetscBool isshell;
882: PetscFunctionBegin;
884: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
885: PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
886: *restriction = dm->ops->createrestriction;
887: PetscFunctionReturn(PETSC_SUCCESS);
888: }
890: /*@C
891: DMShellSetCreateInjection - Set the routine used to create the injection operator
893: Logically Collective on dm
895: Input Parameters:
896: + dm - the shell DM
897: - inject - the routine to create the injection
899: Level: advanced
901: .seealso: `DMShellSetCreateInterpolation()`, `DMCreateInjection()`, `DMShellGetCreateInjection()`, `DMShellSetContext()`, `DMShellGetContext()`
902: @*/
903: PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM, DM, Mat *))
904: {
905: PetscBool isshell;
907: PetscFunctionBegin;
909: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
910: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
911: dm->ops->createinjection = inject;
912: PetscFunctionReturn(PETSC_SUCCESS);
913: }
915: /*@C
916: DMShellGetCreateInjection - Get the routine used to create the injection operator
918: Logically Collective on dm
920: Input Parameter:
921: . dm - the shell DM
923: Output Parameter:
924: . inject - the routine to create the injection
926: Level: advanced
928: .seealso: `DMShellGetCreateInterpolation()`, `DMCreateInjection()`, `DMShellSetContext()`, `DMShellGetContext()`
929: @*/
930: PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM, DM, Mat *))
931: {
932: PetscBool isshell;
934: PetscFunctionBegin;
936: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
937: PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
938: *inject = dm->ops->createinjection;
939: PetscFunctionReturn(PETSC_SUCCESS);
940: }
942: /*@C
943: DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM
945: Logically Collective on dm
947: Input Parameters:
948: + dm - the shell DM
949: - decomp - the routine to create the decomposition
951: Level: advanced
953: .seealso: `DMCreateFieldDecomposition()`, `DMShellSetContext()`, `DMShellGetContext()`
954: @*/
955: PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM, PetscInt *, char ***, IS **, DM **))
956: {
957: PetscBool isshell;
959: PetscFunctionBegin;
961: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
962: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
963: dm->ops->createfielddecomposition = decomp;
964: PetscFunctionReturn(PETSC_SUCCESS);
965: }
967: /*@C
968: DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM
970: Logically Collective on dm
972: Input Parameters:
973: + dm - the shell DM
974: - decomp - the routine to create the decomposition
976: Level: advanced
978: .seealso: `DMCreateDomainDecomposition()`, `DMShellSetContext()`, `DMShellGetContext()`
979: @*/
980: PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM, PetscInt *, char ***, IS **, IS **, DM **))
981: {
982: PetscBool isshell;
984: PetscFunctionBegin;
986: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
987: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
988: dm->ops->createdomaindecomposition = decomp;
989: PetscFunctionReturn(PETSC_SUCCESS);
990: }
992: /*@C
993: DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM
995: Logically Collective on dm
997: Input Parameters:
998: + dm - the shell DM
999: - scatter - the routine to create the scatters
1001: Level: advanced
1003: .seealso: `DMCreateDomainDecompositionScatters()`, `DMShellSetContext()`, `DMShellGetContext()`
1004: @*/
1005: PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **))
1006: {
1007: PetscBool isshell;
1009: PetscFunctionBegin;
1011: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
1012: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
1013: dm->ops->createddscatters = scatter;
1014: PetscFunctionReturn(PETSC_SUCCESS);
1015: }
1017: /*@C
1018: DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM
1020: Logically Collective on dm
1022: Input Parameters:
1023: + dm - the shell DM
1024: - subdm - the routine to create the decomposition
1026: Level: advanced
1028: .seealso: `DMCreateSubDM()`, `DMShellGetCreateSubDM()`, `DMShellSetContext()`, `DMShellGetContext()`
1029: @*/
1030: PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM, PetscInt, const PetscInt[], IS *, DM *))
1031: {
1032: PetscBool isshell;
1034: PetscFunctionBegin;
1036: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
1037: if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
1038: dm->ops->createsubdm = subdm;
1039: PetscFunctionReturn(PETSC_SUCCESS);
1040: }
1042: /*@C
1043: DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM
1045: Logically Collective on dm
1047: Input Parameter:
1048: . dm - the shell DM
1050: Output Parameter:
1051: . subdm - the routine to create the decomposition
1053: Level: advanced
1055: .seealso: `DMCreateSubDM()`, `DMShellSetCreateSubDM()`, `DMShellSetContext()`, `DMShellGetContext()`
1056: @*/
1057: PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM, PetscInt, const PetscInt[], IS *, DM *))
1058: {
1059: PetscBool isshell;
1061: PetscFunctionBegin;
1063: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
1064: PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
1065: *subdm = dm->ops->createsubdm;
1066: PetscFunctionReturn(PETSC_SUCCESS);
1067: }
1069: static PetscErrorCode DMDestroy_Shell(DM dm)
1070: {
1071: DM_Shell *shell = (DM_Shell *)dm->data;
1073: PetscFunctionBegin;
1074: if (shell->destroyctx) PetscCallBack("Destroy Context", (*shell->destroyctx)(shell->ctx));
1075: PetscCall(MatDestroy(&shell->A));
1076: PetscCall(VecDestroy(&shell->Xglobal));
1077: PetscCall(VecDestroy(&shell->Xlocal));
1078: PetscCall(VecScatterDestroy(&shell->gtol));
1079: PetscCall(VecScatterDestroy(&shell->ltog));
1080: PetscCall(VecScatterDestroy(&shell->ltol));
1081: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
1082: PetscCall(PetscFree(shell));
1083: PetscFunctionReturn(PETSC_SUCCESS);
1084: }
1086: static PetscErrorCode DMView_Shell(DM dm, PetscViewer v)
1087: {
1088: DM_Shell *shell = (DM_Shell *)dm->data;
1090: PetscFunctionBegin;
1091: if (shell->Xglobal) PetscCall(VecView(shell->Xglobal, v));
1092: PetscFunctionReturn(PETSC_SUCCESS);
1093: }
1095: static PetscErrorCode DMLoad_Shell(DM dm, PetscViewer v)
1096: {
1097: DM_Shell *shell = (DM_Shell *)dm->data;
1099: PetscFunctionBegin;
1100: PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &shell->Xglobal));
1101: PetscCall(VecLoad(shell->Xglobal, v));
1102: PetscFunctionReturn(PETSC_SUCCESS);
1103: }
1105: PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1106: {
1107: PetscFunctionBegin;
1108: if (subdm) PetscCall(DMShellCreate(PetscObjectComm((PetscObject)dm), subdm));
1109: PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm));
1110: PetscFunctionReturn(PETSC_SUCCESS);
1111: }
1113: PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
1114: {
1115: DM_Shell *shell;
1117: PetscFunctionBegin;
1118: PetscCall(PetscNew(&shell));
1119: dm->data = shell;
1121: dm->ops->destroy = DMDestroy_Shell;
1122: dm->ops->createglobalvector = DMCreateGlobalVector_Shell;
1123: dm->ops->createlocalvector = DMCreateLocalVector_Shell;
1124: dm->ops->creatematrix = DMCreateMatrix_Shell;
1125: dm->ops->view = DMView_Shell;
1126: dm->ops->load = DMLoad_Shell;
1127: dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell;
1128: dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell;
1129: dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell;
1130: dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell;
1131: dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell;
1132: dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell;
1133: dm->ops->createsubdm = DMCreateSubDM_Shell;
1134: PetscCall(DMSetMatType(dm, MATDENSE));
1135: PetscFunctionReturn(PETSC_SUCCESS);
1136: }
1138: /*@
1139: DMShellCreate - Creates a shell DM object, used to manage user-defined problem data
1141: Collective
1143: Input Parameter:
1144: . comm - the processors that will share the global vector
1146: Output Parameters:
1147: . shell - the shell DM
1149: Level: advanced
1151: .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMShellSetContext()`, `DMShellGetContext()`
1152: @*/
1153: PetscErrorCode DMShellCreate(MPI_Comm comm, DM *dm)
1154: {
1155: PetscFunctionBegin;
1157: PetscCall(DMCreate(comm, dm));
1158: PetscCall(DMSetType(*dm, DMSHELL));
1159: PetscCall(DMSetUp(*dm));
1160: PetscFunctionReturn(PETSC_SUCCESS);
1161: }