Actual source code: alpha.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: /*
  2:   Code for timestepping with implicit generalized-\alpha method
  3:   for first order systems.
  4: */
  5: #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/

  7: typedef PetscErrorCode (*TSAlphaAdaptFunction)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);

  9: typedef struct {
 10:   Vec       X0,Xa,X1;
 11:   Vec       V0,Va,V1;
 12:   Vec       R,E;
 13:   PetscReal Alpha_m;
 14:   PetscReal Alpha_f;
 15:   PetscReal Gamma;
 16:   PetscReal stage_time;
 17:   PetscReal shift;

 19:   TSAlphaAdaptFunction adapt;
 20:   void                 *adaptctx;
 21:   PetscReal            rtol;
 22:   PetscReal            atol;
 23:   PetscReal            rho;
 24:   PetscReal            scale_min;
 25:   PetscReal            scale_max;
 26:   PetscReal            dt_min;
 27:   PetscReal            dt_max;
 28: } TS_Alpha;

 32: static PetscErrorCode TSStep_Alpha(TS ts)
 33: {
 34:   TS_Alpha            *th = (TS_Alpha*)ts->data;
 35:   PetscInt            its,lits,reject;
 36:   PetscReal           next_time_step;
 37:   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
 38:   PetscErrorCode      ierr;

 41:   if (ts->steps == 0) {
 42:     VecSet(th->V0,0.0);
 43:   } else {
 44:     VecCopy(th->V1,th->V0);
 45:   }
 46:   VecCopy(ts->vec_sol,th->X0);
 47:   next_time_step = ts->time_step;
 48:   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
 49:     ts->time_step  = next_time_step;
 50:     th->stage_time = ts->ptime + th->Alpha_f*ts->time_step;
 51:     th->shift      = th->Alpha_m/(th->Alpha_f*th->Gamma*ts->time_step);
 52:     TSPreStep(ts);
 53:     TSPreStage(ts,th->stage_time);
 54:     /* predictor */
 55:     VecCopy(th->X0,th->X1);
 56:     /* solve R(X,V) = 0 */
 57:     SNESSolve(ts->snes,NULL,th->X1);
 58:     /* V1 = (1-1/Gamma)*V0 + 1/(Gamma*dT)*(X1-X0) */
 59:     VecWAXPY(th->V1,-1,th->X0,th->X1);
 60:     VecAXPBY(th->V1,1-1/th->Gamma,1/(th->Gamma*ts->time_step),th->V0);
 61:     TSPostStage(ts,th->stage_time,0,&(th->V1));
 62:     /* nonlinear solve convergence */
 63:     SNESGetConvergedReason(ts->snes,&snesreason);
 64:     if (snesreason < 0 && !th->adapt) break;
 65:     SNESGetIterationNumber(ts->snes,&its);
 66:     SNESGetLinearSolveIterations(ts->snes,&lits);
 67:     ts->snes_its += its; ts->ksp_its += lits;
 68:     PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);
 69:     /* time step adaptativity */
 70:     if (!th->adapt) break;
 71:     else {
 72:       PetscReal t1     = ts->ptime + ts->time_step;
 73:       PetscBool stepok = (reject==0) ? PETSC_TRUE : PETSC_FALSE;
 74:       th->adapt(ts,t1,th->X1,th->V1,&next_time_step,&stepok,th->adaptctx);
 75:       PetscInfo5(ts,"Step %D (t=%g,dt=%g) %s, next dt=%g\n",ts->steps,(double)ts->ptime,(double)ts->time_step,stepok?"accepted":"rejected",(double)next_time_step);
 76:       if (stepok) break;
 77:     }
 78:   }
 79:   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
 80:     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
 81:     PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
 82:     return(0);
 83:   }
 84:   if (reject >= ts->max_reject) {
 85:     ts->reason = TS_DIVERGED_STEP_REJECTED;
 86:     PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);
 87:     return(0);
 88:   }
 89:   VecCopy(th->X1,ts->vec_sol);
 90:   ts->ptime    += ts->time_step;
 91:   ts->time_step = next_time_step;
 92:   ts->steps++;
 93:   return(0);
 94: }

 98: static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X)
 99: {
100:   TS_Alpha       *th = (TS_Alpha*)ts->data;
101:   PetscReal      dt  = t - ts->ptime;

105:   VecCopy(ts->vec_sol,X);
106:   VecAXPY(X,th->Gamma*dt,th->V1);
107:   VecAXPY(X,(1-th->Gamma)*dt,th->V0);
108:   return(0);
109: }

