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: } 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 Parameters:
 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 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: {
 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 Parameters:
 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 Parameters:
 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 Parameters:
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 Parameters:
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 Parameters:
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 Parameter:
273: .  dm - shell DM

275:    Output Parameter:
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:   *(void**)ctx = shell->ctx;
293:   return(0);
294: }

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

299:    Collective

301:    Input Parameters:
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 Parameters:
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 Parameters:
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: /*@
400:   DMShellGetGlobalVector - Returns the template global vector associated with the DMShell, or NULL if it was not set

402:    Not collective

404:    Input Parameters:
405: +  dm - shell DM
406: -  X - template vector

408:    Level: advanced

410: .seealso: DMShellSetGlobalVector(), DMShellSetCreateGlobalVector(), DMCreateGlobalVector()
411: @*/
412: PetscErrorCode DMShellGetGlobalVector(DM dm, Vec *X)
413: {
414:   DM_Shell      *shell = (DM_Shell *) dm->data;
415:   PetscBool      isshell;

421:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
422:   if (!isshell) return(0);
423:   *X = shell->Xglobal;
424:   return(0);
425: }

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

430:    Logically Collective

432:    Input Parameters:
433: +  dm - the shell DM
434: -  func - the creation routine

436:    Level: advanced

438: .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
439: @*/
440: PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
441: {
444:   dm->ops->createglobalvector = func;
445:   return(0);
446: }

448: /*@
449:    DMShellSetLocalVector - sets a template local vector associated with the DMShell

451:    Logically Collective on dm

453:    Input Parameters:
454: +  dm - shell DM
455: -  X - template vector

457:    Level: advanced

459: .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector()
460: @*/
461: PetscErrorCode DMShellSetLocalVector(DM dm,Vec X)
462: {
463:   DM_Shell       *shell = (DM_Shell*)dm->data;
465:   PetscBool      isshell;
466:   DM             vdm;

471:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
472:   if (!isshell) return(0);
473:   VecGetDM(X,&vdm);
474:   /*
475:       if the vector proposed as the new base global vector for the DM is a DM vector associated
476:       with the same DM then the current base global vector for the DM is ok and if we replace it with the new one
477:       we get a circular dependency that prevents the DM from being destroy when it should be.
478:       This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided
479:       DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries
480:       to set its input vector (which is associated with the DM) as the base global vector.
481:       Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien
482:       for pointing out the problem.
483:    */
484:   if (vdm == dm) return(0);
485:   PetscObjectReference((PetscObject)X);
486:   VecDestroy(&shell->Xlocal);
487:   shell->Xlocal = X;
488:   return(0);
489: }

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

494:    Logically Collective

496:    Input Parameters:
497: +  dm - the shell DM
498: -  func - the creation routine

500:    Level: advanced

502: .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext()
503: @*/
504: PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*))
505: {
508:   dm->ops->createlocalvector = func;
509:   return(0);
510: }

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

515:    Logically Collective on dm

517:    Input Parameters
518: +  dm - the shell DM
519: .  begin - the routine that begins the global to local scatter
520: -  end - the routine that ends the global to local scatter

522:    Notes:
523:     If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then
524:    DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers

526:    Level: advanced

528: .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
529: @*/
530: PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec))
531: {
534:   dm->ops->globaltolocalbegin = begin;
535:   dm->ops->globaltolocalend = end;
536:   return(0);
537: }

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

542:    Logically Collective on dm

544:    Input Parameters
545: +  dm - the shell DM
546: .  begin - the routine that begins the local to global scatter
547: -  end - the routine that ends the local to global scatter

549:    Notes:
550:     If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then
551:    DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers

553:    Level: advanced

555: .seealso: DMShellSetGlobalToLocal()
556: @*/
557: PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec))
558: {
561:   dm->ops->localtoglobalbegin = begin;
562:   dm->ops->localtoglobalend = end;
563:   return(0);
564: }

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

569:    Logically Collective on dm

571:    Input Parameters
572: +  dm - the shell DM
573: .  begin - the routine that begins the local to local scatter
574: -  end - the routine that ends the local to local scatter

576:    Notes:
577:     If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then
578:    DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers

580:    Level: advanced

582: .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
583: @*/
584: PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec))
585: {
588:   dm->ops->localtolocalbegin = begin;
589:   dm->ops->localtolocalend = end;
590:   return(0);
591: }

593: /*@
594:    DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication

596:    Logically Collective on dm

598:    Input Parameters
599: +  dm - the shell DM
600: -  gtol - the global to local VecScatter context

602:    Level: advanced

604: .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell()
605: @*/
606: PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol)
607: {
608:   DM_Shell       *shell = (DM_Shell*)dm->data;

614:   PetscObjectReference((PetscObject)gtol);
615:   VecScatterDestroy(&shell->gtol);
616:   shell->gtol = gtol;
617:   return(0);
618: }

