Actual source code: fas.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1: /* Defines the basic SNES object */
  2:  #include <../src/snes/impls/fas/fasimpls.h>

  4: const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","SNESFASType","SNES_FAS",0};

  6: static PetscErrorCode SNESReset_FAS(SNES snes)
  7: {
  8:   SNES_FAS       *fas = (SNES_FAS*)snes->data;

 12:   SNESDestroy(&fas->smoothu);
 13:   SNESDestroy(&fas->smoothd);
 14:   MatDestroy(&fas->inject);
 15:   MatDestroy(&fas->interpolate);
 16:   MatDestroy(&fas->restrct);
 17:   VecDestroy(&fas->rscale);
 18:   VecDestroy(&fas->Xg);
 19:   VecDestroy(&fas->Fg);
 20:   if (fas->next) {SNESReset(fas->next);}
 21:   return(0);
 22: }

 24: static PetscErrorCode SNESDestroy_FAS(SNES snes)
 25: {
 26:   SNES_FAS       *fas = (SNES_FAS*)snes->data;

 30:   /* recursively resets and then destroys */
 31:   SNESReset_FAS(snes);
 32:   SNESDestroy(&fas->next);
 33:   PetscFree(fas);
 34:   return(0);
 35: }

 37: static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth)
 38: {
 39:   SNESLineSearch linesearch;
 40:   SNESLineSearch slinesearch;
 41:   void           *lsprectx,*lspostctx;
 42:   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
 43:   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);

 47:   if (!snes->linesearch) return(0);
 48:   SNESGetLineSearch(snes,&linesearch);
 49:   SNESGetLineSearch(smooth,&slinesearch);
 50:   SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
 51:   SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
 52:   SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);
 53:   SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);
 54:   PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);
 55:   return(0);
 56: }

 58: static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth)
 59: {
 60:   SNES_FAS       *fas = (SNES_FAS*) snes->data;

 64:   PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth);
 65:   SNESSetFromOptions(smooth);
 66:   SNESFASSetUpLineSearch_Private(snes, smooth);

 68:   PetscObjectReference((PetscObject)snes->vec_sol);
 69:   PetscObjectReference((PetscObject)snes->vec_sol_update);
 70:   PetscObjectReference((PetscObject)snes->vec_func);
 71:   smooth->vec_sol        = snes->vec_sol;
 72:   smooth->vec_sol_update = snes->vec_sol_update;
 73:   smooth->vec_func       = snes->vec_func;

 75:   if (fas->eventsmoothsetup) {PetscLogEventBegin(fas->eventsmoothsetup,smooth,0,0,0);}
 76:   SNESSetUp(smooth);
 77:   if (fas->eventsmoothsetup) {PetscLogEventEnd(fas->eventsmoothsetup,smooth,0,0,0);}
 78:   return(0);
 79: }

 81: static PetscErrorCode SNESSetUp_FAS(SNES snes)
 82: {
 83:   SNES_FAS       *fas = (SNES_FAS*) snes->data;
 85:   PetscInt       dm_levels;
 86:   SNES           next;
 87:   PetscBool      isFine, hasCreateRestriction, hasCreateInjection;

 90:   SNESFASCycleIsFine(snes, &isFine);
 91:   if (fas->usedmfornumberoflevels && isFine) {
 92:     DMGetRefineLevel(snes->dm,&dm_levels);
 93:     dm_levels++;
 94:     if (dm_levels > fas->levels) {
 95:       /* reset the number of levels */
 96:       SNESFASSetLevels(snes,dm_levels,NULL);
 97:       SNESSetFromOptions(snes);
 98:     }
 99:   }
100:   SNESFASCycleGetCorrection(snes, &next);
101:   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */

103:   SNESSetWorkVecs(snes, 2); /* work vectors used for intergrid transfers */

105:   /* set up the smoothers if they haven't already been set up */
106:   if (!fas->smoothd) {
107:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
108:   }

110:   if (snes->dm) {
111:     /* set the smoother DMs properly */
112:     if (fas->smoothu) {SNESSetDM(fas->smoothu, snes->dm);}
113:     SNESSetDM(fas->smoothd, snes->dm);
114:     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
115:     if (next) {
116:       /* for now -- assume the DM and the evaluation functions have been set externally */
117:       if (!next->dm) {
118:         DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);
119:         SNESSetDM(next, next->dm);
120:       }
121:       /* set the interpolation and restriction from the DM */
122:       if (!fas->interpolate) {
123:         DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);
124:         if (!fas->restrct) {
125:           DMHasCreateRestriction(next->dm, &hasCreateRestriction);
126:           /* DM can create restrictions, use that */
127:           if (hasCreateRestriction) {
128:             DMCreateRestriction(next->dm, snes->dm, &fas->restrct);
129:           } else {
130:             PetscObjectReference((PetscObject)fas->interpolate);
131:             fas->restrct = fas->interpolate;
132:           }
133:         }
134:       }
135:       /* set the injection from the DM */
136:       if (!fas->inject) {
137:         DMHasCreateInjection(next->dm, &hasCreateInjection);
138:         if (hasCreateInjection) {
139:           DMCreateInjection(next->dm, snes->dm, &fas->inject);
140:         }
141:       }
142:     }
143:   }