111: /*------------------------------------------------------------*/
114: static PetscErrorCode TSReset_Alpha(TS ts)
115: {
116:   TS_Alpha       *th = (TS_Alpha*)ts->data;

120:   VecDestroy(&th->X0);
121:   VecDestroy(&th->Xa);
122:   VecDestroy(&th->X1);
123:   VecDestroy(&th->V0);
124:   VecDestroy(&th->Va);
125:   VecDestroy(&th->V1);
126:   VecDestroy(&th->E);
127:   return(0);
128: }

132: static PetscErrorCode TSDestroy_Alpha(TS ts)
133: {

137:   TSReset_Alpha(ts);
138:   PetscFree(ts->data);

140:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",NULL);
141:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",NULL);
142:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",NULL);
143:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetAdapt_C",NULL);
144:   return(0);
145: }

149: static PetscErrorCode SNESTSFormFunction_Alpha(SNES snes,Vec x,Vec y,TS ts)
150: {
151:   TS_Alpha       *th = (TS_Alpha*)ts->data;
152:   Vec            X0  = th->X0, V0 = th->V0;
153:   Vec            X1  = x, V1 = th->V1, R = y;

157:   /* V1 = (1-1/Gamma)*V0 + 1/(Gamma*dT)*(X1-X0) */
158:   VecWAXPY(V1,-1,X0,X1);
159:   VecAXPBY(V1,1-1/th->Gamma,1/(th->Gamma*ts->time_step),V0);
160:   /* Xa = X0 + Alpha_f*(X1-X0) */
161:   VecWAXPY(th->Xa,-1,X0,X1);
162:   VecAYPX(th->Xa,th->Alpha_f,X0);
163:   /* Va = V0 + Alpha_m*(V1-V0) */
164:   VecWAXPY(th->Va,-1,V0,V1);
165:   VecAYPX(th->Va,th->Alpha_m,V0);
166:   /* F = Function(ta,Xa,Va) */
167:   TSComputeIFunction(ts,th->stage_time,th->Xa,th->Va,R,PETSC_FALSE);
168:   VecScale(R,1/th->Alpha_f);
169:   return(0);
170: }

174: static PetscErrorCode SNESTSFormJacobian_Alpha(SNES snes,Vec x,Mat A,Mat B,TS ts)
175: {
176:   TS_Alpha       *th = (TS_Alpha*)ts->data;

180:   /* A,B = Jacobian(ta,Xa,Va) */
181:   TSComputeIJacobian(ts,th->stage_time,th->Xa,th->Va,th->shift,A,B,PETSC_FALSE);
182:   return(0);
183: }

187: static PetscErrorCode TSSetUp_Alpha(TS ts)
188: {
189:   TS_Alpha       *th = (TS_Alpha*)ts->data;

193:   VecDuplicate(ts->vec_sol,&th->X0);
194:   VecDuplicate(ts->vec_sol,&th->Xa);
195:   VecDuplicate(ts->vec_sol,&th->X1);
196:   VecDuplicate(ts->vec_sol,&th->V0);
197:   VecDuplicate(ts->vec_sol,&th->Va);
198:   VecDuplicate(ts->vec_sol,&th->V1);
199:   return(0);
200: }

