Actual source code: dmshell.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  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:   MatDuplicate(A,MAT_SHARE_NONZERO_PATTERN,J);
198:   MatSetDM(*J,dm);
199:   return(0);
200: }

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

211:   *gvec = NULL;
212:   X     = shell->Xglobal;
213:   if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()");
214:   /* Need to create a copy in order to attach the DM to the vector */
215:   VecDuplicate(X,gvec);
216:   VecZeroEntries(*gvec);
217:   VecSetDM(*gvec,dm);
218:   return(0);
219: }

221: PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec)
222: {
224:   DM_Shell       *shell = (DM_Shell*)dm->data;
225:   Vec            X;

230:   *gvec = NULL;
231:   X     = shell->Xlocal;
232:   if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()");
233:   /* Need to create a copy in order to attach the DM to the vector */
234:   VecDuplicate(X,gvec);
235:   VecZeroEntries(*gvec);
236:   VecSetDM(*gvec,dm);
237:   return(0);
238: }

240: /*@
241:    DMShellSetContext - set some data to be usable by this DM

243:    Collective

245:    Input Arguments:
246: +  dm - shell DM
247: -  ctx - the context

249:    Level: advanced

251: .seealso: DMCreateMatrix(), DMShellGetContext()
252: @*/
253: PetscErrorCode DMShellSetContext(DM dm,void *ctx)
254: {
255:   DM_Shell       *shell = (DM_Shell*)dm->data;
257:   PetscBool      isshell;

261:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
262:   if (!isshell) return(0);
263:   shell->ctx = ctx;
264:   return(0);
265: }

267: /*@
268:    DMShellGetContext - Returns the user-provided context associated to the DM

270:    Collective

272:    Input Argument:
273: .  dm - shell DM

275:    Output Argument:
276: .  ctx - the context

278:    Level: advanced

280: .seealso: DMCreateMatrix(), DMShellSetContext()
281: @*/
282: PetscErrorCode DMShellGetContext(DM dm,void **ctx)
283: {
284:   DM_Shell       *shell = (DM_Shell*)dm->data;
286:   PetscBool      isshell;

290:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
291:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
292:   *ctx = shell->ctx;
293:   return(0);
294: }

296: /*@
297:    DMShellSetMatrix - sets a template matrix associated with the DMShell

299:    Collective

301:    Input Arguments:
302: +  dm - shell DM
303: -  J - template matrix

305:    Level: advanced

307:    Developer Notes:
308:     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.

310: .seealso: DMCreateMatrix(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
311: @*/
312: PetscErrorCode DMShellSetMatrix(DM dm,Mat J)
313: {
314:   DM_Shell       *shell = (DM_Shell*)dm->data;
316:   PetscBool      isshell;
317:   DM             mdm;

322:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
323:   if (!isshell) return(0);
324:   if (J == shell->A) return(0);
325:   MatGetDM(J,&mdm);
326:   PetscObjectReference((PetscObject)J);
327:   MatDestroy(&shell->A);
328:   if (mdm == dm) {
329:     MatDuplicate(J,MAT_SHARE_NONZERO_PATTERN,&shell->A);
330:     MatSetDM(shell->A,NULL);
331:   } else shell->A = J;
332:   return(0);
333: }

335: /*@C
336:    DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM

338:    Logically Collective on dm

340:    Input Arguments:
341: +  dm - the shell DM
342: -  func - the function to create a matrix

344:    Level: advanced

346: .seealso: DMCreateMatrix(), DMShellSetMatrix(), DMShellSetContext(), DMShellGetContext()
347: @*/
348: PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*))
349: {
352:   dm->ops->creatematrix = func;
353:   return(0);
354: }

