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:   PetscCtx           ctx;
 13:   PetscCtxDestroyFn *destroyctx;
 14: } DM_Shell;

 16: /*@
 17:   DMGlobalToLocalBeginDefaultShell - Uses the GlobalToLocal `VecScatter` context set by the user to begin a global to local scatter

 19:   Collective

 21:   Input Parameters:
 22: + dm   - `DMSHELL`
 23: . g    - global vector
 24: . mode - `InsertMode`
 25: - l    - local vector

 27:   Level: advanced

 29:   Note:
 30:   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.

 32: .seealso: `DM`, `DMSHELL`, `DMGlobalToLocalEndDefaultShell()`
 33: @*/
 34: PetscErrorCode DMGlobalToLocalBeginDefaultShell(DM dm, Vec g, InsertMode mode, Vec l)
 35: {
 36:   DM_Shell *shell = (DM_Shell *)dm->data;

 38:   PetscFunctionBegin;
 39:   PetscCheck(shell->gtol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
 40:   PetscCall(VecScatterBegin(shell->gtol, g, l, mode, SCATTER_FORWARD));
 41:   PetscFunctionReturn(PETSC_SUCCESS);
 42: }

 44: /*@
 45:   DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal `VecScatter` context set by the user to end a global to local scatter
 46:   Collective

 48:   Input Parameters:
 49: + dm   - `DMSHELL`
 50: . g    - global vector
 51: . mode - `InsertMode`
 52: - l    - local vector

 54:   Level: advanced

 56: .seealso: `DM`, `DMSHELL`, `DMGlobalToLocalBeginDefaultShell()`
 57: @*/
 58: PetscErrorCode DMGlobalToLocalEndDefaultShell(DM dm, Vec g, InsertMode mode, Vec l)
 59: {
 60:   DM_Shell *shell = (DM_Shell *)dm->data;

 62:   PetscFunctionBegin;
 63:   PetscCheck(shell->gtol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
 64:   PetscCall(VecScatterEnd(shell->gtol, g, l, mode, SCATTER_FORWARD));
 65:   PetscFunctionReturn(PETSC_SUCCESS);
 66: }

 68: /*@
 69:   DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal `VecScatter` context set by the user to begin a local to global scatter
 70:   Collective

 72:   Input Parameters:
 73: + dm   - `DMSHELL`
 74: . l    - local vector
 75: . mode - `InsertMode`
 76: - g    - global vector

 78:   Level: advanced

 80:   Note:
 81:   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.

 83: .seealso: `DM`, `DMSHELL`, `DMLocalToGlobalEndDefaultShell()`
 84: @*/
 85: PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm, Vec l, InsertMode mode, Vec g)
 86: {
 87:   DM_Shell *shell = (DM_Shell *)dm->data;

 89:   PetscFunctionBegin;
 90:   PetscCheck(shell->ltog, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
 91:   PetscCall(VecScatterBegin(shell->ltog, l, g, mode, SCATTER_FORWARD));
 92:   PetscFunctionReturn(PETSC_SUCCESS);
 93: }

 95: /*@
 96:   DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal `VecScatter` context set by the user to end a local to global scatter
 97:   Collective

 99:   Input Parameters:
100: + dm   - `DMSHELL`
101: . l    - local vector
102: . mode - `InsertMode`
103: - g    - global vector

105:   Level: advanced

107: .seealso: `DM`, `DMSHELL`, `DMLocalToGlobalBeginDefaultShell()`
108: @*/
109: PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm, Vec l, InsertMode mode, Vec g)
110: {
111:   DM_Shell *shell = (DM_Shell *)dm->data;

113:   PetscFunctionBegin;
114:   PetscCheck(shell->ltog, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()");
115:   PetscCall(VecScatterEnd(shell->ltog, l, g, mode, SCATTER_FORWARD));
116:   PetscFunctionReturn(PETSC_SUCCESS);
117: }

119: /*@
120:   DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal `VecScatter` context set by the user to begin a local to local scatter
121:   Collective

123:   Input Parameters:
124: + dm   - `DMSHELL`
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:
134:   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.

136: .seealso: `DM`, `DMSHELL`, `DMLocalToLocalEndDefaultShell()`
137: @*/
138: PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm, Vec g, InsertMode mode, Vec l)
139: {
140:   DM_Shell *shell = (DM_Shell *)dm->data;

142:   PetscFunctionBegin;
143:   PetscCheck(shell->ltol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()");
144:   PetscCall(VecScatterBegin(shell->ltol, g, l, mode, SCATTER_FORWARD));
145:   PetscFunctionReturn(PETSC_SUCCESS);
146: }

148: /*@
149:   DMLocalToLocalEndDefaultShell - Uses the LocalToLocal `VecScatter` context set by the user to end a local to local scatter
150:   Collective

152:   Input Parameters:
153: + dm   - `DMSHELL`
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: `DM`, `DMSHELL`, `DMLocalToLocalBeginDefaultShell()`
163: @*/
164: PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm, Vec g, InsertMode mode, Vec l)
165: {
166:   DM_Shell *shell = (DM_Shell *)dm->data;

168:   PetscFunctionBegin;
169:   PetscCheck(shell->ltol, ((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()");
170:   PetscCall(VecScatterEnd(shell->ltol, g, l, mode, SCATTER_FORWARD));
171:   PetscFunctionReturn(PETSC_SUCCESS);
172: }

174: static PetscErrorCode DMCreateMatrix_Shell(DM dm, Mat *J)
175: {
176:   DM_Shell *shell = (DM_Shell *)dm->data;
177:   Mat       A;

179:   PetscFunctionBegin;
181:   PetscAssertPointer(J, 2);
182:   if (!shell->A) {
183:     PetscInt m, M;

185:     PetscCheck(shell->Xglobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector");
186:     PetscCall(PetscInfo(dm, "Naively creating matrix using global vector distribution without preallocation\n"));
187:     PetscCall(VecGetSize(shell->Xglobal, &M));
188:     PetscCall(VecGetLocalSize(shell->Xglobal, &m));
189:     PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), &shell->A));
190:     PetscCall(MatSetSizes(shell->A, m, m, M, M));
191:     PetscCall(MatSetType(shell->A, dm->mattype));
192:     PetscCall(MatSetUp(shell->A));
193:   }
194:   A = shell->A;
195:   PetscCall(MatDuplicate(A, MAT_SHARE_NONZERO_PATTERN, J));
196:   PetscCall(MatSetDM(*J, dm));
197:   PetscFunctionReturn(PETSC_SUCCESS);
198: }

200: static PetscErrorCode DMCreateGlobalVector_Shell(DM dm, Vec *gvec)
201: {
202:   DM_Shell *shell = (DM_Shell *)dm->data;
203:   Vec       X;

205:   PetscFunctionBegin;
207:   PetscAssertPointer(gvec, 2);
208:   *gvec = NULL;
209:   X     = shell->Xglobal;
210:   PetscCheck(X, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()");
211:   /* Need to create a copy in order to attach the DM to the vector */
212:   PetscCall(VecDuplicate(X, gvec));
213:   PetscCall(VecZeroEntries(*gvec));
214:   PetscCall(VecSetDM(*gvec, dm));
215:   PetscFunctionReturn(PETSC_SUCCESS);
216: }

218: static PetscErrorCode DMCreateLocalVector_Shell(DM dm, Vec *gvec)
219: {
220:   DM_Shell *shell = (DM_Shell *)dm->data;
221:   Vec       X;

223:   PetscFunctionBegin;
225:   PetscAssertPointer(gvec, 2);
226:   *gvec = NULL;
227:   X     = shell->Xlocal;
228:   PetscCheck(X, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()");
229:   /* Need to create a copy in order to attach the DM to the vector */
230:   PetscCall(VecDuplicate(X, gvec));
231:   PetscCall(VecZeroEntries(*gvec));
232:   PetscCall(VecSetDM(*gvec, dm));
233:   PetscFunctionReturn(PETSC_SUCCESS);
234: }

236: /*@C
237:   DMShellSetDestroyContext - set a function that destroys the context provided with `DMShellSetContext()`

239:   Collective

241:   Input Parameters:
242: + dm         - the `DM` to attach the `destroyctx()` function to
243: - destroyctx - the function that destroys the context

245:   Level: advanced

247: .seealso: `DM`, `DMSHELL`, `DMShellSetContext()`, `DMShellGetContext()`
248: @*/
249: PetscErrorCode DMShellSetDestroyContext(DM dm, PetscCtxDestroyFn *destroyctx)
250: {
251:   DM_Shell *shell = (DM_Shell *)dm->data;
252:   PetscBool isshell;

254:   PetscFunctionBegin;
256:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
257:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
258:   shell->destroyctx = destroyctx;
259:   PetscFunctionReturn(PETSC_SUCCESS);
260: }

262: /*@
263:   DMShellSetContext - set some data to be usable by this `DMSHELL`

265:   Collective

267:   Input Parameters:
268: + dm  - `DMSHELL`
269: - ctx - the context

271:   Level: advanced

273: .seealso: `DM`, `DMSHELL`, `DMCreateMatrix()`, `DMShellGetContext()`
274: @*/
275: PetscErrorCode DMShellSetContext(DM dm, PetscCtx ctx)
276: {
277:   DM_Shell *shell = (DM_Shell *)dm->data;
278:   PetscBool isshell;

280:   PetscFunctionBegin;
282:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
283:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
284:   shell->ctx = ctx;
285:   PetscFunctionReturn(PETSC_SUCCESS);
286: }

288: /*@
289:   DMShellGetContext - Returns the user-provided context associated to the `DMSHELL`

291:   Collective

293:   Input Parameter:
294: . dm - `DMSHELL`

296:   Output Parameter:
297: . ctx - the context

299:   Level: advanced

301:   Fortran Notes:
302:   This only works when the context is a Fortran derived type or a `PetscObject`. Declare `ctx` with
303: .vb
304:   type(tUsertype), pointer :: ctx
305: .ve

307: .seealso: `DM`, `DMSHELL`, `DMCreateMatrix()`, `DMShellSetContext()`
308: @*/
309: PetscErrorCode DMShellGetContext(DM dm, PetscCtxRt ctx)
310: {
311:   DM_Shell *shell = (DM_Shell *)dm->data;
312:   PetscBool isshell;

314:   PetscFunctionBegin;
316:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
317:   PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
318:   *(void **)ctx = shell->ctx;
319:   PetscFunctionReturn(PETSC_SUCCESS);
320: }

322: /*@
323:   DMShellSetMatrix - sets a template matrix associated with the `DMSHELL`

325:   Collective

327:   Input Parameters:
328: + dm - `DMSHELL`
329: - J  - template matrix

331:   Level: advanced

333:   Developer Notes:
334:   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.

336: .seealso: `DM`, `DMSHELL`, `DMCreateMatrix()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
337: @*/
338: PetscErrorCode DMShellSetMatrix(DM dm, Mat J)
339: {
340:   DM_Shell *shell = (DM_Shell *)dm->data;
341:   PetscBool isshell;
342:   DM        mdm;

344:   PetscFunctionBegin;
347:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
348:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
349:   if (J == shell->A) PetscFunctionReturn(PETSC_SUCCESS);
350:   PetscCall(MatGetDM(J, &mdm));
351:   PetscCall(PetscObjectReference((PetscObject)J));
352:   PetscCall(MatDestroy(&shell->A));
353:   if (mdm == dm) {
354:     PetscCall(MatDuplicate(J, MAT_SHARE_NONZERO_PATTERN, &shell->A));
355:     PetscCall(MatSetDM(shell->A, NULL));
356:   } else shell->A = J;
357:   PetscFunctionReturn(PETSC_SUCCESS);
358: }

360: /*@C
361:   DMShellSetCreateMatrix - sets the routine to create a matrix associated with the `DMSHELL`

363:   Logically Collective

365:   Input Parameters:
366: + dm   - the `DMSHELL`
367: - func - the function to create a matrix

369:   Level: advanced

371: .seealso: `DM`, `DMSHELL`, `DMCreateMatrix()`, `DMShellSetMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
372: @*/
373: PetscErrorCode DMShellSetCreateMatrix(DM dm, PetscErrorCode (*func)(DM, Mat *))
374: {
375:   PetscFunctionBegin;
377:   dm->ops->creatematrix = func;
378:   PetscFunctionReturn(PETSC_SUCCESS);
379: }

381: /*@
382:   DMShellSetGlobalVector - sets a template global vector associated with the `DMSHELL`

384:   Logically Collective

386:   Input Parameters:
387: + dm - `DMSHELL`
388: - X  - template vector

390:   Level: advanced

392: .seealso: `DM`, `DMSHELL`, `DMCreateGlobalVector()`, `DMShellSetMatrix()`, `DMShellSetCreateGlobalVector()`
393: @*/
394: PetscErrorCode DMShellSetGlobalVector(DM dm, Vec X)
395: {
396:   DM_Shell *shell = (DM_Shell *)dm->data;
397:   PetscBool isshell;
398:   DM        vdm;

400:   PetscFunctionBegin;
403:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
404:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
405:   PetscCall(VecGetDM(X, &vdm));
406:   /*
407:       if the vector proposed as the new base global vector for the DM is a DM vector associated
408:       with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
409:       we get a circular dependency that prevents the DM from being destroy when it should be.
410:       This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
411:       DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
412:       to set its input vector (which is associated with the DM) as the base global vector.
413:       Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
414:       for pointing out the problem.
415:    */
416:   if (vdm == dm) PetscFunctionReturn(PETSC_SUCCESS);
417:   PetscCall(PetscObjectReference((PetscObject)X));
418:   PetscCall(VecDestroy(&shell->Xglobal));
419:   shell->Xglobal = X;
420:   PetscCall(DMClearGlobalVectors(dm));
421:   PetscCall(DMClearNamedGlobalVectors(dm));
422:   PetscFunctionReturn(PETSC_SUCCESS);
423: }

425: /*@
426:   DMShellGetGlobalVector - Returns the template global vector associated with the `DMSHELL`, or `NULL` if it was not set

428:   Not Collective

430:   Input Parameters:
431: + dm - `DMSHELL`
432: - X  - template vector

434:   Level: advanced

436: .seealso: `DM`, `DMSHELL`, `DMShellSetGlobalVector()`, `DMShellSetCreateGlobalVector()`, `DMCreateGlobalVector()`
437: @*/
438: PetscErrorCode DMShellGetGlobalVector(DM dm, Vec *X)
439: {
440:   DM_Shell *shell = (DM_Shell *)dm->data;
441:   PetscBool isshell;

443:   PetscFunctionBegin;
445:   PetscAssertPointer(X, 2);
446:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
447:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
448:   *X = shell->Xglobal;
449:   PetscFunctionReturn(PETSC_SUCCESS);
450: }

452: /*@C
453:   DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the `DMSHELL`

455:   Logically Collective

457:   Input Parameters:
458: + dm   - the `DMSHELL`
459: - func - the creation routine

461:   Level: advanced

463: .seealso: `DM`, `DMSHELL`, `DMShellSetGlobalVector()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
464: @*/
465: PetscErrorCode DMShellSetCreateGlobalVector(DM dm, PetscErrorCode (*func)(DM, Vec *))
466: {
467:   PetscFunctionBegin;
469:   dm->ops->createglobalvector = func;
470:   PetscFunctionReturn(PETSC_SUCCESS);
471: }

473: /*@
474:   DMShellSetLocalVector - sets a template local vector associated with the `DMSHELL`

476:   Logically Collective

478:   Input Parameters:
479: + dm - `DMSHELL`
480: - X  - template vector

482:   Level: advanced

484: .seealso: `DM`, `DMSHELL`, `DMCreateLocalVector()`, `DMShellSetMatrix()`, `DMShellSetCreateLocalVector()`
485: @*/
486: PetscErrorCode DMShellSetLocalVector(DM dm, Vec X)
487: {
488:   DM_Shell *shell = (DM_Shell *)dm->data;
489:   PetscBool isshell;
490:   DM        vdm;

492:   PetscFunctionBegin;
495:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
496:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
497:   PetscCall(VecGetDM(X, &vdm));
498:   /*
499:       if the vector proposed as the new base global vector for the DM is a DM vector associated
500:       with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
501:       we get a circular dependency that prevents the DM from being destroy when it should be.
502:       This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
503:       DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
504:       to set its input vector (which is associated with the DM) as the base global vector.
505:       Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
506:       for pointing out the problem.
507:    */
508:   if (vdm == dm) PetscFunctionReturn(PETSC_SUCCESS);
509:   PetscCall(PetscObjectReference((PetscObject)X));
510:   PetscCall(VecDestroy(&shell->Xlocal));
511:   shell->Xlocal = X;
512:   PetscCall(DMClearLocalVectors(dm));
513:   PetscCall(DMClearNamedLocalVectors(dm));
514:   PetscFunctionReturn(PETSC_SUCCESS);
515: }

517: /*@C
518:   DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the `DMSHELL`

520:   Logically Collective

522:   Input Parameters:
523: + dm   - the `DMSHELL`
524: - func - the creation routine

526:   Level: advanced

528: .seealso: `DM`, `DMSHELL`, `DMShellSetLocalVector()`, `DMShellSetCreateMatrix()`, `DMShellSetContext()`, `DMShellGetContext()`
529: @*/
530: PetscErrorCode DMShellSetCreateLocalVector(DM dm, PetscErrorCode (*func)(DM, Vec *))
531: {
532:   PetscFunctionBegin;
534:   dm->ops->createlocalvector = func;
535:   PetscFunctionReturn(PETSC_SUCCESS);
536: }

538: /*@C
539:   DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter

541:   Logically Collective

543:   Input Parameters:
544: + dm    - the `DMSHELL`
545: . begin - the routine that begins the global to local scatter
546: - end   - the routine that ends the global to local scatter

548:   Level: advanced

550:   Note:
551:   If these functions are not provided but `DMShellSetGlobalToLocalVecScatter()` is called then
552:   `DMGlobalToLocalBeginDefaultShell()`/`DMGlobalToLocalEndDefaultShell()` are used to perform the transfers

554: .seealso: `DM`, `DMSHELL`, `DMShellSetLocalToGlobal()`, `DMGlobalToLocalBeginDefaultShell()`, `DMGlobalToLocalEndDefaultShell()`
555: @*/
556: PetscErrorCode DMShellSetGlobalToLocal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec))
557: {
558:   PetscFunctionBegin;
560:   dm->ops->globaltolocalbegin = begin;
561:   dm->ops->globaltolocalend   = end;
562:   PetscFunctionReturn(PETSC_SUCCESS);
563: }

565: /*@C
566:   DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter

568:   Logically Collective

570:   Input Parameters:
571: + dm    - the `DMSHELL`
572: . begin - the routine that begins the local to global scatter
573: - end   - the routine that ends the local to global scatter

575:   Level: advanced

577:   Note:
578:   If these functions are not provided but `DMShellSetLocalToGlobalVecScatter()` is called then
579:   `DMLocalToGlobalBeginDefaultShell()`/`DMLocalToGlobalEndDefaultShell()` are used to perform the transfers

581: .seealso: `DM`, `DMSHELL`, `DMShellSetGlobalToLocal()`
582: @*/
583: PetscErrorCode DMShellSetLocalToGlobal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec))
584: {
585:   PetscFunctionBegin;
587:   dm->ops->localtoglobalbegin = begin;
588:   dm->ops->localtoglobalend   = end;
589:   PetscFunctionReturn(PETSC_SUCCESS);
590: }