145:   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
146:   if (fas->galerkin) {
147:     if (next) {
148:       SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next);
149:     }
150:     if (fas->smoothd && fas->level != fas->levels - 1) {
151:       SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes);
152:     }
153:     if (fas->smoothu && fas->level != fas->levels - 1) {
154:       SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes);
155:     }
156:   }

158:   /* sets the down (pre) smoother's default norm and sets it from options */
159:   if (fas->smoothd) {
160:     if (fas->level == 0 && fas->levels != 1) {
161:       SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);
162:     } else {
163:       SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);
164:     }
165:     SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd);
166:   }

168:   /* sets the up (post) smoother's default norm and sets it from options */
169:   if (fas->smoothu) {
170:     if (fas->level != fas->levels - 1) {
171:       SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);
172:     } else {
173:       SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);
174:     }
175:     SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu);
176:   }

178:   if (next) {
179:     /* gotta set up the solution vector for this to work */
180:     if (!next->vec_sol) {SNESFASCreateCoarseVec(snes,&next->vec_sol);}
181:     if (!next->vec_rhs) {SNESFASCreateCoarseVec(snes,&next->vec_rhs);}
182:     PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);
183:     SNESFASSetUpLineSearch_Private(snes, next);
184:     SNESSetUp(next);
185:   }

187:   /* setup FAS work vectors */
188:   if (fas->galerkin) {
189:     VecDuplicate(snes->vec_sol, &fas->Xg);
190:     VecDuplicate(snes->vec_sol, &fas->Fg);
191:   }
192:   return(0);
193: }

195: static PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes)
196: {
197:   SNES_FAS       *fas   = (SNES_FAS*) snes->data;
198:   PetscInt       levels = 1;
199:   PetscBool      flg    = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE;
201:   SNESFASType    fastype;
202:   const char     *optionsprefix;
203:   SNESLineSearch linesearch;
204:   PetscInt       m, n_up, n_down;
205:   SNES           next;
206:   PetscBool      isFine;

209:   SNESFASCycleIsFine(snes, &isFine);
210:   PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");

212:   /* number of levels -- only process most options on the finest level */
213:   if (isFine) {
214:     PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);
215:     if (!flg && snes->dm) {
216:       DMGetRefineLevel(snes->dm,&levels);
217:       levels++;
218:       fas->usedmfornumberoflevels = PETSC_TRUE;
219:     }
220:     SNESFASSetLevels(snes, levels, NULL);
221:     fastype = fas->fastype;
222:     PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);
223:     if (flg) {
224:       SNESFASSetType(snes, fastype);
225:     }

227:     SNESGetOptionsPrefix(snes, &optionsprefix);
228:     PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);
229:     if (flg) {
230:       SNESFASSetCycles(snes, m);
231:     }
232:     PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);
233:     if (flg) {
234:       SNESFASSetContinuation(snes,continuationflg);
235:     }

237:     PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);
238:     if (flg) {
239:       SNESFASSetGalerkin(snes, galerkinflg);
240:     }

