Actual source code: snespatch.c

petsc-master 2019-06-15
Report Typos and Errors
  1: /*
  2:       Defines a SNES that can consist of a collection of SNESes on patches of the domain
  3: */
  4:  #include <petsc/private/snesimpl.h>
  5: #include <petsc/private/pcpatchimpl.h> /* We need internal access to PCPatch right now, until that part is moved to Plex */
  6:  #include <petscsf.h>

  8: typedef struct {
  9:   PC pc; /* The linear patch preconditioner */
 10: } SNES_Patch;

 12: static PetscErrorCode SNESPatchComputeResidual_Private(SNES snes, Vec x, Vec F, void *ctx)
 13: {
 14:   PC                pc      = (PC) ctx;
 15:   PC_PATCH          *pcpatch = (PC_PATCH *) pc->data;
 16:   PetscInt          pt, size, i;
 17:   const PetscInt    *indices;
 18:   const PetscScalar *X;
 19:   PetscScalar       *XWithAll;
 20:   PetscErrorCode    ierr;


 24:   /* scatter from x to patch->patchStateWithAll[pt] */
 25:   pt = pcpatch->currentPatch;
 26:   ISGetSize(pcpatch->dofMappingWithoutToWithAll[pt], &size);

 28:   ISGetIndices(pcpatch->dofMappingWithoutToWithAll[pt], &indices);
 29:   VecGetArrayRead(x, &X);
 30:   VecGetArray(pcpatch->patchStateWithAll[pt], &XWithAll);

 32:   for (i = 0; i < size; ++i) {
 33:     XWithAll[indices[i]] = X[i];
 34:   }

 36:   VecRestoreArray(pcpatch->patchStateWithAll[pt], &XWithAll);
 37:   VecRestoreArrayRead(x, &X);
 38:   ISRestoreIndices(pcpatch->dofMappingWithoutToWithAll[pt], &indices);

 40:   PCPatchComputeFunction_Internal(pc, pcpatch->patchStateWithAll[pt], F, pt);
 41:   return(0);
 42: }

 44: static PetscErrorCode SNESPatchComputeJacobian_Private(SNES snes, Vec x, Mat J, Mat M, void *ctx)
 45: {
 46:   PC                pc      = (PC) ctx;
 47:   PC_PATCH          *pcpatch = (PC_PATCH *) pc->data;
 48:   PetscInt          pt, size, i;
 49:   const PetscInt    *indices;
 50:   const PetscScalar *X;
 51:   PetscScalar       *XWithAll;
 52:   PetscErrorCode    ierr;

 55:   /* scatter from x to patch->patchStateWithAll[pt] */
 56:   pt = pcpatch->currentPatch;
 57:   ISGetSize(pcpatch->dofMappingWithoutToWithAll[pt], &size);

 59:   ISGetIndices(pcpatch->dofMappingWithoutToWithAll[pt], &indices);
 60:   VecGetArrayRead(x, &X);
 61:   VecGetArray(pcpatch->patchStateWithAll[pt], &XWithAll);

 63:   for (i = 0; i < size; ++i) {
 64:     XWithAll[indices[i]] = X[i];
 65:   }

 67:   VecRestoreArray(pcpatch->patchStateWithAll[pt], &XWithAll);
 68:   VecRestoreArrayRead(x, &X);
 69:   ISRestoreIndices(pcpatch->dofMappingWithoutToWithAll[pt], &indices);

 71:   PCPatchComputeOperator_Internal(pc, pcpatch->patchStateWithAll[pt], M, pcpatch->currentPatch, PETSC_FALSE);
 72:   return(0);
 73: }

 75: static PetscErrorCode PCSetUp_PATCH_Nonlinear(PC pc)
 76: {
 77:   PC_PATCH       *patch = (PC_PATCH *) pc->data;
 78:   const char     *prefix;
 79:   PetscInt       i, pStart, dof;

 83:   if (!pc->setupcalled) {
 84:     PetscMalloc1(patch->npatch, &patch->solver);
 85:     PCGetOptionsPrefix(pc, &prefix);
 86:     PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);
 87:     for (i = 0; i < patch->npatch; ++i) {
 88:       SNES snes;

 90:       SNESCreate(PETSC_COMM_SELF, &snes);
 91:       SNESSetOptionsPrefix(snes, prefix);
 92:       SNESAppendOptionsPrefix(snes, "sub_");
 93:       PetscObjectIncrementTabLevel((PetscObject) snes, (PetscObject) pc, 2);
 94:       PetscLogObjectParent((PetscObject) pc, (PetscObject) snes);
 95:       patch->solver[i] = (PetscObject) snes;
 96:     }

 98:     PetscMalloc1(patch->npatch, &patch->patchResidual);
 99:     PetscMalloc1(patch->npatch, &patch->patchState);
100:     PetscMalloc1(patch->npatch, &patch->patchStateWithAll);
101:     for (i = 0; i < patch->npatch; ++i) {
102:       VecDuplicate(patch->patchRHS[i], &patch->patchResidual[i]);
103:       VecDuplicate(patch->patchUpdate[i], &patch->patchState[i]);

105:       PetscSectionGetDof(patch->gtolCountsWithAll, i+pStart, &dof);
106:       VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchStateWithAll[i]);
107:       VecSetUp(patch->patchStateWithAll[i]);
108:     }
109:     VecDuplicate(patch->localUpdate, &patch->localState);
110:   }
111:   for (i = 0; i < patch->npatch; ++i) {
112:     SNES snes = (SNES) patch->solver[i];

114:     SNESSetFunction(snes, patch->patchResidual[i], SNESPatchComputeResidual_Private, pc);
115:     SNESSetJacobian(snes, patch->mat[i], patch->mat[i], SNESPatchComputeJacobian_Private, pc);
116:   }
117:   if (!pc->setupcalled && patch->optionsSet) for (i = 0; i < patch->npatch; ++i) {SNESSetFromOptions((SNES) patch->solver[i]);}
118:   return(0);
119: }