620: /*@
621:    DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication

623:    Logically Collective on dm

625:    Input Parameters
626: +  dm - the shell DM
627: -  ltog - the local to global VecScatter context

629:    Level: advanced

631: .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell()
632: @*/
633: PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog)
634: {
635:   DM_Shell       *shell = (DM_Shell*)dm->data;

641:   PetscObjectReference((PetscObject)ltog);
642:   VecScatterDestroy(&shell->ltog);
643:   shell->ltog = ltog;
644:   return(0);
645: }

647: /*@
648:    DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication

650:    Logically Collective on dm

652:    Input Parameters
653: +  dm - the shell DM
654: -  ltol - the local to local VecScatter context

656:    Level: advanced

658: .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell()
659: @*/
660: PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol)
661: {
662:   DM_Shell       *shell = (DM_Shell*)dm->data;

668:   PetscObjectReference((PetscObject)ltol);
669:   VecScatterDestroy(&shell->ltol);
670:   shell->ltol = ltol;
671:   return(0);
672: }

674: /*@C
675:    DMShellSetCoarsen - Set the routine used to coarsen the shell DM

677:    Logically Collective on dm

679:    Input Parameters
680: +  dm - the shell DM
681: -  coarsen - the routine that coarsens the DM

683:    Level: advanced

685: .seealso: DMShellSetRefine(), DMCoarsen(), DMShellGetCoarsen(), DMShellSetContext(), DMShellGetContext()
686: @*/
687: PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*))
688: {
690:   PetscBool      isshell;

694:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
695:   if (!isshell) return(0);
696:   dm->ops->coarsen = coarsen;
697:   return(0);
698: }

700: /*@C
701:    DMShellGetCoarsen - Get the routine used to coarsen the shell DM

703:    Logically Collective on dm

705:    Input Parameter:
706: .  dm - the shell DM

708:    Output Parameter:
709: .  coarsen - the routine that coarsens the DM

711:    Level: advanced

713: .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
714: @*/
715: PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM,MPI_Comm,DM*))
716: {
718:   PetscBool      isshell;

722:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
723:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
724:   *coarsen = dm->ops->coarsen;
725:   return(0);
726: }

728: /*@C
729:    DMShellSetRefine - Set the routine used to refine the shell DM

731:    Logically Collective on dm

733:    Input Parameters
734: +  dm - the shell DM
735: -  refine - the routine that refines the DM

737:    Level: advanced

739: .seealso: DMShellSetCoarsen(), DMRefine(), DMShellGetRefine(), DMShellSetContext(), DMShellGetContext()
740: @*/
741: PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*))
742: {
744:   PetscBool      isshell;

748:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
749:   if (!isshell) return(0);
750:   dm->ops->refine = refine;
751:   return(0);
752: }

754: /*@C
755:    DMShellGetRefine - Get the routine used to refine the shell DM

757:    Logically Collective on dm

759:    Input Parameter:
760: .  dm - the shell DM

762:    Output Parameter:
763: .  refine - the routine that refines the DM

765:    Level: advanced

767: .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine()
768: @*/
769: PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM,MPI_Comm,DM*))
770: {
772:   PetscBool      isshell;

776:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
777:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
778:   *refine = dm->ops->refine;
779:   return(0);
780: }

782: /*@C
783:    DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator

785:    Logically Collective on dm

787:    Input Parameters
788: +  dm - the shell DM
789: -  interp - the routine to create the interpolation

791:    Level: advanced

793: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
794: @*/
795: PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*))
796: {
798:   PetscBool      isshell;

802:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
803:   if (!isshell) return(0);
804:   dm->ops->createinterpolation = interp;
805:   return(0);
806: }

808: /*@C
809:    DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator

811:    Logically Collective on dm

813:    Input Parameter:
814: .  dm - the shell DM

816:    Output Parameter:
817: .  interp - the routine to create the interpolation

819:    Level: advanced

821: .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
822: @*/
823: PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM,DM,Mat*,Vec*))
824: {
826:   PetscBool      isshell;

830:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
831:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
832:   *interp = dm->ops->createinterpolation;
833:   return(0);
834: }

836: /*@C
837:    DMShellSetCreateRestriction - Set the routine used to create the restriction operator

839:    Logically Collective on dm

841:    Input Parameters
842: +  dm - the shell DM
843: -  striction- the routine to create the restriction

845:    Level: advanced

847: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext()
848: @*/
849: PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*))
850: {
852:   PetscBool      isshell;

856:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
857:   if (!isshell) return(0);
858:   dm->ops->createrestriction = restriction;
859:   return(0);
860: }