242:     if (fas->fastype == SNES_FAS_FULL) {
243:       PetscOptionsBool("-snes_fas_full_downsweep","Smooth on the initial down sweep for full FAS cycles","SNESFASFullSetDownSweep",fas->full_downsweep,&fas->full_downsweep,&flg);
244:       if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);}
245:     }

247:     PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);

249:     PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);

251:     {
252:       PetscViewer viewer;
253:       PetscViewerFormat format;
254:       PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_fas_monitor",&viewer,&format,&monflg);
255:       if (monflg) {
256:         PetscViewerAndFormat *vf;
257:         PetscViewerAndFormatCreate(viewer,format,&vf);
258:         PetscObjectDereference((PetscObject)viewer);
259:         SNESFASSetMonitor(snes,vf,PETSC_TRUE);
260:       }
261:     }
262:     flg    = PETSC_FALSE;
263:     monflg = PETSC_TRUE;
264:     PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);
265:     if (flg) {SNESFASSetLog(snes,monflg);}
266:   }

268:   PetscOptionsTail();

270:   /* setup from the determined types if there is no pointwise procedure or smoother defined */
271:   if (upflg) {
272:     SNESFASSetNumberSmoothUp(snes,n_up);
273:   }
274:   if (downflg) {
275:     SNESFASSetNumberSmoothDown(snes,n_down);
276:   }

278:   /* set up the default line search for coarse grid corrections */
279:   if (fas->fastype == SNES_FAS_ADDITIVE) {
280:     if (!snes->linesearch) {
281:       SNESGetLineSearch(snes, &linesearch);
282:       SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
283:     }
284:   }

286:   /* recursive option setting for the smoothers */
287:   SNESFASCycleGetCorrection(snes, &next);
288:   if (next) {SNESSetFromOptions(next);}
289:   return(0);
290: }

292:  #include <petscdraw.h>
293: static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
294: {
295:   SNES_FAS       *fas = (SNES_FAS*) snes->data;
296:   PetscBool      isFine,iascii,isdraw;
297:   PetscInt       i;
299:   SNES           smoothu, smoothd, levelsnes;

302:   SNESFASCycleIsFine(snes, &isFine);
303:   if (isFine) {
304:     PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
305:     PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
306:     if (iascii) {
307:       PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);
308:       if (fas->galerkin) {
309:         PetscViewerASCIIPrintf(viewer,"  Using Galerkin computed coarse grid function evaluation\n");
310:       } else {
311:         PetscViewerASCIIPrintf(viewer,"  Not using Galerkin computed coarse grid function evaluation\n");
312:       }
313:       for (i=0; i<fas->levels; i++) {
314:         SNESFASGetCycleSNES(snes, i, &levelsnes);
315:         SNESFASCycleGetSmootherUp(levelsnes, &smoothu);
316:         SNESFASCycleGetSmootherDown(levelsnes, &smoothd);
317:         if (!i) {
318:           PetscViewerASCIIPrintf(viewer,"  Coarse grid solver -- level %D -------------------------------\n",i);
319:         } else {
320:           PetscViewerASCIIPrintf(viewer,"  Down solver (pre-smoother) on level %D -------------------------------\n",i);
321:         }
322:         PetscViewerASCIIPushTab(viewer);
323:         if (smoothd) {
324:           SNESView(smoothd,viewer);
325:         } else {
326:           PetscViewerASCIIPrintf(viewer,"Not yet available\n");
327:         }
328:         PetscViewerASCIIPopTab(viewer);
329:         if (i && (smoothd == smoothu)) {
330:           PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) same as down solver (pre-smoother)\n");
331:         } else if (i) {
332:           PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) on level %D -------------------------------\n",i);
333:           PetscViewerASCIIPushTab(viewer);
334:           if (smoothu) {
335:             SNESView(smoothu,viewer);
336:           } else {
337:             PetscViewerASCIIPrintf(viewer,"Not yet available\n");
338:           }
339:           PetscViewerASCIIPopTab(viewer);
340:         }
341:       }
342:     } else if (isdraw) {
343:       PetscDraw draw;
344:       PetscReal x,w,y,bottom,th,wth;
345:       SNES_FAS  *curfas = fas;
346:       PetscViewerDrawGetDraw(viewer,0,&draw);
347:       PetscDrawGetCurrentPoint(draw,&x,&y);
348:       PetscDrawStringGetSize(draw,&wth,&th);
349:       bottom = y - th;
350:       while (curfas) {
351:         if (!curfas->smoothu) {
352:           PetscDrawPushCurrentPoint(draw,x,bottom);
353:           if (curfas->smoothd) SNESView(curfas->smoothd,viewer);
354:           PetscDrawPopCurrentPoint(draw);
355:         } else {
356:           w    = 0.5*PetscMin(1.0-x,x);
357:           PetscDrawPushCurrentPoint(draw,x-w,bottom);
358:           if (curfas->smoothd) SNESView(curfas->smoothd,viewer);
359:           PetscDrawPopCurrentPoint(draw);
360:           PetscDrawPushCurrentPoint(draw,x+w,bottom);
361:           if (curfas->smoothu) SNESView(curfas->smoothu,viewer);
362:           PetscDrawPopCurrentPoint(draw);
363:         }
364:         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
365:         bottom -= 5*th;
366:         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
367:         else curfas = NULL;
368:       }
369:     }
370:   }
371:   return(0);
372: }

