Actual source code: dmsnes.c

  1: #include <petsc/private/snesimpl.h>
  2: #include <petsc/private/dmimpl.h>

  4: static PetscErrorCode DMSNESDestroy(DMSNES *kdm)
  5: {

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

 17: PetscErrorCode DMSNESLoad(DMSNES kdm,PetscViewer viewer)
 18: {

 22:   PetscViewerBinaryRead(viewer,&kdm->ops->computefunction,1,NULL,PETSC_FUNCTION);
 23:   PetscViewerBinaryRead(viewer,&kdm->ops->computejacobian,1,NULL,PETSC_FUNCTION);
 24:   return(0);
 25: }

 27: PetscErrorCode DMSNESView(DMSNES kdm,PetscViewer viewer)
 28: {
 30:   PetscBool      isascii,isbinary;

 33:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 34:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
 35:   if (isascii) {
 36: #if defined(PETSC_SERIALIZE_FUNCTIONS) && defined(PETSC_SERIALIZE_FUNCTIONS_VIEW)
 37:     const char *fname;

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

 63: static PetscErrorCode DMSNESCreate(MPI_Comm comm,DMSNES *kdm)
 64: {

 68:   SNESInitializePackage();
 69:   PetscHeaderCreate(*kdm, DMSNES_CLASSID,  "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView);
 70:   return(0);
 71: }

 73: /* Attaches the DMSNES to the coarse level.
 74:  * Under what conditions should we copy versus duplicate?
 75:  */
 76: static PetscErrorCode DMCoarsenHook_DMSNES(DM dm,DM dmc,void *ctx)
 77: {

 81:   DMCopyDMSNES(dm,dmc);
 82:   return(0);
 83: }

 85: /* This could restrict auxiliary information to the coarse level.
 86:  */
 87: static PetscErrorCode DMRestrictHook_DMSNES(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
 88: {

 91:   return(0);
 92: }

 94: /* Attaches the DMSNES to the subdomain. */
 95: static PetscErrorCode DMSubDomainHook_DMSNES(DM dm,DM subdm,void *ctx)
 96: {

100:   DMCopyDMSNES(dm,subdm);
101:   return(0);
102: }

104: /* This could restrict auxiliary information to the coarse level.
105:  */
106: static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
107: {

110:   return(0);
111: }

113: static PetscErrorCode DMRefineHook_DMSNES(DM dm,DM dmf,void *ctx)
114: {

118:   DMCopyDMSNES(dm,dmf);
119:   return(0);
120: }

122: /* This could restrict auxiliary information to the coarse level.
123:  */
124: static PetscErrorCode DMInterpolateHook_DMSNES(DM dm,Mat Interp,DM dmf,void *ctx)
125: {

128:   return(0);
129: }

131: /*@C
132:    DMSNESCopy - copies the information in a DMSNES to another DMSNES

134:    Not Collective

136:    Input Parameters:
137: +  kdm - Original DMSNES
138: -  nkdm - DMSNES to receive the data, should have been created with DMSNESCreate()

140:    Level: developer

142: .seealso: DMSNESCreate(), DMSNESDestroy()
143: @*/
144: PetscErrorCode DMSNESCopy(DMSNES kdm,DMSNES nkdm)
145: {

151:   nkdm->ops->computefunction  = kdm->ops->computefunction;
152:   nkdm->ops->computejacobian  = kdm->ops->computejacobian;
153:   nkdm->ops->computegs        = kdm->ops->computegs;
154:   nkdm->ops->computeobjective = kdm->ops->computeobjective;
155:   nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
156:   nkdm->ops->computepfunction = kdm->ops->computepfunction;
157:   nkdm->ops->destroy          = kdm->ops->destroy;
158:   nkdm->ops->duplicate        = kdm->ops->duplicate;

160:   nkdm->functionctx  = kdm->functionctx;
161:   nkdm->gsctx        = kdm->gsctx;
162:   nkdm->pctx         = kdm->pctx;
163:   nkdm->jacobianctx  = kdm->jacobianctx;
164:   nkdm->objectivectx = kdm->objectivectx;
165:   nkdm->originaldm   = kdm->originaldm;

167:   /*
168:   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
169:   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
170:   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
171:   */

173:   /* implementation specific copy hooks */
174:   if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
175:   return(0);
176: }

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

181:    Not Collective

183:    Input Parameter:
184: .  dm - DM to be used with SNES

186:    Output Parameter:
187: .  snesdm - private DMSNES context

189:    Level: developer

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

194: .seealso: DMGetDMSNESWrite()
195: @*/
196: PetscErrorCode DMGetDMSNES(DM dm,DMSNES *snesdm)
197: {

202:   *snesdm = (DMSNES) dm->dmsnes;
203:   if (!*snesdm) {
204:     PetscInfo(dm,"Creating new DMSNES\n");
205:     DMSNESCreate(PetscObjectComm((PetscObject)dm),snesdm);

207:     dm->dmsnes            = (PetscObject) *snesdm;
208:     (*snesdm)->originaldm = dm;
209:     DMCoarsenHookAdd(dm,DMCoarsenHook_DMSNES,DMRestrictHook_DMSNES,NULL);
210:     DMRefineHookAdd(dm,DMRefineHook_DMSNES,DMInterpolateHook_DMSNES,NULL);
211:     DMSubDomainHookAdd(dm,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
212:   }
213:   return(0);
214: }

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

219:    Not Collective

221:    Input Parameter:
222: .  dm - DM to be used with SNES

224:    Output Parameter:
225: .  snesdm - private DMSNES context

227:    Level: developer

229: .seealso: DMGetDMSNES()
230: @*/
231: PetscErrorCode DMGetDMSNESWrite(DM dm,DMSNES *snesdm)
232: {
234:   DMSNES         sdm;

238:   DMGetDMSNES(dm,&sdm);
239:   if (!sdm->originaldm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMSNES has a NULL originaldm");
240:   if (sdm->originaldm != dm) {  /* Copy on write */
241:     DMSNES oldsdm = sdm;
242:     PetscInfo(dm,"Copying DMSNES due to write\n");
243:     DMSNESCreate(PetscObjectComm((PetscObject)dm),&sdm);
244:     DMSNESCopy(oldsdm,sdm);
245:     DMSNESDestroy((DMSNES*)&dm->dmsnes);
246:     dm->dmsnes = (PetscObject)sdm;
247:     sdm->originaldm = dm;
248:   }
249:   *snesdm = sdm;
250:   return(0);
251: }

253: /*@C
254:    DMCopyDMSNES - copies a DM context to a new DM

256:    Logically Collective

258:    Input Parameters:
259: +  dmsrc - DM to obtain context from
260: -  dmdest - DM to add context to

262:    Level: developer

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

267: .seealso: DMGetDMSNES(), SNESSetDM()
268: @*/
269: PetscErrorCode DMCopyDMSNES(DM dmsrc,DM dmdest)
270: {

276:   if (!dmdest->dmsnes) {DMSNESCreate(PetscObjectComm((PetscObject) dmdest), (DMSNES *) &dmdest->dmsnes);}
277:   DMSNESCopy((DMSNES) dmsrc->dmsnes, (DMSNES) dmdest->dmsnes);
278:   DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);
279:   DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);
280:   DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
281:   return(0);
282: }

284: /*@C
285:    DMSNESSetFunction - set SNES residual evaluation function

287:    Not Collective

289:    Input Parameters:
290: +  dm - DM to be used with SNES
291: .  f - residual evaluation function; see SNESFunction for details
292: -  ctx - context for residual evaluation

294:    Level: advanced

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

301: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction
302: @*/
303: PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
304: {
306:   DMSNES         sdm;

310:   if (f || ctx) {
311:     DMGetDMSNESWrite(dm,&sdm);
312:   }
313:   if (f) sdm->ops->computefunction = f;
314:   if (ctx) sdm->functionctx = ctx;
315:   return(0);
316: }

318: /*@C
319:    DMSNESSetMFFunction - set SNES residual evaluation function used in applying the matrix-free Jacobian with -snes_mf_operator

321:    Logically Collective on dm

323:    Input Parameters:
324: +  dm - DM to be used with SNES
325: -  f - residual evaluation function; see SNESFunction for details

327:    Level: advanced

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

338:   if (f || ctx) {
339:     DMGetDMSNESWrite(dm,&sdm);
340:   }
341:   if (f) sdm->ops->computemffunction = f;
342:   if (ctx) sdm->mffunctionctx = ctx;
343:   return(0);
344: }

346: /*@C
347:    DMSNESGetFunction - get SNES residual evaluation function

349:    Not Collective

351:    Input Parameter:
352: .  dm - DM to be used with SNES

354:    Output Parameters:
355: +  f - residual evaluation function; see SNESFunction for details
356: -  ctx - context for residual evaluation

358:    Level: advanced

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

364: .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
365: @*/
366: PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
367: {
369:   DMSNES         sdm;

373:   DMGetDMSNES(dm,&sdm);
374:   if (f) *f = sdm->ops->computefunction;
375:   if (ctx) *ctx = sdm->functionctx;
376:   return(0);
377: }

379: /*@C
380:    DMSNESSetObjective - set SNES objective evaluation function

382:    Not Collective

384:    Input Parameters:
385: +  dm - DM to be used with SNES
386: .  obj - objective evaluation function; see SNESObjectiveFunction for details
387: -  ctx - context for residual evaluation

389:    Level: advanced

391: .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
392: @*/
393: PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
394: {
396:   DMSNES         sdm;

400:   if (obj || ctx) {
401:     DMGetDMSNESWrite(dm,&sdm);
402:   }
403:   if (obj) sdm->ops->computeobjective = obj;
404:   if (ctx) sdm->objectivectx = ctx;
405:   return(0);
406: }

408: /*@C
409:    DMSNESGetObjective - get SNES objective evaluation function

411:    Not Collective

413:    Input Parameter:
414: .  dm - DM to be used with SNES

416:    Output Parameters:
417: +  obj- residual evaluation function; see SNESObjectiveFunction for details
418: -  ctx - context for residual evaluation

420:    Level: advanced

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

426: .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
427: @*/
428: PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
429: {
431:   DMSNES         sdm;

435:   DMGetDMSNES(dm,&sdm);
436:   if (obj) *obj = sdm->ops->computeobjective;
437:   if (ctx) *ctx = sdm->objectivectx;
438:   return(0);
439: }

441: /*@C
442:    DMSNESSetNGS - set SNES Gauss-Seidel relaxation function

444:    Not Collective

446:    Input Parameters:
447: +  dm - DM to be used with SNES
448: .  f  - relaxation function, see SNESGSFunction
449: -  ctx - context for residual evaluation

451:    Level: advanced

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

458: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
459: @*/
460: PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
461: {
463:   DMSNES         sdm;

467:   if (f || ctx) {
468:     DMGetDMSNESWrite(dm,&sdm);
469:   }
470:   if (f) sdm->ops->computegs = f;
471:   if (ctx) sdm->gsctx = ctx;
472:   return(0);
473: }

475: /*@C
476:    DMSNESGetNGS - get SNES Gauss-Seidel relaxation function

478:    Not Collective

480:    Input Parameter:
481: .  dm - DM to be used with SNES

483:    Output Parameters:
484: +  f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
485: -  ctx - context for residual evaluation

487:    Level: advanced

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

494: .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction
495: @*/
496: PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
497: {
499:   DMSNES         sdm;

503:   DMGetDMSNES(dm,&sdm);
504:   if (f) *f = sdm->ops->computegs;
505:   if (ctx) *ctx = sdm->gsctx;
506:   return(0);
507: }

509: /*@C
510:    DMSNESSetJacobian - set SNES Jacobian evaluation function

512:    Not Collective

514:    Input Parameters:
515: +  dm - DM to be used with SNES
516: .  J - Jacobian evaluation function
517: -  ctx - context for residual evaluation

519:    Level: advanced

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

526: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
527: @*/
528: PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
529: {
531:   DMSNES         sdm;

535:   if (J || ctx) {
536:     DMGetDMSNESWrite(dm,&sdm);
537:   }
538:   if (J) sdm->ops->computejacobian = J;
539:   if (ctx) sdm->jacobianctx = ctx;
540:   return(0);
541: }

543: /*@C
544:    DMSNESGetJacobian - get SNES Jacobian evaluation function

546:    Not Collective

548:    Input Parameter:
549: .  dm - DM to be used with SNES

551:    Output Parameters:
552: +  J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
553: -  ctx - context for residual evaluation

555:    Level: advanced

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

562: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
563: @*/
564: PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
565: {
567:   DMSNES         sdm;

571:   DMGetDMSNES(dm,&sdm);
572:   if (J) *J = sdm->ops->computejacobian;
573:   if (ctx) *ctx = sdm->jacobianctx;
574:   return(0);
575: }

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

580:    Not Collective

582:    Input Parameters:
583: +  dm - DM to be used with SNES
584: .  b - RHS evaluation function
585: .  J - Picard matrix evaluation function
586: -  ctx - context for residual evaluation

588:    Level: advanced

590: .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
591: @*/
592: PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
593: {
595:   DMSNES         sdm;

599:   DMGetDMSNES(dm,&sdm);
600:   if (b) sdm->ops->computepfunction = b;
601:   if (J) sdm->ops->computepjacobian = J;
602:   if (ctx) sdm->pctx = ctx;
603:   return(0);
604: }

606: /*@C
607:    DMSNESGetPicard - get SNES Picard iteration evaluation functions

609:    Not Collective

611:    Input Parameter:
612: .  dm - DM to be used with SNES

614:    Output Parameters:
615: +  b - RHS evaluation function; see SNESFunction for details
616: .  J  - RHS evaluation function; see SNESJacobianFunction for detailsa
617: -  ctx - context for residual evaluation

619:    Level: advanced

621: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
622: @*/
623: PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
624: {
626:   DMSNES         sdm;

630:   DMGetDMSNES(dm,&sdm);
631:   if (b) *b = sdm->ops->computepfunction;
632:   if (J) *J = sdm->ops->computepjacobian;
633:   if (ctx) *ctx = sdm->pctx;
634:   return(0);
635: }