592: /*@C
593:   DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter

595:   Logically Collective

597:   Input Parameters:
598: + dm    - the `DMSHELL`
599: . begin - the routine that begins the local to local scatter
600: - end   - the routine that ends the local to local scatter

602:   Level: advanced

604:   Note:
605:   If these functions are not provided but `DMShellSetLocalToLocalVecScatter()` is called then
606:   `DMLocalToLocalBeginDefaultShell()`/`DMLocalToLocalEndDefaultShell()` are used to perform the transfers

608: .seealso: `DM`, `DMSHELL`, `DMShellSetGlobalToLocal()`, `DMLocalToLocalBeginDefaultShell()`, `DMLocalToLocalEndDefaultShell()`
609: @*/
610: PetscErrorCode DMShellSetLocalToLocal(DM dm, PetscErrorCode (*begin)(DM, Vec, InsertMode, Vec), PetscErrorCode (*end)(DM, Vec, InsertMode, Vec))
611: {
612:   PetscFunctionBegin;
614:   dm->ops->localtolocalbegin = begin;
615:   dm->ops->localtolocalend   = end;
616:   PetscFunctionReturn(PETSC_SUCCESS);
617: }

619: /*@
620:   DMShellSetGlobalToLocalVecScatter - Sets a `VecScatter` context for global to local communication

622:   Logically Collective

624:   Input Parameters:
625: + dm   - the `DMSHELL`
626: - gtol - the global to local `VecScatter` context

628:   Level: advanced

630: .seealso: `DM`, `DMSHELL`, `DMShellSetGlobalToLocal()`, `DMGlobalToLocalBeginDefaultShell()`, `DMGlobalToLocalEndDefaultShell()`
631: @*/
632: PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
633: {
634:   DM_Shell *shell = (DM_Shell *)dm->data;

636:   PetscFunctionBegin;
639:   PetscCall(PetscObjectReference((PetscObject)gtol));
640:   PetscCall(VecScatterDestroy(&shell->gtol));
641:   shell->gtol = gtol;
642:   PetscFunctionReturn(PETSC_SUCCESS);
643: }