374: /*
375: Defines the action of the downsmoother
376:  */
377: static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
378: {
379:   PetscErrorCode      ierr;
380:   SNESConvergedReason reason;
381:   Vec                 FPC;
382:   SNES                smoothd;
383:   PetscBool           flg;
384:   SNES_FAS            *fas = (SNES_FAS*) snes->data;

387:   SNESFASCycleGetSmootherDown(snes, &smoothd);
388:   SNESSetInitialFunction(smoothd, F);
389:   if (fas->eventsmoothsolve) {PetscLogEventBegin(fas->eventsmoothsolve,smoothd,B,X,0);}
390:   SNESSolve(smoothd, B, X);
391:   if (fas->eventsmoothsolve) {PetscLogEventEnd(fas->eventsmoothsolve,smoothd,B,X,0);}
392:   /* check convergence reason for the smoother */
393:   SNESGetConvergedReason(smoothd,&reason);
394:   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
395:     snes->reason = SNES_DIVERGED_INNER;
396:     return(0);
397:   }

399:   SNESGetFunction(smoothd, &FPC, NULL, NULL);
400:   SNESGetAlwaysComputesFinalResidual(smoothd, &flg);
401:   if (!flg) {
402:     SNESComputeFunction(smoothd, X, FPC);
403:   }
404:   VecCopy(FPC, F);
405:   if (fnorm) {VecNorm(F,NORM_2,fnorm);}
406:   return(0);
407: }


410: /*
411: Defines the action of the upsmoother
412:  */
413: static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
414: {
415:   PetscErrorCode      ierr;
416:   SNESConvergedReason reason;
417:   Vec                 FPC;
418:   SNES                smoothu;
419:   PetscBool           flg;
420:   SNES_FAS            *fas = (SNES_FAS*) snes->data;

423:   SNESFASCycleGetSmootherUp(snes, &smoothu);
424:   if (fas->eventsmoothsolve) {PetscLogEventBegin(fas->eventsmoothsolve,smoothu,0,0,0);}
425:   SNESSolve(smoothu, B, X);
426:   if (fas->eventsmoothsolve) {PetscLogEventEnd(fas->eventsmoothsolve,smoothu,0,0,0);}
427:   /* check convergence reason for the smoother */
428:   SNESGetConvergedReason(smoothu,&reason);
429:   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
430:     snes->reason = SNES_DIVERGED_INNER;
431:     return(0);
432:   }
433:   SNESGetFunction(smoothu, &FPC, NULL, NULL);
434:   SNESGetAlwaysComputesFinalResidual(smoothu, &flg);
435:   if (!flg) {
436:     SNESComputeFunction(smoothu, X, FPC);
437:   }
438:   VecCopy(FPC, F);
439:   if (fnorm) {VecNorm(F,NORM_2,fnorm);}
440:   return(0);
441: }