121: static PetscErrorCode PCApply_PATCH_Nonlinear(PC pc, PetscInt i, Vec patchRHS, Vec patchUpdate)
122: {
123:   PC_PATCH      *patch = (PC_PATCH *) pc->data;
124:   PetscInt       pStart;

128:   patch->currentPatch = i;
129:   PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);

131:   /* Scatter the overlapped global state to our patch state vector */
132:   PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);
133:   PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localState, patch->patchState[i], INSERT_VALUES, SCATTER_FORWARD, SCATTER_INTERIOR);
134:   PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localState, patch->patchStateWithAll[i], INSERT_VALUES, SCATTER_FORWARD, SCATTER_WITHALL);

136:   /* Set initial guess to be current state*/
137:   VecCopy(patch->patchState[i], patchUpdate);
138:   /* Solve for new state */
139:   SNESSolve((SNES) patch->solver[i], patchRHS, patchUpdate);
140:   /* To compute update, subtract off previous state */
141:   VecAXPY(patchUpdate, -1.0, patch->patchState[i]);

143:   PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);
144:   return(0);
145: }

147: static PetscErrorCode PCReset_PATCH_Nonlinear(PC pc)
148: {
149:   PC_PATCH      *patch = (PC_PATCH *) pc->data;
150:   PetscInt       i;

154:   if (patch->solver) {
155:     for (i = 0; i < patch->npatch; ++i) {SNESReset((SNES) patch->solver[i]);}
156:   }

158:   if (patch->patchResidual) {
159:     for (i = 0; i < patch->npatch; ++i) {VecDestroy(&patch->patchResidual[i]);}
160:     PetscFree(patch->patchResidual);
161:   }

163:   if (patch->patchState) {
164:     for (i = 0; i < patch->npatch; ++i) {VecDestroy(&patch->patchState[i]);}
165:     PetscFree(patch->patchState);
166:   }

168:   if (patch->patchStateWithAll) {
169:     for (i = 0; i < patch->npatch; ++i) {VecDestroy(&patch->patchStateWithAll[i]);}
170:     PetscFree(patch->patchStateWithAll);
171:   }

173:   VecDestroy(&patch->localState);
174:   return(0);
175: }

177: static PetscErrorCode PCDestroy_PATCH_Nonlinear(PC pc)
178: {
179:   PC_PATCH      *patch = (PC_PATCH *) pc->data;
180:   PetscInt       i;

184:   if (patch->solver) {
185:     for (i = 0; i < patch->npatch; ++i) {SNESDestroy((SNES *) &patch->solver[i]);}
186:     PetscFree(patch->solver);
187:   }
188:   return(0);
189: }

191: static PetscErrorCode PCUpdateMultiplicative_PATCH_Nonlinear(PC pc, PetscInt i, PetscInt pStart)
192: {
193:   PC_PATCH      *patch = (PC_PATCH *) pc->data;

197:   PCPatch_ScatterLocal_Private(pc, i + pStart, patch->patchUpdate[i], patch->localState, ADD_VALUES, SCATTER_REVERSE, SCATTER_INTERIOR);
198:   return(0);
199: }