645: /*@
646:   DMShellSetLocalToGlobalVecScatter - Sets a` VecScatter` context for local to global communication

648:   Logically Collective

650:   Input Parameters:
651: + dm   - the `DMSHELL`
652: - ltog - the local to global `VecScatter` context

654:   Level: advanced

656: .seealso: `DM`, `DMSHELL`, `DMShellSetLocalToGlobal()`, `DMLocalToGlobalBeginDefaultShell()`, `DMLocalToGlobalEndDefaultShell()`
657: @*/
658: PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
659: {
660:   DM_Shell *shell = (DM_Shell *)dm->data;

662:   PetscFunctionBegin;
665:   PetscCall(PetscObjectReference((PetscObject)ltog));
666:   PetscCall(VecScatterDestroy(&shell->ltog));
667:   shell->ltog = ltog;
668:   PetscFunctionReturn(PETSC_SUCCESS);
669: }

671: /*@
672:   DMShellSetLocalToLocalVecScatter - Sets a `VecScatter` context for local to local communication

674:   Logically Collective

676:   Input Parameters:
677: + dm   - the `DMSHELL`
678: - ltol - the local to local `VecScatter` context

680:   Level: advanced

682: .seealso: `DM`, `DMSHELL`, `DMShellSetLocalToLocal()`, `DMLocalToLocalBeginDefaultShell()`, `DMLocalToLocalEndDefaultShell()`
683: @*/
684: PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
685: {
686:   DM_Shell *shell = (DM_Shell *)dm->data;

688:   PetscFunctionBegin;
691:   PetscCall(PetscObjectReference((PetscObject)ltol));
692:   PetscCall(VecScatterDestroy(&shell->ltol));
693:   shell->ltol = ltol;
694:   PetscFunctionReturn(PETSC_SUCCESS);
695: }