443: /*@
444:    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level

446:    Collective

448:    Input Arguments:
449: .  snes - SNESFAS

451:    Output Arguments:
452: .  Xcoarse - vector on level one coarser than snes

454:    Level: developer

456: .seealso: SNESFASSetRestriction(), SNESFASRestrict()
457: @*/
458: PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
459: {
461:   SNES_FAS       *fas;

466:   fas = (SNES_FAS*)snes->data;
467:   if (fas->rscale) {
468:     VecDuplicate(fas->rscale,Xcoarse);
469:   } else if (fas->inject) {
470:     MatCreateVecs(fas->inject,Xcoarse,NULL);
471:   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");
472:   return(0);
473: }

475: /*@
476:    SNESFASRestrict - restrict a Vec to the next coarser level

478:    Collective

480:    Input Arguments:
481: +  fine - SNES from which to restrict
482: -  Xfine - vector to restrict

484:    Output Arguments:
485: .  Xcoarse - result of restriction

487:    Level: developer

489: .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
490: @*/
491: PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
492: {
494:   SNES_FAS       *fas;

500:   fas = (SNES_FAS*)fine->data;
501:   if (fas->inject) {
502:     MatRestrict(fas->inject,Xfine,Xcoarse);
503:   } else {
504:     MatRestrict(fas->restrct,Xfine,Xcoarse);
505:     VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);
506:   }
507:   return(0);
508: }

510: /*

512: Performs the FAS coarse correction as:

514: fine problem:   F(x) = b
515: coarse problem: F^c(x^c) = b^c

517: b^c = F^c(Rx) - R(F(x) - b)

519:  */
520: PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
521: {
522:   PetscErrorCode      ierr;
523:   Vec                 X_c, Xo_c, F_c, B_c;
524:   SNESConvergedReason reason;
525:   SNES                next;
526:   Mat                 restrct, interpolate;
527:   SNES_FAS            *fasc;

530:   SNESFASCycleGetCorrection(snes, &next);
531:   if (next) {
532:     fasc = (SNES_FAS*)next->data;

534:     SNESFASCycleGetRestriction(snes, &restrct);
535:     SNESFASCycleGetInterpolation(snes, &interpolate);

537:     X_c  = next->vec_sol;
538:     Xo_c = next->work[0];
539:     F_c  = next->vec_func;
540:     B_c  = next->vec_rhs;

542:     if (fasc->eventinterprestrict) {PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);}
543:     SNESFASRestrict(snes, X, Xo_c);
544:     /* restrict the defect: R(F(x) - b) */
545:     MatRestrict(restrct, F, B_c);
546:     if (fasc->eventinterprestrict) {PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);}

548:     if (fasc->eventresidual) {PetscLogEventBegin(fasc->eventresidual,next,0,0,0);}
549:     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
550:     SNESComputeFunction(next, Xo_c, F_c);
551:     if (fasc->eventresidual) {PetscLogEventEnd(fasc->eventresidual,next,0,0,0);}

553:     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
554:     VecCopy(B_c, X_c);
555:     VecCopy(F_c, B_c);
556:     VecCopy(X_c, F_c);
557:     /* set initial guess of the coarse problem to the projected fine solution */
558:     VecCopy(Xo_c, X_c);

560:     /* recurse to the next level */
561:     SNESSetInitialFunction(next, F_c);
562:     SNESSolve(next, B_c, X_c);
563:     SNESGetConvergedReason(next,&reason);
564:     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
565:       snes->reason = SNES_DIVERGED_INNER;
566:       return(0);
567:     }
568:     /* correct as x <- x + I(x^c - Rx)*/
569:     VecAXPY(X_c, -1.0, Xo_c);

571:     if (fasc->eventinterprestrict) {PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);}
572:     MatInterpolateAdd(interpolate, X_c, X, X_new);
573:     PetscObjectSetName((PetscObject) X_c, "Coarse correction");
574:     VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view");
575:     PetscObjectSetName((PetscObject) X_new, "Updated Fine solution");
576:     VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view");
577:     if (fasc->eventinterprestrict) {PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);}
578:   }
579:   return(0);
580: }