201: static PetscErrorCode SNESSetUp_Patch(SNES snes)
202: {
203:   SNES_Patch    *patch = (SNES_Patch *) snes->data;
204:   DM             dm;
205:   Mat            dummy;
206:   Vec            F;
207:   PetscInt       n, N;

211:   SNESGetDM(snes, &dm);
212:   PCSetDM(patch->pc, dm);
213:   SNESGetFunction(snes, &F, NULL, NULL);
214:   VecGetLocalSize(F, &n);
215:   VecGetSize(F, &N);
216:   MatCreateShell(PetscObjectComm((PetscObject) snes), n, n, N, N, (void *) snes, &dummy);
217:   PCSetOperators(patch->pc, dummy, dummy);
218:   MatDestroy(&dummy);
219:   PCSetUp(patch->pc);
220:   /* allocate workspace */
221:   return(0);
222: }

224: static PetscErrorCode SNESReset_Patch(SNES snes)
225: {
226:   SNES_Patch    *patch = (SNES_Patch *) snes->data;

230:   PCReset(patch->pc);
231:   return(0);
232: }

234: static PetscErrorCode SNESDestroy_Patch(SNES snes)
235: {
236:   SNES_Patch    *patch = (SNES_Patch *) snes->data;

240:   SNESReset_Patch(snes);
241:   PCDestroy(&patch->pc);
242:   PetscFree(snes->data);
243:   return(0);
244: }

246: static PetscErrorCode SNESSetFromOptions_Patch(PetscOptionItems *PetscOptionsObject, SNES snes)
247: {
248:   SNES_Patch    *patch = (SNES_Patch *) snes->data;
249:   const char    *prefix;

253:   PetscObjectGetOptionsPrefix((PetscObject)snes, &prefix);
254:   PetscObjectSetOptionsPrefix((PetscObject)patch->pc, prefix);
255:   PCSetFromOptions(patch->pc);
256:   return(0);
257: }

259: static PetscErrorCode SNESView_Patch(SNES snes,PetscViewer viewer)
260: {
261:   SNES_Patch    *patch = (SNES_Patch *) snes->data;
262:   PetscBool      iascii;

266:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
267:   if (iascii) {
268:     PetscViewerASCIIPrintf(viewer,"SNESPATCH\n");
269:   }
270:   PetscViewerASCIIPushTab(viewer);
271:   PCView(patch->pc, viewer);
272:   PetscViewerASCIIPopTab(viewer);
273:   return(0);
274: }

