Actual source code: alpha.c

  1: /*
  2:   Code for timestepping with implicit generalized-\alpha method
  3:   for first order systems.
  4: */
  5: #include <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:     /* predictor */
 53:     VecCopy(th->X0,th->X1);
 54:     /* solve R(X,V) = 0 */
 55:     SNESSolve(ts->snes,PETSC_NULL,th->X1);
 56:     /* V1 = (1-1/Gamma)*V0 + 1/(Gamma*dT)*(X1-X0) */
 57:     VecWAXPY(th->V1,-1,th->X0,th->X1);
 58:     VecAXPBY(th->V1,1-1/th->Gamma,1/(th->Gamma*ts->time_step),th->V0);
 59:     /* nonlinear solve convergence */
 60:     SNESGetConvergedReason(ts->snes,&snesreason);
 61:     if (snesreason < 0 && !th->adapt) break;
 62:     SNESGetIterationNumber(ts->snes,&its);
 63:     SNESGetLinearSolveIterations(ts->snes,&lits);
 64:     ts->nonlinear_its += its; ts->linear_its += lits;
 65:     PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);
 66:     /* time step adaptativity */
 67:     if (!th->adapt) break;
 68:     else {
 69:       PetscReal t1 = ts->ptime + ts->time_step;
 70:       PetscBool stepok = (reject==0) ? PETSC_TRUE : PETSC_FALSE;
 71:       th->adapt(ts,t1,th->X1,th->V1,&next_time_step,&stepok,th->adaptctx);
 72:       PetscInfo5(ts,"Step %D (t=%G,dt=%G) %s, next dt=%G\n",ts->steps,ts->ptime,ts->time_step,stepok?"accepted":"rejected",next_time_step);
 73:       if (stepok) break;
 74:     }
 75:   }
 76:   if (snesreason < 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
 77:     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
 78:     PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
 79:     return(0);
 80:   }
 81:   if (reject >= ts->max_reject) {
 82:     ts->reason = TS_DIVERGED_STEP_REJECTED;
 83:     PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);
 84:     return(0);
 85:   }
 86:   VecCopy(th->X1,ts->vec_sol);
 87:   ts->ptime += ts->time_step;
 88:   ts->time_step = next_time_step;
 89:   ts->steps++;
 90:   return(0);
 91: }

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

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

