Actual source code: dmsnes.c

petsc-3.4.5 2014-06-29
  1: #include <petsc-private/snesimpl.h>   /*I "petscsnes.h" I*/
  2: #include <petsc-private/dmimpl.h>     /*I "petscdm.h" I*/

  6: static PetscErrorCode DMSNESDestroy(DMSNES *kdm)
  7: {

 11:   if (!*kdm) return(0);
 13:   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; return(0);}
 14:   if ((*kdm)->ops->destroy) {((*kdm)->ops->destroy)(*kdm);}
 15:   PetscHeaderDestroy(kdm);
 16:   return(0);
 17: }

 21: PetscErrorCode DMSNESLoad(DMSNES kdm,PetscViewer viewer)
 22: {

 26:   PetscViewerBinaryRead(viewer,&kdm->ops->computefunction,1,PETSC_FUNCTION);
 27:   PetscViewerBinaryRead(viewer,&kdm->ops->computejacobian,1,PETSC_FUNCTION);
 28:   return(0);
 29: }

 33: PetscErrorCode DMSNESView(DMSNES kdm,PetscViewer viewer)
 34: {
 36:   PetscBool      isascii,isbinary;

 39:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 40:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
 41:   if (isascii) {
 42: #if defined(PETSC_SERIALIZE_FUNCTIONS)
 43:     const char *fname;

 45:     PetscFPTFind(kdm->ops->computefunction,&fname);
 46:     if (fname) {
 47:       PetscViewerASCIIPrintf(viewer,"Function used by SNES: %s\n",fname);
 48:     }
 49:     PetscFPTFind(kdm->ops->computejacobian,&fname);
 50:     if (fname) {
 51:       PetscViewerASCIIPrintf(viewer,"Jacobian function used by SNES: %s\n",fname);
 52:     }
 53: #endif
 54:   } else if (isbinary) {
 55:     struct {
 56:       PetscErrorCode (*func)(SNES,Vec,Vec,void*);
 57:       PetscErrorCode (*jac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 58:     } funcstruct = {kdm->ops->computefunction,kdm->ops->computejacobian};
 59:     PetscViewerBinaryWrite(viewer,&funcstruct,2,PETSC_FUNCTION,PETSC_FALSE);
 60:   }
 61:   return(0);
 62: }

 66: static PetscErrorCode DMSNESCreate(MPI_Comm comm,DMSNES *kdm)
 67: {

 71: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
 72:   SNESInitializePackage();
 73: #endif
 74:   PetscHeaderCreate(*kdm, _p_DMSNES, struct _DMSNESOps, DMSNES_CLASSID,  "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView);
 75:   PetscMemzero((*kdm)->ops, sizeof(struct _DMSNESOps));
 76:   return(0);
 77: }

 81: /* Attaches the DMSNES to the coarse level.
 82:  * Under what conditions should we copy versus duplicate?
 83:  */
 84: static PetscErrorCode DMCoarsenHook_DMSNES(DM dm,DM dmc,void *ctx)
 85: {

 89:   DMCopyDMSNES(dm,dmc);
 90:   return(0);
 91: }

 95: /* This could restrict auxiliary information to the coarse level.
 96:  */
 97: static PetscErrorCode DMRestrictHook_DMSNES(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
 98: {

101:   return(0);
102: }

106: /* Attaches the DMSNES to the subdomain. */
107: static PetscErrorCode DMSubDomainHook_DMSNES(DM dm,DM subdm,void *ctx)
108: {

112:   DMCopyDMSNES(dm,subdm);
113:   return(0);
114: }

118: /* This could restrict auxiliary information to the coarse level.
119:  */
120: static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
121: {

124:   return(0);
125: }

129: static PetscErrorCode DMRefineHook_DMSNES(DM dm,DM dmf,void *ctx)
130: {

134:   DMCopyDMSNES(dm,dmf);
135:   return(0);
136: }

140: /* This could restrict auxiliary information to the coarse level.
141:  */
142: static PetscErrorCode DMInterpolateHook_DMSNES(DM dm,Mat Interp,DM dmf,void *ctx)
143: {

146:   return(0);
147: }

151: /*@C
152:    DMSNESCopy - copies the information in a DMSNES to another DMSNES

154:    Not Collective

156:    Input Argument:
157: +  kdm - Original DMSNES
158: -  nkdm - DMSNES to receive the data, should have been created with DMSNESCreate()

160:    Level: developer

162: .seealso: DMSNESCreate(), DMSNESDestroy()
163: @*/
164: PetscErrorCode DMSNESCopy(DMSNES kdm,DMSNES nkdm)
165: {

171:   nkdm->ops->computefunction  = kdm->ops->computefunction;
172:   nkdm->ops->computejacobian  = kdm->ops->computejacobian;
173:   nkdm->ops->computegs        = kdm->ops->computegs;
174:   nkdm->ops->computeobjective = kdm->ops->computeobjective;
175:   nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
176:   nkdm->ops->computepfunction = kdm->ops->computepfunction;
177:   nkdm->ops->destroy          = kdm->ops->destroy;
178:   nkdm->ops->duplicate        = kdm->ops->duplicate;

180:   nkdm->functionctx  = kdm->functionctx;
181:   nkdm->gsctx        = kdm->gsctx;
182:   nkdm->pctx         = kdm->pctx;
183:   nkdm->jacobianctx  = kdm->jacobianctx;
184:   nkdm->objectivectx = kdm->objectivectx;
185:   nkdm->data         = kdm->data;

187:   /*
188:   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
189:   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
190:   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
191:   */

193:   /* implementation specific copy hooks */
194:   if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
195:   return(0);
196: }

200: /*@C
201:    DMGetDMSNES - get read-only private DMSNES context from a DM

203:    Not Collective

205:    Input Argument:
206: .  dm - DM to be used with SNES

208:    Output Argument:
209: .  snesdm - private DMSNES context

211:    Level: developer

213:    Notes:
214:    Use DMGetDMSNESWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.

216: .seealso: DMGetDMSNESWrite()
217: @*/
218: PetscErrorCode DMGetDMSNES(DM dm,DMSNES *snesdm)
219: {

224:   *snesdm = (DMSNES) dm->dmsnes;
225:   if (!*snesdm) {
226:     PetscInfo(dm,"Creating new DMSNES\n");
227:     DMSNESCreate(PetscObjectComm((PetscObject)dm),snesdm);

229:     dm->dmsnes = (PetscObject) *snesdm;

231:     DMCoarsenHookAdd(dm,DMCoarsenHook_DMSNES,DMRestrictHook_DMSNES,NULL);
232:     DMRefineHookAdd(dm,DMRefineHook_DMSNES,DMInterpolateHook_DMSNES,NULL);
233:     DMSubDomainHookAdd(dm,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
234:   }
235:   return(0);
236: }

240: /*@C
241:    DMGetDMSNESWrite - get write access to private DMSNES context from a DM

243:    Not Collective

245:    Input Argument:
246: .  dm - DM to be used with SNES

248:    Output Argument:
249: .  snesdm - private DMSNES context

251:    Level: developer

253: .seealso: DMGetDMSNES()
254: @*/
255: PetscErrorCode DMGetDMSNESWrite(DM dm,DMSNES *snesdm)
256: {
258:   DMSNES         sdm;

262:   DMGetDMSNES(dm,&sdm);
263:   if (!sdm->originaldm) sdm->originaldm = dm;
264:   if (sdm->originaldm != dm) {  /* Copy on write */
265:     DMSNES oldsdm = sdm;
266:     PetscInfo(dm,"Copying DMSNES due to write\n");
267:     DMSNESCreate(PetscObjectComm((PetscObject)dm),&sdm);
268:     DMSNESCopy(oldsdm,sdm);
269:     DMSNESDestroy((DMSNES*)&dm->dmsnes);
270:     dm->dmsnes = (PetscObject)sdm;
271:   }
272:   *snesdm = sdm;
273:   return(0);
274: }

278: /*@C
279:    DMCopyDMSNES - copies a DM context to a new DM

281:    Logically Collective

283:    Input Arguments:
284: +  dmsrc - DM to obtain context from
285: -  dmdest - DM to add context to

287:    Level: developer

289:    Note:
290:    The context is copied by reference. This function does not ensure that a context exists.

292: .seealso: DMGetDMSNES(), SNESSetDM()
293: @*/
294: PetscErrorCode DMCopyDMSNES(DM dmsrc,DM dmdest)
295: {

301:   DMSNESDestroy((DMSNES*)&dmdest->dmsnes);

303:   dmdest->dmsnes = dmsrc->dmsnes;

305:   PetscObjectReference(dmdest->dmsnes);
306:   DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);
307:   DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);
308:   DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
309:   return(0);
310: }

314: /*@C
315:    DMSNESSetFunction - set SNES residual evaluation function

317:    Not Collective

319:    Input Arguments:
320: +  dm - DM to be used with SNES
321: .  SNESFunction - residual evaluation function
322: -  ctx - context for residual evaluation

324:    Level: advanced

326:    Note:
327:    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
328:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
329:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

331: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction
332: @*/
333: PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*SNESFunction)(SNES,Vec,Vec,void*),void *ctx)
334: {
336:   DMSNES         sdm;

340:   if (SNESFunction || ctx) {
341:     DMGetDMSNESWrite(dm,&sdm);
342:   }
343:   if (SNESFunction) sdm->ops->computefunction = SNESFunction;
344:   if (ctx) sdm->functionctx = ctx;
345:   return(0);
346: }