204: static PetscErrorCode TSSetFromOptions_Alpha(TS ts)
205: {
206:   TS_Alpha       *th = (TS_Alpha*)ts->data;

210:   PetscOptionsHead("Alpha ODE solver options");
211:   {
212:     PetscBool flag, adapt = PETSC_FALSE;
213:     PetscReal radius = 1.0;
214:     PetscOptionsReal("-ts_alpha_radius","spectral radius","TSAlphaSetRadius",radius,&radius,&flag);
215:     if (flag) { TSAlphaSetRadius(ts,radius); }
216:     PetscOptionsReal("-ts_alpha_alpha_m","algoritmic parameter alpha_m","TSAlphaSetParams",th->Alpha_m,&th->Alpha_m,NULL);
217:     PetscOptionsReal("-ts_alpha_alpha_f","algoritmic parameter alpha_f","TSAlphaSetParams",th->Alpha_f,&th->Alpha_f,NULL);
218:     PetscOptionsReal("-ts_alpha_gamma","algoritmic parameter gamma","TSAlphaSetParams",th->Gamma,&th->Gamma,NULL);
219:     TSAlphaSetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma);

221:     PetscOptionsBool("-ts_alpha_adapt","default time step adaptativity","TSAlphaSetAdapt",adapt,&adapt,&flag);
222:     if (flag) { TSAlphaSetAdapt(ts,adapt ? TSAlphaAdaptDefault : NULL,NULL); }
223:     PetscOptionsReal("-ts_alpha_adapt_rtol","relative tolerance for dt adaptativity","",th->rtol,&th->rtol,NULL);
224:     PetscOptionsReal("-ts_alpha_adapt_atol","absolute tolerance for dt adaptativity","",th->atol,&th->atol,NULL);
225:     PetscOptionsReal("-ts_alpha_adapt_min","minimum dt scale","",th->scale_min,&th->scale_min,NULL);
226:     PetscOptionsReal("-ts_alpha_adapt_max","maximum dt scale","",th->scale_max,&th->scale_max,NULL);
227:     PetscOptionsReal("-ts_alpha_adapt_dt_min","minimum dt","",th->dt_min,&th->dt_min,NULL);
228:     PetscOptionsReal("-ts_alpha_adapt_dt_max","maximum dt","",th->dt_max,&th->dt_max,NULL);
229:     SNESSetFromOptions(ts->snes);
230:   }
231:   PetscOptionsTail();
232:   return(0);
233: }

237: static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer)
238: {
239:   TS_Alpha       *th = (TS_Alpha*)ts->data;
240:   PetscBool      iascii;

244:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
245:   if (iascii) {
246:     PetscViewerASCIIPrintf(viewer,"  Alpha_m=%g, Alpha_f=%g, Gamma=%g\n",(double)th->Alpha_m,(double)th->Alpha_f,(double)th->Gamma);
247:   }
248:   SNESView(ts->snes,viewer);
249:   return(0);
250: }

252: /*------------------------------------------------------------*/

256: PetscErrorCode  TSAlphaSetRadius_Alpha(TS ts,PetscReal radius)
257: {
258:   TS_Alpha *th = (TS_Alpha*)ts->data;

261:   if (radius < 0 || radius > 1) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius);
262:   th->Alpha_m = 0.5*(3-radius)/(1+radius);
263:   th->Alpha_f = 1/(1+radius);
264:   th->Gamma   = 0.5 + th->Alpha_m - th->Alpha_f;
265:   return(0);
266: }

270: PetscErrorCode  TSAlphaSetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)
271: {
272:   TS_Alpha *th = (TS_Alpha*)ts->data;

275:   th->Alpha_m = alpha_m;
276:   th->Alpha_f = alpha_f;
277:   th->Gamma   = gamma;
278:   return(0);
279: }

283: PetscErrorCode  TSAlphaGetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma)
284: {
285:   TS_Alpha *th = (TS_Alpha*)ts->data;

288:   if (alpha_m) *alpha_m = th->Alpha_m;
289:   if (alpha_f) *alpha_f = th->Alpha_f;
290:   if (gamma)   *gamma   = th->Gamma;
291:   return(0);
292: }

296: PetscErrorCode  TSAlphaSetAdapt_Alpha(TS ts,TSAlphaAdaptFunction adapt,void *ctx)
297: {
298:   TS_Alpha *th = (TS_Alpha*)ts->data;

301:   th->adapt    = adapt;
302:   th->adaptctx = ctx;
303:   return(0);
304: }