276: static PetscErrorCode SNESSolve_Patch(SNES snes)
277: {
278:   SNES_Patch        *patch = (SNES_Patch *) snes->data;
279:   PC_PATCH          *pcpatch = (PC_PATCH *) patch->pc->data;
280:   SNESLineSearch    ls;
281:   Vec               rhs, update, state, residual;
282:   const PetscScalar *globalState  = NULL;
283:   PetscScalar       *localState   = NULL;
284:   PetscInt          its = 0;
285:   PetscReal         xnorm = 0.0, ynorm = 0.0, fnorm = 0.0;
286:   PetscErrorCode    ierr;

289:   SNESGetSolution(snes, &state);
290:   SNESGetSolutionUpdate(snes, &update);
291:   SNESGetRhs(snes, &rhs);

293:   SNESGetFunction(snes, &residual, NULL, NULL);
294:   SNESGetLineSearch(snes, &ls);

296:   SNESSetConvergedReason(snes, SNES_CONVERGED_ITERATING);
297:   VecSet(update, 0.0);
298:   SNESComputeFunction(snes, state, residual);

300:   VecNorm(state, NORM_2, &xnorm);
301:   VecNorm(residual, NORM_2, &fnorm);
302:   snes->ttol = fnorm*snes->rtol;

304:   if (snes->ops->converged) {
305:     (*snes->ops->converged)(snes,its,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
306:   } else {
307:     SNESConvergedSkip(snes,its,xnorm,ynorm,fnorm,&snes->reason,0);
308:   }
309:   SNESLogConvergenceHistory(snes, fnorm, 0); /* should we count lits from the patches? */
310:   SNESMonitor(snes, its, fnorm);

312:   /* The main solver loop */
313:   for (its = 0; its < snes->max_its; its++) {

315:     SNESSetIterationNumber(snes, its);

317:     /* Scatter state vector to overlapped vector on all patches.
318:        The vector pcpatch->localState is scattered to each patch
319:        in PCApply_PATCH_Nonlinear. */
320:     VecGetArrayRead(state, &globalState);
321:     VecGetArray(pcpatch->localState, &localState);
322:     PetscSFBcastBegin(pcpatch->defaultSF, MPIU_SCALAR, globalState, localState);
323:     PetscSFBcastEnd(pcpatch->defaultSF, MPIU_SCALAR, globalState, localState);
324:     VecRestoreArray(pcpatch->localState, &localState);
325:     VecRestoreArrayRead(state, &globalState);

327:     /* The looping over patches happens here */
328:     PCApply(patch->pc, rhs, update);

330:     /* Apply a line search. This will often be basic with
331:        damping = 1/(max number of patches a dof can be in),
332:        but not always */
333:     VecScale(update, -1.0);
334:     SNESLineSearchApply(ls, state, residual, &fnorm, update);

336:     VecNorm(state, NORM_2, &xnorm);
337:     VecNorm(update, NORM_2, &ynorm);

339:     if (snes->ops->converged) {
340:       (*snes->ops->converged)(snes,its,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
341:     } else {
342:       SNESConvergedSkip(snes,its,xnorm,ynorm,fnorm,&snes->reason,0);
343:     }
344:     SNESLogConvergenceHistory(snes, fnorm, 0); /* FIXME: should we count lits? */
345:     SNESMonitor(snes, its, fnorm);
346:   }

348:   if (its == snes->max_its) { SNESSetConvergedReason(snes, SNES_DIVERGED_MAX_IT); }
349:   return(0);
350: }

352: /*MC
353:   SNESPATCH - Solve a nonlinear problem by composing together many nonlinear solvers on patches

355:   Level: intermediate

357: .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
358:            PCPATCH

360:    References:
361: .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", SIAM Review, 57(4), 2015

363: M*/
364: PETSC_EXTERN PetscErrorCode SNESCreate_Patch(SNES snes)
365: {
367:   SNES_Patch     *patch;
368:   PC_PATCH       *patchpc;
369:   SNESLineSearch linesearch;

372:   PetscNewLog(snes, &patch);

374:   snes->ops->solve          = SNESSolve_Patch;
375:   snes->ops->setup          = SNESSetUp_Patch;
376:   snes->ops->reset          = SNESReset_Patch;
377:   snes->ops->destroy        = SNESDestroy_Patch;
378:   snes->ops->setfromoptions = SNESSetFromOptions_Patch;
379:   snes->ops->view           = SNESView_Patch;

381:   SNESGetLineSearch(snes,&linesearch);
382:   SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
383:   snes->usesksp        = PETSC_FALSE;

385:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

387:   snes->data = (void *) patch;
388:   PCCreate(PetscObjectComm((PetscObject) snes), &patch->pc);
389:   PCSetType(patch->pc, PCPATCH);

391:   patchpc = (PC_PATCH*) patch->pc->data;
392:   patchpc->classname = "snes";
393:   patchpc->isNonlinear = PETSC_TRUE;

395:   patchpc->setupsolver   = PCSetUp_PATCH_Nonlinear;
396:   patchpc->applysolver   = PCApply_PATCH_Nonlinear;
397:   patchpc->resetsolver   = PCReset_PATCH_Nonlinear;
398:   patchpc->destroysolver = PCDestroy_PATCH_Nonlinear;
399:   patchpc->updatemultiplicative = PCUpdateMultiplicative_PATCH_Nonlinear;

401:   return(0);
402: }

404: PetscErrorCode SNESPatchSetDiscretisationInfo(SNES snes, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap,
405:                                             const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
406: {
407:   SNES_Patch     *patch = (SNES_Patch *) snes->data;
409:   DM             dm;

412:   SNESGetDM(snes, &dm);
413:   if (!dm) SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch SNES\n");
414:   PCSetDM(patch->pc, dm);
415:   PCPatchSetDiscretisationInfo(patch->pc, nsubspaces, dms, bs, nodesPerCell, cellNodeMap, subspaceOffsets, numGhostBcs, ghostBcNodes, numGlobalBcs, globalBcNodes);
416:   return(0);
417: }

419: PetscErrorCode SNESPatchSetComputeOperator(SNES snes, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx)
420: {
421:   SNES_Patch    *patch = (SNES_Patch *) snes->data;

425:   PCPatchSetComputeOperator(patch->pc, func, ctx);
426:   return(0);
427: }

429: PetscErrorCode SNESPatchSetComputeFunction(SNES snes, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx)
430: {
431:   SNES_Patch    *patch = (SNES_Patch *) snes->data;

435:   PCPatchSetComputeFunction(patch->pc, func, ctx);
436:   return(0);
437: }

439: PetscErrorCode SNESPatchSetConstructType(SNES snes, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx)
440: {
441:   SNES_Patch    *patch = (SNES_Patch *) snes->data;

445:   PCPatchSetConstructType(patch->pc, ctype, func, ctx);
446:   return(0);
447: }

449: PetscErrorCode SNESPatchSetCellNumbering(SNES snes, PetscSection cellNumbering)
450: {
451:   SNES_Patch    *patch = (SNES_Patch *) snes->data;

455:   PCPatchSetCellNumbering(patch->pc, cellNumbering);
456:   return(0);
457: }