356: /*@
357:    DMShellSetGlobalVector - sets a template global vector associated with the DMShell

359:    Logically Collective on dm

361:    Input Arguments:
362: +  dm - shell DM
363: -  X - template vector

365:    Level: advanced

367: .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector()
368: @*/
369: PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X)
370: {
371:   DM_Shell       *shell = (DM_Shell*)dm->data;
373:   PetscBool      isshell;
374:   DM             vdm;

379:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
380:   if (!isshell) return(0);
381:   VecGetDM(X,&vdm);
382:   /*
383:       if the vector proposed as the new base global vector for the DM is a DM vector associated
384:       with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
385:       we get a circular dependency that prevents the DM from being destroy when it should be.
386:       This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
387:       DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
388:       to set its input vector (which is associated with the DM) as the base global vector.
389:       Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
390:       for pointing out the problem.
391:    */
392:   if (vdm == dm) return(0);
393:   PetscObjectReference((PetscObject)X);
394:   VecDestroy(&shell->Xglobal);
395:   shell->Xglobal = X;
396:   return(0);
397: }

399: /*@C
400:    DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM

402:    Logically Collective

404:    Input Arguments:
405: +  dm - the shell DM
406: -  func - the creation routine

408:    Level: advanced

410: .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
411: @*/
412: PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
413: {
416:   dm->ops->createglobalvector = func;
417:   return(0);
418: }

420: /*@
421:    DMShellSetLocalVector - sets a template local vector associated with the DMShell

423:    Logically Collective on dm

425:    Input Arguments:
426: +  dm - shell DM
427: -  X - template vector

429:    Level: advanced

431: .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector()
432: @*/
433: PetscErrorCode DMShellSetLocalVector(DM dm,Vec X)
434: {
435:   DM_Shell       *shell = (DM_Shell*)dm->data;
437:   PetscBool      isshell;
438:   DM             vdm;

443:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
444:   if (!isshell) return(0);
445:   VecGetDM(X,&vdm);
446:   /*
447:       if the vector proposed as the new base global vector for the DM is a DM vector associated
448:       with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
449:       we get a circular dependency that prevents the DM from being destroy when it should be.
450:       This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
451:       DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
452:       to set its input vector (which is associated with the DM) as the base global vector.
453:       Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
454:       for pointing out the problem.
455:    */
456:   if (vdm == dm) return(0);
457:   PetscObjectReference((PetscObject)X);
458:   VecDestroy(&shell->Xlocal);
459:   shell->Xlocal = X;
460:   return(0);
461: }

463: /*@C
464:    DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM

466:    Logically Collective

468:    Input Arguments:
469: +  dm - the shell DM
470: -  func - the creation routine

472:    Level: advanced

474: .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
475: @*/
476: PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
477: {
480:   dm->ops->createlocalvector = func;
481:   return(0);
482: }

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

487:    Logically Collective on dm

489:    Input Arguments
490: +  dm - the shell DM
491: .  begin - the routine that begins the global to local scatter
492: -  end - the routine that ends the global to local scatter

494:    Notes:
495:     If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then
496:    DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers

498:    Level: advanced

500: .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
501: @*/
502: PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
505:   dm->ops->globaltolocalbegin = begin;
506:   dm->ops->globaltolocalend = end;
507:   return(0);
508: }

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

513:    Logically Collective on dm

515:    Input Arguments
516: +  dm - the shell DM
517: .  begin - the routine that begins the local to global scatter
518: -  end - the routine that ends the local to global scatter

520:    Notes:
521:     If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then
522:    DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers

524:    Level: advanced

526: .seealso: DMShellSetGlobalToLocal()
527: @*/
528: PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
531:   dm->ops->localtoglobalbegin = begin;
532:   dm->ops->localtoglobalend = end;
533:   return(0);
534: }

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

539:    Logically Collective on dm

541:    Input Arguments
542: +  dm - the shell DM
543: .  begin - the routine that begins the local to local scatter
544: -  end - the routine that ends the local to local scatter

546:    Notes:
547:     If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then
548:    DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers

550:    Level: advanced

552: .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
553: @*/
554: PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) {
557:   dm->ops->localtolocalbegin = begin;
558:   dm->ops->localtolocalend = end;
559:   return(0);
560: }