582: /*

584: The additive cycle looks like:

586: xhat = x
587: xhat = dS(x, b)
588: x = coarsecorrection(xhat, b_d)
589: x = x + nu*(xhat - x);
590: (optional) x = uS(x, b)

592: With the coarse RHS (defect correction) as below.

594:  */
595: static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
596: {
597:   Vec                  F, B, Xhat;
598:   Vec                  X_c, Xo_c, F_c, B_c;
599:   PetscErrorCode       ierr;
600:   SNESConvergedReason  reason;
601:   PetscReal            xnorm, fnorm, ynorm;
602:   SNESLineSearchReason lsresult;
603:   SNES                 next;
604:   Mat                  restrct, interpolate;
605:   SNES_FAS             *fas = (SNES_FAS*)snes->data,*fasc;

608:   SNESFASCycleGetCorrection(snes, &next);
609:   F    = snes->vec_func;
610:   B    = snes->vec_rhs;
611:   Xhat = snes->work[1];
612:   VecCopy(X, Xhat);
613:   /* recurse first */
614:   if (next) {
615:     fasc = (SNES_FAS*)next->data;
616:     SNESFASCycleGetRestriction(snes, &restrct);
617:     SNESFASCycleGetInterpolation(snes, &interpolate);
618:     if (fas->eventresidual) {PetscLogEventBegin(fas->eventresidual,snes,0,0,0);}
619:     SNESComputeFunction(snes, Xhat, F);
620:     if (fas->eventresidual) {PetscLogEventEnd(fas->eventresidual,snes,0,0,0);}
621:     VecNorm(F, NORM_2, &fnorm);
622:     X_c  = next->vec_sol;
623:     Xo_c = next->work[0];
624:     F_c  = next->vec_func;
625:     B_c  = next->vec_rhs;

627:     SNESFASRestrict(snes,Xhat,Xo_c);
628:     /* restrict the defect */
629:     MatRestrict(restrct, F, B_c);

631:     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
632:     if (fasc->eventresidual) {PetscLogEventBegin(fasc->eventresidual,next,0,0,0);}
633:     SNESComputeFunction(next, Xo_c, F_c);
634:     if (fasc->eventresidual) {PetscLogEventEnd(fasc->eventresidual,next,0,0,0);}
635:     VecCopy(B_c, X_c);
636:     VecCopy(F_c, B_c);
637:     VecCopy(X_c, F_c);
638:     /* set initial guess of the coarse problem to the projected fine solution */
639:     VecCopy(Xo_c, X_c);

641:     /* recurse */
642:     SNESSetInitialFunction(next, F_c);
643:     SNESSolve(next, B_c, X_c);

645:     /* smooth on this level */
646:     SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);

648:     SNESGetConvergedReason(next,&reason);
649:     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
650:       snes->reason = SNES_DIVERGED_INNER;
651:       return(0);
652:     }

654:     /* correct as x <- x + I(x^c - Rx)*/
655:     VecAYPX(X_c, -1.0, Xo_c);
656:     MatInterpolate(interpolate, X_c, Xhat);

658:     /* additive correction of the coarse direction*/
659:     SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);
660:     SNESLineSearchGetReason(snes->linesearch, &lsresult);
661:     SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);
662:     if (lsresult) {
663:       if (++snes->numFailures >= snes->maxFailures) {
664:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
665:         return(0);
666:       }
667:     }
668:   } else {
669:     SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);
670:   }
671:   return(0);
672: }

674: /*

676: Defines the FAS cycle as:

678: fine problem: F(x) = b
679: coarse problem: F^c(x) = b^c

681: b^c = F^c(Rx) - R(F(x) - b)

683: correction:

685: x = x + I(x^c - Rx)

687:  */
688: static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
689: {

692:   Vec            F,B;
693:   SNES           next;

696:   F = snes->vec_func;
697:   B = snes->vec_rhs;
698:   /* pre-smooth -- just update using the pre-smoother */
699:   SNESFASCycleGetCorrection(snes, &next);
700:   SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);
701:   if (next) {
702:     SNESFASCoarseCorrection(snes, X, F, X);
703:     SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);
704:   }
705:   return(0);
706: }

708: static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
709: {
710:   SNES           next;
711:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
712:   PetscBool      isFine;

716:   /* pre-smooth -- just update using the pre-smoother */
717:   SNESFASCycleIsFine(snes,&isFine);
718:   SNESFASCycleGetCorrection(snes,&next);
719:   fas->full_stage = 0;
720:   if (next) {SNESFASCycleSetupPhase_Full(next);}
721:   return(0);
722: }