350: /*@C
351:    DMSNESGetFunction - get SNES residual evaluation function

353:    Not Collective

355:    Input Argument:
356: .  dm - DM to be used with SNES

358:    Output Arguments:
359: +  SNESFunction - residual evaluation function
360: -  ctx - context for residual evaluation

362:    Level: advanced

364:    Note:
365:    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
366:    associated with the DM.

368: .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
369: @*/
370: PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**SNESFunction)(SNES,Vec,Vec,void*),void **ctx)
371: {
373:   DMSNES         sdm;

377:   DMGetDMSNES(dm,&sdm);
378:   if (SNESFunction) *SNESFunction = sdm->ops->computefunction;
379:   if (ctx) *ctx = sdm->functionctx;
380:   return(0);
381: }

385: /*@C
386:    DMSNESSetObjective - set SNES objective evaluation function

388:    Not Collective

390:    Input Arguments:
391: +  dm - DM to be used with SNES
392: .  SNESObjectiveFunction - residual evaluation function
393: -  ctx - context for residual evaluation

395:    Level: advanced

397: .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
398: @*/
399: PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*SNESObjectiveFunction)(SNES,Vec,PetscReal*,void*),void *ctx)
400: {
402:   DMSNES         sdm;

406:   if (SNESObjectiveFunction || ctx) {
407:     DMGetDMSNESWrite(dm,&sdm);
408:   }
409:   if (SNESObjectiveFunction) sdm->ops->computeobjective = SNESObjectiveFunction;
410:   if (ctx) sdm->objectivectx = ctx;
411:   return(0);
412: }