562: /*@
563:    DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication

565:    Logically Collective on dm

567:    Input Arguments
568: +  dm - the shell DM
569: -  gtol - the global to local VecScatter context

571:    Level: advanced

573: .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
574: @*/
575: PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
576: {
577:   DM_Shell       *shell = (DM_Shell*)dm->data;

583:   PetscObjectReference((PetscObject)gtol);
584:   VecScatterDestroy(&shell->gtol);
585:   shell->gtol = gtol;
586:   return(0);
587: }

589: /*@
590:    DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication

592:    Logically Collective on dm

594:    Input Arguments
595: +  dm - the shell DM
596: -  ltog - the local to global VecScatter context

598:    Level: advanced

600: .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell()
601: @*/
602: PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
603: {
604:   DM_Shell       *shell = (DM_Shell*)dm->data;

610:   PetscObjectReference((PetscObject)ltog);
611:   VecScatterDestroy(&shell->ltog);
612:   shell->ltog = ltog;
613:   return(0);
614: }

616: /*@
617:    DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication

619:    Logically Collective on dm

621:    Input Arguments
622: +  dm - the shell DM
623: -  ltol - the local to local VecScatter context

625:    Level: advanced

627: .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
628: @*/
629: PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
630: {
631:   DM_Shell       *shell = (DM_Shell*)dm->data;

637:   PetscObjectReference((PetscObject)ltol);
638:   VecScatterDestroy(&shell->ltol);
639:   shell->ltol = ltol;
640:   return(0);
641: }

643: /*@C
644:    DMShellSetCoarsen - Set the routine used to coarsen the shell DM

646:    Logically Collective on dm

648:    Input Arguments
649: +  dm - the shell DM
650: -  coarsen - the routine that coarsens the DM

652:    Level: advanced

654: .seealso: DMShellSetRefine(), DMCoarsen(), DMShellGetCoarsen(), DMShellSetContext(), DMShellGetContext()
655: @*/
656: PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*))
657: {
659:   PetscBool      isshell;

663:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
664:   if (!isshell) return(0);
665:   dm->ops->coarsen = coarsen;
666:   return(0);
667: }

669: /*@C
670:    DMShellGetCoarsen - Get the routine used to coarsen the shell DM

672:    Logically Collective on dm

674:    Input Argument:
675: .  dm - the shell DM

677:    Output Argument:
678: .  coarsen - the routine that coarsens the DM

680:    Level: advanced

682: .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
683: @*/
684: PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM,MPI_Comm,DM*))
685: {
687:   PetscBool      isshell;

691:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
692:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
693:   *coarsen = dm->ops->coarsen;
694:   return(0);
695: }

697: /*@C
698:    DMShellSetRefine - Set the routine used to refine the shell DM

700:    Logically Collective on dm

702:    Input Arguments
703: +  dm - the shell DM
704: -  refine - the routine that refines the DM

706:    Level: advanced

708: .seealso: DMShellSetCoarsen(), DMRefine(), DMShellGetRefine(), DMShellSetContext(), DMShellGetContext()
709: @*/
710: PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*))
711: {
713:   PetscBool      isshell;

717:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
718:   if (!isshell) return(0);
719:   dm->ops->refine = refine;
720:   return(0);
721: }

723: /*@C
724:    DMShellGetRefine - Get the routine used to refine the shell DM

726:    Logically Collective on dm

728:    Input Argument:
729: .  dm - the shell DM

731:    Output Argument:
732: .  refine - the routine that refines the DM

734:    Level: advanced

736: .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
737: @*/
738: PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM,MPI_Comm,DM*))
739: {
741:   PetscBool      isshell;

745:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
746:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
747:   *refine = dm->ops->refine;
748:   return(0);
749: }

751: /*@C
752:    DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator

754:    Logically Collective on dm

756:    Input Arguments
757: +  dm - the shell DM
758: -  interp - the routine to create the interpolation

760:    Level: advanced

762: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
763: @*/
764: PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*))
765: {
767:   PetscBool      isshell;

771:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
772:   if (!isshell) return(0);
773:   dm->ops->createinterpolation = interp;
774:   return(0);
775: }