724: static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
725: {
727:   Vec            F,B;
728:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
729:   PetscBool      isFine;
730:   SNES           next;

733:   F = snes->vec_func;
734:   B = snes->vec_rhs;
735:   SNESFASCycleIsFine(snes,&isFine);
736:   SNESFASCycleGetCorrection(snes,&next);

738:   if (isFine) {
739:     SNESFASCycleSetupPhase_Full(snes);
740:   }

742:   if (fas->full_stage == 0) {
743:     /* downsweep */
744:     if (next) {
745:       if (fas->level != 1) next->max_its += 1;
746:       if (fas->full_downsweep) {SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);}
747:       fas->full_downsweep = PETSC_TRUE;
748:       SNESFASCoarseCorrection(snes,X,F,X);
749:       SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);
750:       if (fas->level != 1) next->max_its -= 1;
751:     } else {
752:       /* The smoother on the coarse level is the coarse solver */
753:       SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);
754:     }
755:     fas->full_stage = 1;
756:   } else if (fas->full_stage == 1) {
757:     if (snes->iter == 0) {SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);}
758:     if (next) {
759:       SNESFASCoarseCorrection(snes,X,F,X);
760:       SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);
761:     }
762:   }
763:   /* final v-cycle */
764:   if (isFine) {
765:     if (next) {
766:       SNESFASCoarseCorrection(snes,X,F,X);
767:       SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);
768:     }
769:   }
770:   return(0);
771: }

773: static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
774: {
776:   Vec            F,B;
777:   SNES           next;

780:   F = snes->vec_func;
781:   B = snes->vec_rhs;
782:   SNESFASCycleGetCorrection(snes,&next);
783:   if (next) {
784:     SNESFASCoarseCorrection(snes,X,F,X);
785:     SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);
786:   } else {
787:     SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);
788:   }
789:   return(0);
790: }

792: PetscBool SNEScite = PETSC_FALSE;
793: const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
794:                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
795:                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
796:                             "  year = 2013,\n"
797:                             "  type = Preprint,\n"
798:                             "  number = {ANL/MCS-P2010-0112},\n"
799:                             "  institution = {Argonne National Laboratory}\n}\n";

801: static PetscErrorCode SNESSolve_FAS(SNES snes)
802: {
804:   PetscInt       i;
805:   Vec            X, F;
806:   PetscReal      fnorm;
807:   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
808:   DM             dm;
809:   PetscBool      isFine;

812:   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);

814:   PetscCitationsRegister(SNESCitation,&SNEScite);
815:   snes->reason = SNES_CONVERGED_ITERATING;
816:   X            = snes->vec_sol;
817:   F            = snes->vec_func;

819:   SNESFASCycleIsFine(snes, &isFine);
820:   /* norm setup */
821:   PetscObjectSAWsTakeAccess((PetscObject)snes);
822:   snes->iter = 0;
823:   snes->norm = 0;
824:   PetscObjectSAWsGrantAccess((PetscObject)snes);
825:   if (!snes->vec_func_init_set) {
826:     if (fas->eventresidual) {PetscLogEventBegin(fas->eventresidual,snes,0,0,0);}
827:     SNESComputeFunction(snes,X,F);
828:     if (fas->eventresidual) {PetscLogEventEnd(fas->eventresidual,snes,0,0,0);}
829:   } else snes->vec_func_init_set = PETSC_FALSE;

831:   VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F||  */
832:   SNESCheckFunctionNorm(snes,fnorm);
833:   PetscObjectSAWsTakeAccess((PetscObject)snes);
834:   snes->norm = fnorm;
835:   PetscObjectSAWsGrantAccess((PetscObject)snes);
836:   SNESLogConvergenceHistory(snes,fnorm,0);
837:   SNESMonitor(snes,snes->iter,fnorm);

839:   /* test convergence */
840:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
841:   if (snes->reason) return(0);


844:   if (isFine) {
845:     /* propagate scale-dependent data up the hierarchy */
846:     SNESGetDM(snes,&dm);
847:     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
848:       DM dmcoarse;
849:       SNESGetDM(ffas->next,&dmcoarse);
850:       DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);
851:       dm   = dmcoarse;
852:     }
853:   }