862: /*@C
863:    DMShellGetCreateRestriction - Get the routine used to create the restriction operator

865:    Logically Collective on dm

867:    Input Parameter:
868: .  dm - the shell DM

870:    Output Parameter:
871: .  restriction - the routine to create the restriction

873:    Level: advanced

875: .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext()
876: @*/
877: PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM,DM,Mat*))
878: {
880:   PetscBool      isshell;

884:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
885:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
886:   *restriction = dm->ops->createrestriction;
887:   return(0);
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: {
906:   PetscBool      isshell;

910:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
911:   if (!isshell) return(0);
912:   dm->ops->createinjection = inject;
913:   return(0);
914: }

916: /*@C
917:    DMShellGetCreateInjection - Get the routine used to create the injection operator

919:    Logically Collective on dm

921:    Input Parameter:
922: .  dm - the shell DM

924:    Output Parameter:
925: .  inject - the routine to create the injection

927:    Level: advanced

929: .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext()
930: @*/
931: PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM,DM,Mat*))
932: {
934:   PetscBool      isshell;

938:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
939:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
940:   *inject = dm->ops->createinjection;
941:   return(0);
942: }

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

947:    Logically Collective on dm

949:    Input Parameters:
950: +  dm - the shell DM
951: -  decomp - the routine to create the decomposition

953:    Level: advanced

955: .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext()
956: @*/
957: PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**))
958: {
960:   PetscBool      isshell;

964:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
965:   if (!isshell) return(0);
966:   dm->ops->createfielddecomposition = decomp;
967:   return(0);
968: }

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

973:    Logically Collective on dm

975:    Input Parameters:
976: +  dm - the shell DM
977: -  decomp - the routine to create the decomposition

979:    Level: advanced

981: .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext()
982: @*/
983: PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**))
984: {
986:   PetscBool      isshell;

990:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
991:   if (!isshell) return(0);
992:   dm->ops->createdomaindecomposition = decomp;
993:   return(0);
994: }

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

999:    Logically Collective on dm

1001:    Input Parameters:
1002: +  dm - the shell DM
1003: -  scatter - the routine to create the scatters

1005:    Level: advanced

1007: .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext()
1008: @*/
1009: PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**))
1010: {
1012:   PetscBool      isshell;

1016:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
1017:   if (!isshell) return(0);
1018:   dm->ops->createddscatters = scatter;
1019:   return(0);
1020: }

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

1025:    Logically Collective on dm

1027:    Input Parameters:
1028: +  dm - the shell DM
1029: -  subdm - the routine to create the decomposition

1031:    Level: advanced

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

1042:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
1043:   if (!isshell) return(0);
1044:   dm->ops->createsubdm = subdm;
1045:   return(0);
1046: }

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

1051:    Logically Collective on dm

1053:    Input Parameter:
1054: .  dm - the shell DM

1056:    Output Parameter:
1057: .  subdm - the routine to create the decomposition

1059:    Level: advanced

1061: .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext()
1062: @*/
1063: PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM,PetscInt,const PetscInt[],IS*,DM*))
1064: {
1066:   PetscBool      isshell;

1070:   PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);
1071:   if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs");
1072:   *subdm = dm->ops->createsubdm;
1073:   return(0);
1074: }

1076: static PetscErrorCode DMDestroy_Shell(DM dm)
1077: {
1079:   DM_Shell       *shell = (DM_Shell*)dm->data;

1082:   MatDestroy(&shell->A);
1083:   VecDestroy(&shell->Xglobal);
1084:   VecDestroy(&shell->Xlocal);
1085:   VecScatterDestroy(&shell->gtol);
1086:   VecScatterDestroy(&shell->ltog);
1087:   VecScatterDestroy(&shell->ltol);
1088:   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
1089:   PetscFree(shell);
1090:   return(0);
1091: }

1093: static PetscErrorCode DMView_Shell(DM dm,PetscViewer v)
1094: {
1096:   DM_Shell       *shell = (DM_Shell*)dm->data;

1099:   VecView(shell->Xglobal,v);
1100:   return(0);
1101: }

1103: static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v)
1104: {
1106:   DM_Shell       *shell = (DM_Shell*)dm->data;

1109:   VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);
1110:   VecLoad(shell->Xglobal,v);
1111:   return(0);
1112: }

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

1119:   if (subdm) {DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);}
1120:   DMCreateSectionSubDM(dm, numFields, fields, is, subdm);
1121:   return(0);
1122: }

1124: PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm)
1125: {
1127:   DM_Shell       *shell;

1130:   PetscNewLog(dm,&shell);
1131:   dm->data = shell;

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

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

1153:     Collective

1155:     Input Parameter:
1156: .   comm - the processors that will share the global vector

1158:     Output Parameters:
1159: .   shell - the shell DM

1161:     Level: advanced

1163: .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext()
1164: @*/
1165: PetscErrorCode  DMShellCreate(MPI_Comm comm,DM *dm)
1166: {

1171:   DMCreate(comm,dm);
1172:   DMSetType(*dm,DMSHELL);
1173:   DMSetUp(*dm);
1174:   return(0);
1175: }