777: /*@C
778:    DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator

780:    Logically Collective on dm

782:    Input Argument:
783: +  dm - the shell DM

785:    Output Argument:
786: -  interp - the routine to create the interpolation

788:    Level: advanced

790: .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
791: @*/
792: PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM,DM,Mat*,Vec*))
793: {
795:   PetscBool      isshell;

799:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
800:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
801:   *interp = dm->ops->createinterpolation;
802:   return(0);
803: }

805: /*@C
806:    DMShellSetCreateRestriction - Set the routine used to create the restriction operator

808:    Logically Collective on dm

810:    Input Arguments
811: +  dm - the shell DM
812: -  striction- the routine to create the restriction

814:    Level: advanced

816: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
817: @*/
818: PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*))
819: {
821:   PetscBool      isshell;

825:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
826:   if (!isshell) return(0);
827:   dm->ops->createrestriction = restriction;
828:   return(0);
829: }

831: /*@C
832:    DMShellGetCreateRestriction - Get the routine used to create the restriction operator

834:    Logically Collective on dm

836:    Input Argument:
837: +  dm - the shell DM

839:    Output Argument:
840: -  restriction - the routine to create the restriction

842:    Level: advanced

844: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext()
845: @*/
846: PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM,DM,Mat*))
847: {
849:   PetscBool      isshell;

853:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
854:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
855:   *restriction = dm->ops->createrestriction;
856:   return(0);
857: }

859: /*@C
860:    DMShellSetCreateInjection - Set the routine used to create the injection operator

862:    Logically Collective on dm

864:    Input Arguments
865: +  dm - the shell DM
866: -  inject - the routine to create the injection

868:    Level: advanced

870: .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellGetCreateInjection(), DMShellSetContext(), DMShellGetContext()
871: @*/
872: PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*))
873: {
875:   PetscBool      isshell;

879:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
880:   if (!isshell) return(0);
881:   dm->ops->createinjection = inject;
882:   return(0);
883: }

885: /*@C
886:    DMShellGetCreateInjection - Get the routine used to create the injection operator

888:    Logically Collective on dm

890:    Input Argument:
891: +  dm - the shell DM

893:    Output Argument:
894: -  inject - the routine to create the injection

896:    Level: advanced

898: .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext()
899: @*/
900: PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM,DM,Mat*))
901: {
903:   PetscBool      isshell;

907:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
908:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
909:   *inject = dm->ops->createinjection;
910:   return(0);
911: }

913: /*@C
914:    DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM

916:    Logically Collective on dm

918:    Input Arguments
919: +  dm - the shell DM
920: -  decomp - the routine to create the decomposition

922:    Level: advanced

924: .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext()
925: @*/
926: PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**))
927: {
929:   PetscBool      isshell;

933:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
934:   if (!isshell) return(0);
935:   dm->ops->createfielddecomposition = decomp;
936:   return(0);
937: }

939: /*@C
940:    DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM

942:    Logically Collective on dm

944:    Input Arguments
945: +  dm - the shell DM
946: -  decomp - the routine to create the decomposition

948:    Level: advanced

950: .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext()
951: @*/
952: PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**))
953: {
955:   PetscBool      isshell;

959:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
960:   if (!isshell) return(0);
961:   dm->ops->createdomaindecomposition = decomp;
962:   return(0);
963: }

965: /*@C
966:    DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM

968:    Logically Collective on dm

970:    Input Arguments
971: +  dm - the shell DM
972: -  scatter - the routine to create the scatters

974:    Level: advanced

976: .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext()
977: @*/
978: PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**))
979: {
981:   PetscBool      isshell;

985:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
986:   if (!isshell) return(0);
987:   dm->ops->createddscatters = scatter;
988:   return(0);
989: }