697: /*@C
698:   DMShellSetCoarsen - Set the routine used to coarsen the `DMSHELL`

700:   Logically Collective

702:   Input Parameters:
703: + dm      - the `DMSHELL`
704: - coarsen - the routine that coarsens the `DM`

706:   Level: advanced

708: .seealso: `DM`, `DMSHELL`, `DMShellSetRefine()`, `DMCoarsen()`, `DMShellGetCoarsen()`, `DMShellSetContext()`, `DMShellGetContext()`
709: @*/
710: PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *))
711: {
712:   PetscBool isshell;

714:   PetscFunctionBegin;
716:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
717:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
718:   dm->ops->coarsen = coarsen;
719:   PetscFunctionReturn(PETSC_SUCCESS);
720: }

722: /*@C
723:   DMShellGetCoarsen - Get the routine used to coarsen the `DMSHELL`

725:   Logically Collective

727:   Input Parameter:
728: . dm - the `DMSHELL`

730:   Output Parameter:
731: . coarsen - the routine that coarsens the `DM`

733:   Level: advanced

735: .seealso: `DM`, `DMSHELL`, `DMShellSetCoarsen()`, `DMCoarsen()`, `DMShellSetRefine()`, `DMRefine()`
736: @*/
737: PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM, MPI_Comm, DM *))
738: {
739:   PetscBool isshell;

741:   PetscFunctionBegin;
743:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
744:   PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
745:   *coarsen = dm->ops->coarsen;
746:   PetscFunctionReturn(PETSC_SUCCESS);
747: }