416: /*@C
417:    DMSNESGetObjective - get SNES objective evaluation function

419:    Not Collective

421:    Input Argument:
422: .  dm - DM to be used with SNES

424:    Output Arguments:
425: +  SNESObjectiveFunction- residual evaluation function
426: -  ctx - context for residual evaluation

428:    Level: advanced

430:    Note:
431:    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
432:    associated with the DM.

434: .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
435: @*/
436: PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**SNESObjectiveFunction)(SNES,Vec,PetscReal*,void*),void **ctx)
437: {
439:   DMSNES         sdm;

443:   DMGetDMSNES(dm,&sdm);
444:   if (SNESObjectiveFunction) *SNESObjectiveFunction = sdm->ops->computeobjective;
445:   if (ctx) *ctx = sdm->objectivectx;
446:   return(0);
447: }

451: /*@C
452:    DMSNESSetGS - set SNES Gauss-Seidel relaxation function

454:    Not Collective

456:    Input Argument:
457: +  dm - DM to be used with SNES
458: .  SNESGSFunction - relaxation function
459: -  ctx - context for residual evaluation

461:    Level: advanced

463:    Note:
464:    SNESSetGS() is normally used, but it calls this function internally because the user context is actually
465:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
466:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

468: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
469: @*/
470: PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*SNESGSFunction)(SNES,Vec,Vec,void*),void *ctx)
471: {
473:   DMSNES         sdm;

477:   if (SNESGSFunction || ctx) {
478:     DMGetDMSNESWrite(dm,&sdm);
479:   }
480:   if (SNESGSFunction) sdm->ops->computegs = SNESGSFunction;
481:   if (ctx) sdm->gsctx = ctx;
482:   return(0);
483: }