306: /* ------------------------------------------------------------ */
307: /*MC
308:       TSALPHA - DAE solver using the implicit Generalized-Alpha method

310:   Level: beginner

312:   References:
313:   K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha
314:   method for integrating the filtered Navier-Stokes equations with a
315:   stabilized finite element method", Computer Methods in Applied
316:   Mechanics and Engineering, 190, 305-319, 2000.
317:   DOI: 10.1016/S0045-7825(00)00203-6.

319:   J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
320:   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
321:   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.

323: .seealso:  TSCreate(), TS, TSSetType()

325: M*/
328: PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts)
329: {
330:   TS_Alpha       *th;

334:   ts->ops->reset          = TSReset_Alpha;
335:   ts->ops->destroy        = TSDestroy_Alpha;
336:   ts->ops->view           = TSView_Alpha;
337:   ts->ops->setup          = TSSetUp_Alpha;
338:   ts->ops->step           = TSStep_Alpha;
339:   ts->ops->interpolate    = TSInterpolate_Alpha;
340:   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
341:   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
342:   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;

344:   PetscNewLog(ts,&th);
345:   ts->data = (void*)th;

347:   th->Alpha_m = 0.5;
348:   th->Alpha_f = 0.5;
349:   th->Gamma   = 0.5;

351:   th->rtol      = 1e-3;
352:   th->atol      = 1e-3;
353:   th->rho       = 0.9;
354:   th->scale_min = 0.1;
355:   th->scale_max = 5.0;
356:   th->dt_min    = 0.0;
357:   th->dt_max    = PETSC_MAX_REAL;

359:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetAdapt_C",TSAlphaSetAdapt_Alpha);
360:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",TSAlphaSetRadius_Alpha);
361:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",TSAlphaSetParams_Alpha);
362:   PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",TSAlphaGetParams_Alpha);
363:   return(0);
364: }

368: /*@C
369:   TSAlphaSetAdapt - sets the time step adaptativity and acceptance test routine

371:   This function allows to accept/reject a step and select the
372:   next time step to use.

374:   Not Collective

376:   Input Parameter:
377: +  ts - timestepping context
378: .  adapt - user-defined adapt routine
379: -  ctx  - [optional] user-defined context for private data for the
380:          adapt routine (may be NULL)

382:    Calling sequence of adapt:
383: $    adapt (TS ts,PetscReal t,Vec X,Vec Xdot,
384: $            PetscReal *next_dt,PetscBool *accepted,void *ctx);

386:   Level: intermediate

388: @*/
389: PetscErrorCode  TSAlphaSetAdapt(TS ts,TSAlphaAdaptFunction adapt,void *ctx)
390: {

395:   PetscTryMethod(ts,"TSAlphaSetAdapt_C",(TS,TSAlphaAdaptFunction,void*),(ts,adapt,ctx));
396:   return(0);
397: }

401: PetscErrorCode  TSAlphaAdaptDefault(TS ts,PetscReal t,Vec X,Vec Xdot, PetscReal *nextdt,PetscBool *ok,void *ctx)
402: {
403:   TS_Alpha            *th;
404:   SNESConvergedReason snesreason;
405:   PetscReal           dt,normX,normE,Emax,scale;
406:   PetscErrorCode      ierr;

410: #if PETSC_USE_DEBUG
411:   {
412:     PetscBool match;
413:     PetscObjectTypeCompare((PetscObject)ts,TSALPHA,&match);
414:     if (!match) SETERRQ(PetscObjectComm((PetscObject)ts),1,"Only for TSALPHA");
415:   }
416: #endif
417:   th = (TS_Alpha*)ts->data;

419:   SNESGetConvergedReason(ts->snes,&snesreason);
420:   if (snesreason < 0) {
421:     *ok      = PETSC_FALSE;
422:     *nextdt *= th->scale_min;
423:     goto finally;
424:   }

426:   /* first-order aproximation to the local error */
427:   /* E = (X0 + dt*Xdot) - X */
428:   TSGetTimeStep(ts,&dt);
429:   if (!th->E) {VecDuplicate(th->X0,&th->E);}
430:   VecWAXPY(th->E,dt,Xdot,th->X0);
431:   VecAXPY(th->E,-1,X);
432:   VecNorm(th->E,NORM_2,&normE);
433:   /* compute maximum allowable error */
434:   VecNorm(X,NORM_2,&normX);
435:   if (normX == 0) {VecNorm(th->X0,NORM_2,&normX);}
436:   Emax =  th->rtol * normX + th->atol;
437:   /* compute next time step */
438:   if (normE > 0) {
439:     scale = th->rho * PetscRealPart(PetscSqrtScalar((PetscScalar)(Emax/normE)));
440:     scale = PetscMax(scale,th->scale_min);
441:     scale = PetscMin(scale,th->scale_max);
442:     if (!(*ok)) scale = PetscMin(1.0,scale);
443:     *nextdt *= scale;
444:   }
445:   /* accept or reject step */
446:   if (normE <= Emax) *ok = PETSC_TRUE;
447:   else               *ok = PETSC_FALSE;

449: finally:
450:   *nextdt = PetscMax(*nextdt,th->dt_min);
451:   *nextdt = PetscMin(*nextdt,th->dt_max);
452:   return(0);
453: }