749: /*@C
750:   DMShellSetRefine - Set the routine used to refine the `DMSHELL`

752:   Logically Collective

754:   Input Parameters:
755: + dm     - the `DMSHELL`
756: - refine - the routine that refines the `DM`

758:   Level: advanced

760: .seealso: `DM`, `DMSHELL`, `DMShellSetCoarsen()`, `DMRefine()`, `DMShellGetRefine()`, `DMShellSetContext()`, `DMShellGetContext()`
761: @*/
762: PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM, MPI_Comm, DM *))
763: {
764:   PetscBool isshell;

766:   PetscFunctionBegin;
768:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
769:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
770:   dm->ops->refine = refine;
771:   PetscFunctionReturn(PETSC_SUCCESS);
772: }

774: /*@C
775:   DMShellGetRefine - Get the routine used to refine the `DMSHELL`

777:   Logically Collective

779:   Input Parameter:
780: . dm - the `DMSHELL`

782:   Output Parameter:
783: . refine - the routine that refines the `DM`

785:   Level: advanced

787: .seealso: `DM`, `DMSHELL`, `DMShellSetCoarsen()`, `DMCoarsen()`, `DMShellSetRefine()`, `DMRefine()`
788: @*/
789: PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM, MPI_Comm, DM *))
790: {
791:   PetscBool isshell;

793:   PetscFunctionBegin;
795:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
796:   PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
797:   *refine = dm->ops->refine;
798:   PetscFunctionReturn(PETSC_SUCCESS);
799: }