108: /*------------------------------------------------------------*/
111: static PetscErrorCode TSReset_Alpha(TS ts)
112: {
113:   TS_Alpha       *th = (TS_Alpha*)ts->data;
114:   PetscErrorCode  ierr;

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

129: static PetscErrorCode TSDestroy_Alpha(TS ts)
130: {
131:   PetscErrorCode  ierr;

134:   TSReset_Alpha(ts);
135:   PetscFree(ts->data);

137:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaSetRadius_C","",PETSC_NULL);
138:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaSetParams_C","",PETSC_NULL);
139:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaGetParams_C","",PETSC_NULL);
140:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaSetAdapt_C","",PETSC_NULL);
141:   return(0);
142: }

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

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

171: static PetscErrorCode SNESTSFormJacobian_Alpha(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
172: {
173:   TS_Alpha       *th = (TS_Alpha*)ts->data;

177:   /* A,B = Jacobian(ta,Xa,Va) */
178:   TSComputeIJacobian(ts,th->stage_time,th->Xa,th->Va,th->shift,A,B,str,PETSC_FALSE);
179:   return(0);
180: }

184: static PetscErrorCode TSSetUp_Alpha(TS ts)
185: {
186:   TS_Alpha *th = (TS_Alpha*)ts->data;

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

201: static PetscErrorCode TSSetFromOptions_Alpha(TS ts)
202: {
203:   TS_Alpha *th = (TS_Alpha*)ts->data;

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

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

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

241:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
242:   if (iascii) {
243:     PetscViewerASCIIPrintf(viewer,"  Alpha_m=%G, Alpha_f=%G, Gamma=%G\n",th->Alpha_m,th->Alpha_f,th->Gamma);
244:   }
245:   SNESView(ts->snes,viewer);
246:   return(0);
247: }

249: /*------------------------------------------------------------*/


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

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

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

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

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

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

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

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


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

311:   Level: beginner

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

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

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

326: M*/
330: PetscErrorCode  TSCreate_Alpha(TS ts)
331: {
332:   TS_Alpha       *th;

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

346:   PetscNewLog(ts,TS_Alpha,&th);
347:   ts->data = (void*)th;

349:   th->Alpha_m = 0.5;
350:   th->Alpha_f = 0.5;
351:   th->Gamma   = 0.5;

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

361:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaSetAdapt_C","TSAlphaSetAdapt_Alpha",TSAlphaSetAdapt_Alpha);
362:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaSetRadius_C","TSAlphaSetRadius_Alpha",TSAlphaSetRadius_Alpha);
363:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaSetParams_C","TSAlphaSetParams_Alpha",TSAlphaSetParams_Alpha);
364:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSAlphaGetParams_C","TSAlphaGetParams_Alpha",TSAlphaGetParams_Alpha);

366:   return(0);
367: }

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

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

378:   Not Collective

380:   Input Parameter:
381: +  ts - timestepping context
382: .  adapt - user-defined adapt routine
383: -  ctx  - [optional] user-defined context for private data for the
384:          adapt routine (may be PETSC_NULL)

386:    Calling sequence of adapt:
387: $    adapt (TS ts,PetscReal t,Vec X,Vec Xdot,
388: $            PetscReal *next_dt,PetscBool *accepted,void *ctx);

390:   Level: intermediate

392: @*/
393: PetscErrorCode  TSAlphaSetAdapt(TS ts,TSAlphaAdaptFunction adapt,void *ctx)
394: {
398:   PetscTryMethod(ts,"TSAlphaSetAdapt_C",(TS,TSAlphaAdaptFunction,void*),(ts,adapt,ctx));
399:   return(0);
400: }

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

413: #if PETSC_USE_DEBUG
414:   {
415:     PetscBool match;
416:     PetscTypeCompare((PetscObject)ts,TSALPHA,&match);
417:     if (!match) SETERRQ(((PetscObject)ts)->comm,1,"Only for TSALPHA");
418:   }
419: #endif
420:   th = (TS_Alpha*)ts->data;

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

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

455:   finally:
456:   *nextdt = PetscMax(*nextdt,th->dt_min);
457:   *nextdt = PetscMin(*nextdt,th->dt_max);
458:   return(0);
459: }

463: /*@
464:   TSAlphaSetRadius - sets the desired spectral radius of the method
465:                      (i.e. high-frequency numerical damping)

467:   Logically Collective on TS

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

476:   Input Parameter:
477: +  ts - timestepping context
478: -  radius - the desired spectral radius

480:   Options Database:
481: .  -ts_alpha_radius <radius>

483:   Level: intermediate

485: .seealso: TSAlphaSetParams(), TSAlphaGetParams()
486: @*/
487: PetscErrorCode  TSAlphaSetRadius(TS ts,PetscReal radius)
488: {

493:   PetscTryMethod(ts,"TSAlphaSetRadius_C",(TS,PetscReal),(ts,radius));
494:   return(0);
495: }

499: /*@
500:   TSAlphaSetParams - sets the algorithmic parameters for TSALPHA

502:   Not Collective

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

507:   Unconditional stability requires:
508:     \alpha_m >= \alpha_f >= 0.5

510:   Backward Euler method is recovered when:
511:     \alpha_m = \alpha_f = gamma = 1


514:   Input Parameter:
515: +  ts - timestepping context
516: .  \alpha_m - algorithmic paramenter
517: .  \alpha_f - algorithmic paramenter
518: -  \gamma   - algorithmic paramenter

520:    Options Database:
521: +  -ts_alpha_alpha_m <alpha_m>
522: .  -ts_alpha_alpha_f <alpha_f>
523: -  -ts_alpha_gamma <gamma>

525:   Note:
526:   Use of this function is normally only required to hack TSALPHA to
527:   use a modified integration scheme. Users should call
528:   TSAlphaSetRadius() to set the desired spectral radius of the methods
529:   (i.e. high-frequency damping) in order so select optimal values for
530:   these parameters.

532:   Level: advanced

534: .seealso: TSAlphaSetRadius(), TSAlphaGetParams()
535: @*/
536: PetscErrorCode  TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)
537: {

542:   PetscTryMethod(ts,"TSAlphaSetParams_C",(TS,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma));
543:   return(0);
544: }

548: /*@
549:   TSAlphaGetParams - gets the algorithmic parameters for TSALPHA

551:   Not Collective

553:   Input Parameter:
554: +  ts - timestepping context
555: .  \alpha_m - algorithmic parameter
556: .  \alpha_f - algorithmic parameter
557: -  \gamma   - algorithmic parameter

559:   Note:
560:   Use of this function is normally only required to hack TSALPHA to
561:   use a modified integration scheme. Users should call
562:   TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral
563:   radius of the method) in order so select optimal values for these
564:   parameters.

566:   Level: advanced

568: .seealso: TSAlphaSetRadius(), TSAlphaSetParams()
569: @*/
570: PetscErrorCode  TSAlphaGetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma)
571: {

579:   PetscUseMethod(ts,"TSAlphaGetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma));
580:   return(0);
581: }