487: /*@C
488:    DMSNESGetGS - get SNES Gauss-Seidel relaxation function

490:    Not Collective

492:    Input Argument:
493: .  dm - DM to be used with SNES

495:    Output Arguments:
496: +  SNESGSFunction - relaxation function which performs Gauss-Seidel sweeps
497: -  ctx - context for residual evaluation

499:    Level: advanced

501:    Note:
502:    SNESGetGS() is normally used, but it calls this function internally because the user context is actually
503:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
504:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

506: .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESGSFunction
507: @*/
508: PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**SNESGSFunction)(SNES,Vec,Vec,void*),void **ctx)
509: {
511:   DMSNES         sdm;

515:   DMGetDMSNES(dm,&sdm);
516:   if (SNESGSFunction) *SNESGSFunction = sdm->ops->computegs;
517:   if (ctx) *ctx = sdm->gsctx;
518:   return(0);
519: }

523: /*@C
524:    DMSNESSetJacobian - set SNES Jacobian evaluation function

526:    Not Collective

528:    Input Argument:
529: +  dm - DM to be used with SNES
530: .  SNESJacobianFunction - Jacobian evaluation function
531: -  ctx - context for residual evaluation

533:    Level: advanced

535:    Note:
536:    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
537:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
538:    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.

540: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
541: @*/
542: PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
543: {
545:   DMSNES         sdm;

549:   if (SNESJacobianFunction || ctx) {
550:     DMGetDMSNESWrite(dm,&sdm);
551:   }
552:   if (SNESJacobianFunction) sdm->ops->computejacobian = SNESJacobianFunction;
553:   if (ctx) sdm->jacobianctx = ctx;
554:   return(0);
555: }

559: /*@C
560:    DMSNESGetJacobian - get SNES Jacobian evaluation function

562:    Not Collective

564:    Input Argument:
565: .  dm - DM to be used with SNES

567:    Output Arguments:
568: +  SNESJacobianFunction - Jacobian evaluation function
569: -  ctx - context for residual evaluation

571:    Level: advanced

573:    Note:
574:    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
575:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
576:    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.

578: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
579: @*/
580: PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
581: {
583:   DMSNES         sdm;

587:   DMGetDMSNES(dm,&sdm);
588:   if (SNESJacobianFunction) *SNESJacobianFunction = sdm->ops->computejacobian;
589:   if (ctx) *ctx = sdm->jacobianctx;
590:   return(0);
591: }

595: /*@C
596:    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.

598:    Not Collective

600:    Input Argument:
601: +  dm - DM to be used with SNES
602: .  SNESFunction - RHS evaluation function
603: .  SNESJacobianFunction - Picard matrix evaluation function
604: -  ctx - context for residual evaluation

606:    Level: advanced

608: .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
609: @*/
610: PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*SNESFunction)(SNES,Vec,Vec,void*),PetscErrorCode (*SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
611: {
613:   DMSNES         sdm;

617:   DMGetDMSNES(dm,&sdm);
618:   if (SNESFunction) sdm->ops->computepfunction = SNESFunction;
619:   if (SNESJacobianFunction) sdm->ops->computepjacobian = SNESJacobianFunction;
620:   if (ctx) sdm->pctx = ctx;
621:   return(0);
622: }

626: /*@C
627:    DMSNESGetPicard - get SNES Picard iteration evaluation functions

629:    Not Collective

631:    Input Argument:
632: .  dm - DM to be used with SNES

634:    Output Arguments:
635: +  SNESFunction - Jacobian evaluation function;
636: .  SNESJacobianFunction  - RHS evaluation function;
637: -  ctx - context for residual evaluation

639:    Level: advanced

641: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
642: @*/
643: PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**SNESFunction)(SNES,Vec,Vec,void*),PetscErrorCode (**SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
644: {
646:   DMSNES         sdm;

650:   DMGetDMSNES(dm,&sdm);
651:   if (SNESFunction) *SNESFunction = sdm->ops->computepfunction;
652:   if (SNESJacobianFunction) *SNESJacobianFunction = sdm->ops->computepjacobian;
653:   if (ctx) *ctx = sdm->pctx;
654:   return(0);
655: }