801: /*@C
802:   DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator

804:   Logically Collective

806:   Input Parameters:
807: + dm     - the `DMSHELL`
808: - interp - the routine to create the interpolation

810:   Level: advanced

812: .seealso: `DM`, `DMSHELL`, `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateInterpolation()`, `DMShellSetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()`
813: @*/
814: PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM, DM, Mat *, Vec *))
815: {
816:   PetscBool isshell;

818:   PetscFunctionBegin;
820:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
821:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
822:   dm->ops->createinterpolation = interp;
823:   PetscFunctionReturn(PETSC_SUCCESS);
824: }

826: /*@C
827:   DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator

829:   Logically Collective

831:   Input Parameter:
832: . dm - the `DMSHELL`

834:   Output Parameter:
835: . interp - the routine to create the interpolation

837:   Level: advanced

839: .seealso: `DM`, `DMSHELL`, `DMShellGetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()`
840: @*/
841: PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM, DM, Mat *, Vec *))
842: {
843:   PetscBool isshell;

845:   PetscFunctionBegin;
847:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
848:   PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
849:   *interp = dm->ops->createinterpolation;
850:   PetscFunctionReturn(PETSC_SUCCESS);
851: }

853: /*@C
854:   DMShellSetCreateRestriction - Set the routine used to create the restriction operator

856:   Logically Collective

858:   Input Parameters:
859: + dm          - the `DMSHELL`
860: - restriction - the routine to create the restriction

862:   Level: advanced

864: .seealso: `DM`, `DMSHELL`, `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellGetCreateRestriction()`, `DMShellSetContext()`, `DMShellGetContext()`
865: @*/
866: PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM, DM, Mat *))
867: {
868:   PetscBool isshell;

870:   PetscFunctionBegin;
872:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
873:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
874:   dm->ops->createrestriction = restriction;
875:   PetscFunctionReturn(PETSC_SUCCESS);
876: }

878: /*@C
879:   DMShellGetCreateRestriction - Get the routine used to create the restriction operator

881:   Logically Collective

883:   Input Parameter:
884: . dm - the `DMSHELL`

886:   Output Parameter:
887: . restriction - the routine to create the restriction

889:   Level: advanced

891: .seealso: `DM`, `DMSHELL`, `DMShellSetCreateInjection()`, `DMCreateInterpolation()`, `DMShellSetContext()`, `DMShellGetContext()`
892: @*/
893: PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM, DM, Mat *))
894: {
895:   PetscBool isshell;

897:   PetscFunctionBegin;
899:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
900:   PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
901:   *restriction = dm->ops->createrestriction;
902:   PetscFunctionReturn(PETSC_SUCCESS);
903: }