457: /*@
458:   TSAlphaSetRadius - sets the desired spectral radius of the method
459:                      (i.e. high-frequency numerical damping)

461:   Logically Collective on TS

463:   The algorithmic parameters \alpha_m and \alpha_f of the
464:   generalized-\alpha method can be computed in terms of a specified
465:   spectral radius \rho in [0,1] for infinite time step in order to
466:   control high-frequency numerical damping:
467:     alpha_m = 0.5*(3-\rho)/(1+\rho)
468:     alpha_f = 1/(1+\rho)

470:   Input Parameter:
471: +  ts - timestepping context
472: -  radius - the desired spectral radius

474:   Options Database:
475: .  -ts_alpha_radius <radius>

477:   Level: intermediate

479: .seealso: TSAlphaSetParams(), TSAlphaGetParams()
480: @*/
481: PetscErrorCode  TSAlphaSetRadius(TS ts,PetscReal radius)
482: {

487:   PetscTryMethod(ts,"TSAlphaSetRadius_C",(TS,PetscReal),(ts,radius));
488:   return(0);
489: }

493: /*@
494:   TSAlphaSetParams - sets the algorithmic parameters for TSALPHA

496:   Not Collective

498:   Second-order accuracy can be obtained so long as:
499:     \gamma = 0.5 + alpha_m - alpha_f

501:   Unconditional stability requires:
502:     \alpha_m >= \alpha_f >= 0.5

504:   Backward Euler method is recovered when:
505:     \alpha_m = \alpha_f = gamma = 1


508:   Input Parameter:
509: +  ts - timestepping context
510: .  \alpha_m - algorithmic paramenter
511: .  \alpha_f - algorithmic paramenter
512: -  \gamma   - algorithmic paramenter

514:    Options Database:
515: +  -ts_alpha_alpha_m <alpha_m>
516: .  -ts_alpha_alpha_f <alpha_f>
517: -  -ts_alpha_gamma <gamma>

519:   Note:
520:   Use of this function is normally only required to hack TSALPHA to
521:   use a modified integration scheme. Users should call
522:   TSAlphaSetRadius() to set the desired spectral radius of the methods
523:   (i.e. high-frequency damping) in order so select optimal values for
524:   these parameters.

526:   Level: advanced

528: .seealso: TSAlphaSetRadius(), TSAlphaGetParams()
529: @*/
530: PetscErrorCode  TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)
531: {

536:   PetscTryMethod(ts,"TSAlphaSetParams_C",(TS,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma));
537:   return(0);
538: }

542: /*@
543:   TSAlphaGetParams - gets the algorithmic parameters for TSALPHA

545:   Not Collective

547:   Input Parameter:
548: +  ts - timestepping context
549: .  \alpha_m - algorithmic parameter
550: .  \alpha_f - algorithmic parameter
551: -  \gamma   - algorithmic parameter

553:   Note:
554:   Use of this function is normally only required to hack TSALPHA to
555:   use a modified integration scheme. Users should call
556:   TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral
557:   radius of the method) in order so select optimal values for these
558:   parameters.

560:   Level: advanced

562: .seealso: TSAlphaSetRadius(), TSAlphaSetParams()
563: @*/
564: PetscErrorCode  TSAlphaGetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma)
565: {

573:   PetscUseMethod(ts,"TSAlphaGetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma));
574:   return(0);
575: }