991: /*@C
992:    DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM

994:    Logically Collective on dm

996:    Input Arguments
997: +  dm - the shell DM
998: -  subdm - the routine to create the decomposition

1000:    Level: advanced

1002: .seealso: DMCreateSubDM(), DMShellGetCreateSubDM(), DMShellSetContext(), DMShellGetContext()
1003: @*/
1004: PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
1005: {
1007:   PetscBool      isshell;

1011:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
1012:   if (!isshell) return(0);
1013:   dm->ops->createsubdm = subdm;
1014:   return(0);
1015: }

1017: /*@C
1018:    DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM

1020:    Logically Collective on dm

1022:    Input Argument:
1023: +  dm - the shell DM

1025:    Output Argument:
1026: -  subdm - the routine to create the decomposition

1028:    Level: advanced

1030: .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext()
1031: @*/
1032: PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
1033: {
1035:   PetscBool      isshell;

1039:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
1040:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
1041:   *subdm = dm->ops->createsubdm;
1042:   return(0);
1043: }

1045: static PetscErrorCode DMDestroy_Shell(DM dm)
1046: {
1048:   DM_Shell       *shell = (DM_Shell*)dm->data;

1051:   MatDestroy(&shell->A);
1052:   VecDestroy(&shell->Xglobal);
1053:   VecDestroy(&shell->Xlocal);
1054:   VecScatterDestroy(&shell->gtol);
1055:   VecScatterDestroy(&shell->ltog);
1056:   VecScatterDestroy(&shell->ltol);
1057:   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
1058:   PetscFree(shell);
1059:   return(0);
1060: }

1062: static PetscErrorCode DMView_Shell(DM dm,PetscViewer v)
1063: {
1065:   DM_Shell       *shell = (DM_Shell*)dm->data;

1068:   VecView(shell->Xglobal,v);
1069:   return(0);
1070: }

1072: static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v)
1073: {
1075:   DM_Shell       *shell = (DM_Shell*)dm->data;

1078:   VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);
1079:   VecLoad(shell->Xglobal,v);
1080:   return(0);
1081: }

1083: PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1084: {

1088:   if (subdm) {DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);}
1089:   DMCreateSectionSubDM(dm, numFields, fields, is, subdm);
1090:   return(0);
1091: }

1093: PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
1094: {
1096:   DM_Shell       *shell;

1099:   PetscNewLog(dm,&shell);
1100:   dm->data = shell;

1102:   dm->ops->destroy            = DMDestroy_Shell;
1103:   dm->ops->createglobalvector = DMCreateGlobalVector_Shell;
1104:   dm->ops->createlocalvector  = DMCreateLocalVector_Shell;
1105:   dm->ops->creatematrix       = DMCreateMatrix_Shell;
1106:   dm->ops->view               = DMView_Shell;
1107:   dm->ops->load               = DMLoad_Shell;
1108:   dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell;
1109:   dm->ops->globaltolocalend   = DMGlobalToLocalEndDefaultShell;
1110:   dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell;
1111:   dm->ops->localtoglobalend   = DMLocalToGlobalEndDefaultShell;
1112:   dm->ops->localtolocalbegin  = DMLocalToLocalBeginDefaultShell;
1113:   dm->ops->localtolocalend    = DMLocalToLocalEndDefaultShell;
1114:   dm->ops->createsubdm        = DMCreateSubDM_Shell;
1115:   DMSetMatType(dm,MATDENSE);
1116:   return(0);
1117: }

1119: /*@
1120:     DMShellCreate - Creates a shell DM object, used to manage user-defined problem data

1122:     Collective

1124:     Input Parameter:
1125: .   comm - the processors that will share the global vector

1127:     Output Parameters:
1128: .   shell - the shell DM

1130:     Level: advanced

1132: .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext()
1133: @*/
1134: PetscErrorCode  DMShellCreate(MPI_Comm comm,DM *dm)
1135: {

1140:   DMCreate(comm,dm);
1141:   DMSetType(*dm,DMSHELL);
1142:   DMSetUp(*dm);
1143:   return(0);
1144: }