905: /*@C
906:   DMShellSetCreateInjection - Set the routine used to create the injection operator

908:   Logically Collective

910:   Input Parameters:
911: + dm     - the `DMSHELL`
912: - inject - the routine to create the injection

914:   Level: advanced

916: .seealso: `DM`, `DMSHELL`, `DMShellSetCreateInterpolation()`, `DMCreateInjection()`, `DMShellGetCreateInjection()`, `DMShellSetContext()`, `DMShellGetContext()`
917: @*/
918: PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM, DM, Mat *))
919: {
920:   PetscBool isshell;

922:   PetscFunctionBegin;
924:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
925:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
926:   dm->ops->createinjection = inject;
927:   PetscFunctionReturn(PETSC_SUCCESS);
928: }

930: /*@C
931:   DMShellGetCreateInjection - Get the routine used to create the injection operator

933:   Logically Collective

935:   Input Parameter:
936: . dm - the `DMSHELL`

938:   Output Parameter:
939: . inject - the routine to create the injection

941:   Level: advanced

943: .seealso: `DM`, `DMSHELL`, `DMShellGetCreateInterpolation()`, `DMCreateInjection()`, `DMShellSetContext()`, `DMShellGetContext()`
944: @*/
945: PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM, DM, Mat *))
946: {
947:   PetscBool isshell;

949:   PetscFunctionBegin;
951:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
952:   PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
953:   *inject = dm->ops->createinjection;
954:   PetscFunctionReturn(PETSC_SUCCESS);
955: }

957: /*@C
958:   DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the `DMSHELL`

960:   Logically Collective

962:   Input Parameters:
963: + dm     - the `DMSHELL`
964: - decomp - the routine to create the decomposition

966:   Level: advanced

968: .seealso: `DM`, `DMSHELL`, `DMCreateFieldDecomposition()`, `DMShellSetContext()`, `DMShellGetContext()`
969: @*/
970: PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM, PetscInt *, char ***, IS **, DM **))
971: {
972:   PetscBool isshell;

974:   PetscFunctionBegin;
976:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
977:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
978:   dm->ops->createfielddecomposition = decomp;
979:   PetscFunctionReturn(PETSC_SUCCESS);
980: }

982: /*@C
983:   DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the `DMSHELL`

985:   Logically Collective

987:   Input Parameters:
988: + dm     - the `DMSHELL`
989: - decomp - the routine to create the decomposition

991:   Level: advanced

993: .seealso: `DM`, `DMSHELL`, `DMCreateDomainDecomposition()`, `DMShellSetContext()`, `DMShellGetContext()`
994: @*/
995: PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM, PetscInt *, char ***, IS **, IS **, DM **))
996: {
997:   PetscBool isshell;

999:   PetscFunctionBegin;
1001:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
1002:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
1003:   dm->ops->createdomaindecomposition = decomp;
1004:   PetscFunctionReturn(PETSC_SUCCESS);
1005: }

1007: /*@C
1008:   DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a `DMSHELL`

1010:   Logically Collective

1012:   Input Parameters:
1013: + dm      - the `DMSHELL`
1014: - scatter - the routine to create the scatters

1016:   Level: advanced

1018: .seealso: `DM`, `DMSHELL`, `DMCreateDomainDecompositionScatters()`, `DMShellSetContext()`, `DMShellGetContext()`
1019: @*/
1020: PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **))
1021: {
1022:   PetscBool isshell;

1024:   PetscFunctionBegin;
1026:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
1027:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
1028:   dm->ops->createddscatters = scatter;
1029:   PetscFunctionReturn(PETSC_SUCCESS);
1030: }

1032: /*@C
1033:   DMShellSetCreateSubDM - Set the routine used to create a sub `DM` from the `DMSHELL`

1035:   Logically Collective

1037:   Input Parameters:
1038: + dm    - the `DMSHELL`
1039: - subdm - the routine to create the decomposition

1041:   Level: advanced

1043: .seealso: `DM`, `DMSHELL`, `DMCreateSubDM()`, `DMShellGetCreateSubDM()`, `DMShellSetContext()`, `DMShellGetContext()`
1044: @*/
1045: PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM, PetscInt, const PetscInt[], IS *, DM *))
1046: {
1047:   PetscBool isshell;

1049:   PetscFunctionBegin;
1051:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
1052:   if (!isshell) PetscFunctionReturn(PETSC_SUCCESS);
1053:   dm->ops->createsubdm = subdm;
1054:   PetscFunctionReturn(PETSC_SUCCESS);
1055: }