855:   for (i = 0; i < snes->max_its; i++) {
856:     /* Call general purpose update function */
857:     if (snes->ops->update) {
858:       (*snes->ops->update)(snes, snes->iter);
859:     }

861:     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
862:       SNESFASCycle_Multiplicative(snes, X);
863:     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
864:       SNESFASCycle_Additive(snes, X);
865:     } else if (fas->fastype == SNES_FAS_FULL) {
866:       SNESFASCycle_Full(snes, X);
867:     } else if (fas->fastype == SNES_FAS_KASKADE) {
868:       SNESFASCycle_Kaskade(snes, X);
869:     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");

871:     /* check for FAS cycle divergence */
872:     if (snes->reason != SNES_CONVERGED_ITERATING) return(0);

874:     /* Monitor convergence */
875:     PetscObjectSAWsTakeAccess((PetscObject)snes);
876:     snes->iter = i+1;
877:     PetscObjectSAWsGrantAccess((PetscObject)snes);
878:     SNESLogConvergenceHistory(snes,snes->norm,0);
879:     SNESMonitor(snes,snes->iter,snes->norm);
880:     /* Test for convergence */
881:     if (isFine) {
882:       (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);
883:       if (snes->reason) break;
884:     }
885:   }
886:   if (i == snes->max_its) {
887:     PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", i);
888:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
889:   }
890:   return(0);
891: }

893: /*MC

895: SNESFAS - Full Approximation Scheme nonlinear multigrid solver.

897:    The nonlinear problem is solved by correction using coarse versions
898:    of the nonlinear problem.  This problem is perturbed so that a projected
899:    solution of the fine problem elicits no correction from the coarse problem.

901: Options Database:
902: +   -snes_fas_levels -  The number of levels
903: .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
904: .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  Additive or multiplicative cycle
905: .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
906: .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
907: .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
908: .   -snes_fas_monitor -  Monitor progress of all of the levels
909: .   -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
910: .   -fas_levels_snes_ -  SNES options for all smoothers
911: .   -fas_levels_cycle_snes_ -  SNES options for all cycles
912: .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
913: .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
914: -   -fas_coarse_snes_ -  SNES options for the coarsest smoother

916: Notes:
917:    The organization of the FAS solver is slightly different from the organization of PCMG
918:    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
919:    The cycle SNES instance may be used for monitoring convergence on a particular level.

921: Level: beginner

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

927: .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
928: M*/

930: PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
931: {
932:   SNES_FAS       *fas;

936:   snes->ops->destroy        = SNESDestroy_FAS;
937:   snes->ops->setup          = SNESSetUp_FAS;
938:   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
939:   snes->ops->view           = SNESView_FAS;
940:   snes->ops->solve          = SNESSolve_FAS;
941:   snes->ops->reset          = SNESReset_FAS;

943:   snes->usesksp = PETSC_FALSE;
944:   snes->usesnpc = PETSC_FALSE;

946:   if (!snes->tolerancesset) {
947:     snes->max_funcs = 30000;
948:     snes->max_its   = 10000;
949:   }

951:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

953:   PetscNewLog(snes,&fas);

955:   snes->data                  = (void*) fas;
956:   fas->level                  = 0;
957:   fas->levels                 = 1;
958:   fas->n_cycles               = 1;
959:   fas->max_up_it              = 1;
960:   fas->max_down_it            = 1;
961:   fas->smoothu                = NULL;
962:   fas->smoothd                = NULL;
963:   fas->next                   = NULL;
964:   fas->previous               = NULL;
965:   fas->fine                   = snes;
966:   fas->interpolate            = NULL;
967:   fas->restrct                = NULL;
968:   fas->inject                 = NULL;
969:   fas->usedmfornumberoflevels = PETSC_FALSE;
970:   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
971:   fas->full_downsweep         = PETSC_FALSE;

973:   fas->eventsmoothsetup    = 0;
974:   fas->eventsmoothsolve    = 0;
975:   fas->eventresidual       = 0;
976:   fas->eventinterprestrict = 0;
977:   return(0);
978: }