Actual source code: posindep.c
1: /*
2: Code for Timestepping with implicit backwards Euler.
3: */
4: #include src/ts/tsimpl.h
6: typedef struct {
7: Vec update; /* work vector where new solution is formed */
8: Vec func; /* work vector where F(t[i],u[i]) is stored */
9: Vec rhs; /* work vector for RHS; vec_sol/dt */
11: /* information used for Pseudo-timestepping */
13: PetscErrorCode (*dt)(TS,PetscReal*,void*); /* compute next timestep, and related context */
14: void *dtctx;
15: PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscTruth*); /* verify previous timestep and related context */
16: void *verifyctx;
18: PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */
19: PetscReal fnorm_previous;
21: PetscReal dt_increment; /* scaling that dt is incremented each time-step */
22: PetscTruth increment_dt_from_initial_dt;
23: } TS_Pseudo;
25: /* ------------------------------------------------------------------------------*/
29: /*@
30: TSPseudoComputeTimeStep - Computes the next timestep for a currently running
31: pseudo-timestepping process.
33: Collective on TS
35: Input Parameter:
36: . ts - timestep context
38: Output Parameter:
39: . dt - newly computed timestep
41: Level: advanced
43: Notes:
44: The routine to be called here to compute the timestep should be
45: set by calling TSPseudoSetTimeStep().
47: .keywords: timestep, pseudo, compute
49: .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
50: @*/
51: PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
52: {
53: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
57: PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
58: (*pseudo->dt)(ts,dt,pseudo->dtctx);
59: PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
60: return(0);
61: }
64: /* ------------------------------------------------------------------------------*/
67: /*@C
68: TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
70: Collective on TS
72: Input Parameters:
73: + ts - the timestep context
74: . dtctx - unused timestep context
75: - update - latest solution vector
77: Output Parameters:
78: + newdt - the timestep to use for the next step
79: - flag - flag indicating whether the last time step was acceptable
81: Level: advanced
83: Note:
84: This routine always returns a flag of 1, indicating an acceptable
85: timestep.
87: .keywords: timestep, pseudo, default, verify
89: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
90: @*/
91: PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscTruth *flag)
92: {
94: *flag = PETSC_TRUE;
95: return(0);
96: }
101: /*@
102: TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
104: Collective on TS
106: Input Parameters:
107: + ts - timestep context
108: - update - latest solution vector
110: Output Parameters:
111: + dt - newly computed timestep (if it had to shrink)
112: - flag - indicates if current timestep was ok
114: Level: advanced
116: Notes:
117: The routine to be called here to compute the timestep should be
118: set by calling TSPseudoSetVerifyTimeStep().
120: .keywords: timestep, pseudo, verify
122: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
123: @*/
124: PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscTruth *flag)
125: {
126: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
130: if (!pseudo->verify) {*flag = PETSC_TRUE; return(0);}
132: (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);
134: return(0);
135: }
137: /* --------------------------------------------------------------------------------*/
141: static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime)
142: {
143: Vec sol = ts->vec_sol;
145: PetscInt i,max_steps = ts->max_steps,its,lits;
146: PetscTruth ok;
147: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
148: PetscReal current_time_step;
149:
151: *steps = -ts->steps;
153: VecCopy(sol,pseudo->update);
154: for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
155: TSPseudoComputeTimeStep(ts,&ts->time_step);
156: TSMonitor(ts,ts->steps,ts->ptime,sol);
157: current_time_step = ts->time_step;
158: while (PETSC_TRUE) {
159: ts->ptime += current_time_step;
160: SNESSolve(ts->snes,pseudo->update);
161: SNESGetNumberLinearIterations(ts->snes,&lits);
162: SNESGetIterationNumber(ts->snes,&its);
163: ts->nonlinear_its += its; ts->linear_its += lits;
164: TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);
165: if (ok) break;
166: ts->ptime -= current_time_step;
167: current_time_step = ts->time_step;
168: }
169: VecCopy(pseudo->update,sol);
170: ts->steps++;
171: }
172: TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);
173: VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
174: TSMonitor(ts,ts->steps,ts->ptime,sol);
176: *steps += ts->steps;
177: *ptime = ts->ptime;
178: return(0);
179: }
181: /*------------------------------------------------------------*/
184: static PetscErrorCode TSDestroy_Pseudo(TS ts)
185: {
186: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
190: if (pseudo->update) {VecDestroy(pseudo->update);}
191: if (pseudo->func) {VecDestroy(pseudo->func);}
192: if (pseudo->rhs) {VecDestroy(pseudo->rhs);}
193: PetscFree(pseudo);
194: return(0);
195: }
198: /*------------------------------------------------------------*/
200: /*
201: This defines the nonlinear equation that is to be solved with SNES
203: (U^{n+1} - U^{n})/dt - F(U^{n+1})
204: */
207: PetscErrorCode TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
208: {
209: TS ts = (TS) ctx;
210: PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
212: PetscInt i,n;
215: /* apply user provided function */
216: TSComputeRHSFunction(ts,ts->ptime,x,y);
217: /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
218: VecGetArray(ts->vec_sol,&un);
219: VecGetArray(x,&unp1);
220: VecGetArray(y,&Funp1);
221: VecGetLocalSize(x,&n);
222: for (i=0; i<n; i++) {
223: Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
224: }
225: VecRestoreArray(ts->vec_sol,&un);
226: VecRestoreArray(x,&unp1);
227: VecRestoreArray(y,&Funp1);
229: return(0);
230: }
232: /*
233: This constructs the Jacobian needed for SNES
235: J = I/dt - J_{F} where J_{F} is the given Jacobian of F.
236: */
239: PetscErrorCode TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
240: {
241: TS ts = (TS) ctx;
245: /* construct users Jacobian */
246: TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);
248: /* shift and scale Jacobian */
249: TSScaleShiftMatrices(ts,*AA,*BB,*str);
250: return(0);
251: }
256: static PetscErrorCode TSSetUp_Pseudo(TS ts)
257: {
258: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
262: /* SNESSetFromOptions(ts->snes); */
263: VecDuplicate(ts->vec_sol,&pseudo->update);
264: VecDuplicate(ts->vec_sol,&pseudo->func);
265: SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);
266: SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);
267: return(0);
268: }
269: /*------------------------------------------------------------*/
273: PetscErrorCode TSPseudoDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
274: {
275: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
279: (*PetscHelpPrintf)(ts->comm,"TS %D dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);
280: return(0);
281: }
285: static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
286: {
287: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
289: PetscTruth flg;
293: PetscOptionsHead("Pseudo-timestepping options");
294: PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);
295: if (flg) {
296: TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);
297: }
298: PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);
299: if (flg) {
300: TSPseudoIncrementDtFromInitialDt(ts);
301: }
302: PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);
303: PetscOptionsTail();
305: return(0);
306: }
310: static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
311: {
313: return(0);
314: }
316: /* ----------------------------------------------------------------------------- */
319: /*@C
320: TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
321: last timestep.
323: Collective on TS
325: Input Parameters:
326: + ts - timestep context
327: . dt - user-defined function to verify timestep
328: - ctx - [optional] user-defined context for private data
329: for the timestep verification routine (may be PETSC_NULL)
331: Level: advanced
333: Calling sequence of func:
334: . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscTruth *flag);
336: . update - latest solution vector
337: . ctx - [optional] timestep context
338: . newdt - the timestep to use for the next step
339: . flag - flag indicating whether the last time step was acceptable
341: Notes:
342: The routine set here will be called by TSPseudoVerifyTimeStep()
343: during the timestepping process.
345: .keywords: timestep, pseudo, set, verify
347: .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
348: @*/
349: PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscTruth*),void* ctx)
350: {
351: PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscTruth *),void *);
356: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);
357: if (f) {
358: (*f)(ts,dt,ctx);
359: }
360: return(0);
361: }
365: /*@
366: TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
367: dt when using the TSPseudoDefaultTimeStep() routine.
369: Collective on TS
371: Input Parameters:
372: + ts - the timestep context
373: - inc - the scaling factor >= 1.0
375: Options Database Key:
376: $ -ts_pseudo_increment <increment>
378: Level: advanced
380: .keywords: timestep, pseudo, set, increment
382: .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
383: @*/
384: PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
385: {
386: PetscErrorCode ierr,(*f)(TS,PetscReal);
391: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);
392: if (f) {
393: (*f)(ts,inc);
394: }
395: return(0);
396: }
400: /*@
401: TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
402: is computed via the formula
403: $ dt = initial_dt*initial_fnorm/current_fnorm
404: rather than the default update,
405: $ dt = current_dt*previous_fnorm/current_fnorm.
407: Collective on TS
409: Input Parameter:
410: . ts - the timestep context
412: Options Database Key:
413: $ -ts_pseudo_increment_dt_from_initial_dt
415: Level: advanced
417: .keywords: timestep, pseudo, set, increment
419: .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
420: @*/
421: PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
422: {
423: PetscErrorCode ierr,(*f)(TS);
428: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);
429: if (f) {
430: (*f)(ts);
431: }
432: return(0);
433: }
438: /*@C
439: TSPseudoSetTimeStep - Sets the user-defined routine to be
440: called at each pseudo-timestep to update the timestep.
442: Collective on TS
444: Input Parameters:
445: + ts - timestep context
446: . dt - function to compute timestep
447: - ctx - [optional] user-defined context for private data
448: required by the function (may be PETSC_NULL)
450: Level: intermediate
452: Calling sequence of func:
453: . func (TS ts,PetscReal *newdt,void *ctx);
455: . newdt - the newly computed timestep
456: . ctx - [optional] timestep context
458: Notes:
459: The routine set here will be called by TSPseudoComputeTimeStep()
460: during the timestepping process.
462: .keywords: timestep, pseudo, set
464: .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
465: @*/
466: PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
467: {
468: PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *);
473: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);
474: if (f) {
475: (*f)(ts,dt,ctx);
476: }
477: return(0);
478: }
480: /* ----------------------------------------------------------------------------- */
486: PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
487: {
488: TS_Pseudo *pseudo;
491: pseudo = (TS_Pseudo*)ts->data;
492: pseudo->verify = dt;
493: pseudo->verifyctx = ctx;
494: return(0);
495: }
501: PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
502: {
503: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
506: pseudo->dt_increment = inc;
507: return(0);
508: }
514: PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
515: {
516: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
519: pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
520: return(0);
521: }
528: PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
529: {
530: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
533: pseudo->dt = dt;
534: pseudo->dtctx = ctx;
535: return(0);
536: }
539: /* ----------------------------------------------------------------------------- */
544: PetscErrorCode TSCreate_Pseudo(TS ts)
545: {
546: TS_Pseudo *pseudo;
550: ts->ops->destroy = TSDestroy_Pseudo;
551: ts->ops->view = TSView_Pseudo;
553: if (ts->problem_type == TS_LINEAR) {
554: SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
555: }
556: if (!ts->A) {
557: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
558: }
560: ts->ops->setup = TSSetUp_Pseudo;
561: ts->ops->step = TSStep_Pseudo;
562: ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
564: /* create the required nonlinear solver context */
565: SNESCreate(ts->comm,&ts->snes);
567: PetscNew(TS_Pseudo,&pseudo);
568: PetscLogObjectMemory(ts,sizeof(TS_Pseudo));
570: PetscMemzero(pseudo,sizeof(TS_Pseudo));
571: ts->data = (void*)pseudo;
573: pseudo->dt_increment = 1.1;
574: pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
575: pseudo->dt = TSPseudoDefaultTimeStep;
577: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
578: "TSPseudoSetVerifyTimeStep_Pseudo",
579: TSPseudoSetVerifyTimeStep_Pseudo);
580: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
581: "TSPseudoSetTimeStepIncrement_Pseudo",
582: TSPseudoSetTimeStepIncrement_Pseudo);
583: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
584: "TSPseudoIncrementDtFromInitialDt_Pseudo",
585: TSPseudoIncrementDtFromInitialDt_Pseudo);
586: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
587: "TSPseudoSetTimeStep_Pseudo",
588: TSPseudoSetTimeStep_Pseudo);
589: return(0);
590: }
595: /*@C
596: TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
597: Use with TSPseudoSetTimeStep().
599: Collective on TS
601: Input Parameters:
602: . ts - the timestep context
603: . dtctx - unused timestep context
605: Output Parameter:
606: . newdt - the timestep to use for the next step
608: Level: advanced
610: .keywords: timestep, pseudo, default
612: .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
613: @*/
614: PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
615: {
616: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
617: PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
621: TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);
622: VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
623: if (pseudo->initial_fnorm == 0.0) {
624: /* first time through so compute initial function norm */
625: pseudo->initial_fnorm = pseudo->fnorm;
626: fnorm_previous = pseudo->fnorm;
627: }
628: if (pseudo->fnorm == 0.0) {
629: *newdt = 1.e12*inc*ts->time_step;
630: } else if (pseudo->increment_dt_from_initial_dt) {
631: *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
632: } else {
633: *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
634: }
635: pseudo->fnorm_previous = pseudo->fnorm;
636: return(0);
637: }