1057: /*@C
1058:   DMShellGetCreateSubDM - Get the routine used to create a sub DM from the `DMSHELL`

1060:   Logically Collective

1062:   Input Parameter:
1063: . dm - the `DMSHELL`

1065:   Output Parameter:
1066: . subdm - the routine to create the decomposition

1068:   Level: advanced

1070: .seealso: `DM`, `DMSHELL`, `DMCreateSubDM()`, `DMShellSetCreateSubDM()`, `DMShellSetContext()`, `DMShellGetContext()`
1071: @*/
1072: PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM, PetscInt, const PetscInt[], IS *, DM *))
1073: {
1074:   PetscBool isshell;

1076:   PetscFunctionBegin;
1078:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
1079:   PetscCheck(isshell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only use with DMSHELL type DMs");
1080:   *subdm = dm->ops->createsubdm;
1081:   PetscFunctionReturn(PETSC_SUCCESS);
1082: }

1084: static PetscErrorCode DMDestroy_Shell(DM dm)
1085: {
1086:   DM_Shell *shell = (DM_Shell *)dm->data;

1088:   PetscFunctionBegin;
1089:   if (shell->destroyctx) PetscCallBack("Destroy Context", (*shell->destroyctx)(&shell->ctx));
1090:   PetscCall(MatDestroy(&shell->A));
1091:   PetscCall(VecDestroy(&shell->Xglobal));
1092:   PetscCall(VecDestroy(&shell->Xlocal));
1093:   PetscCall(VecScatterDestroy(&shell->gtol));
1094:   PetscCall(VecScatterDestroy(&shell->ltog));
1095:   PetscCall(VecScatterDestroy(&shell->ltol));
1096:   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
1097:   PetscCall(PetscFree(shell));
1098:   PetscFunctionReturn(PETSC_SUCCESS);
1099: }

1101: static PetscErrorCode DMView_Shell(DM dm, PetscViewer v)
1102: {
1103:   DM_Shell *shell = (DM_Shell *)dm->data;

1105:   PetscFunctionBegin;
1106:   if (shell->Xglobal) PetscCall(VecView(shell->Xglobal, v));
1107:   PetscFunctionReturn(PETSC_SUCCESS);
1108: }

1110: static PetscErrorCode DMLoad_Shell(DM dm, PetscViewer v)
1111: {
1112:   DM_Shell *shell = (DM_Shell *)dm->data;

1114:   PetscFunctionBegin;
1115:   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &shell->Xglobal));
1116:   PetscCall(VecLoad(shell->Xglobal, v));
1117:   PetscFunctionReturn(PETSC_SUCCESS);
1118: }

1120: static PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1121: {
1122:   PetscFunctionBegin;
1123:   if (subdm) PetscCall(DMShellCreate(PetscObjectComm((PetscObject)dm), subdm));
1124:   PetscCall(DMCreateSectionSubDM(dm, numFields, fields, NULL, NULL, is, subdm));
1125:   PetscFunctionReturn(PETSC_SUCCESS);
1126: }

1128: PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
1129: {
1130:   DM_Shell *shell;

1132:   PetscFunctionBegin;
1133:   PetscCall(PetscNew(&shell));
1134:   dm->data = shell;

1136:   dm->ops->destroy            = DMDestroy_Shell;
1137:   dm->ops->createglobalvector = DMCreateGlobalVector_Shell;
1138:   dm->ops->createlocalvector  = DMCreateLocalVector_Shell;
1139:   dm->ops->creatematrix       = DMCreateMatrix_Shell;
1140:   dm->ops->view               = DMView_Shell;
1141:   dm->ops->load               = DMLoad_Shell;
1142:   dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell;
1143:   dm->ops->globaltolocalend   = DMGlobalToLocalEndDefaultShell;
1144:   dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell;
1145:   dm->ops->localtoglobalend   = DMLocalToGlobalEndDefaultShell;
1146:   dm->ops->localtolocalbegin  = DMLocalToLocalBeginDefaultShell;
1147:   dm->ops->localtolocalend    = DMLocalToLocalEndDefaultShell;
1148:   dm->ops->createsubdm        = DMCreateSubDM_Shell;
1149:   PetscCall(DMSetMatType(dm, MATDENSE));
1150:   PetscFunctionReturn(PETSC_SUCCESS);
1151: }

1153: /*@
1154:   DMShellCreate - Creates a `DMSHELL` object, used to manage user-defined problem data

1156:   Collective

1158:   Input Parameter:
1159: . comm - the processors that will share the global vector

1161:   Output Parameter:
1162: . dm - the `DMSHELL`

1164:   Level: advanced

1166: .seealso: `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMShellSetContext()`, `DMShellGetContext()`
1167: @*/
1168: PetscErrorCode DMShellCreate(MPI_Comm comm, DM *dm)
1169: {
1170:   PetscFunctionBegin;
1171:   PetscAssertPointer(dm, 2);
1172:   PetscCall(DMCreate(comm, dm));
1173:   PetscCall(DMSetType(*dm, DMSHELL));
1174:   PetscCall(DMSetUp(*dm));
1175:   PetscFunctionReturn(PETSC_SUCCESS);
1176: }