Actual source code: posindep.c
1: /*$Id: posindep.c,v 1.60 2001/09/11 16:34:22 bsmith Exp $*/
2: /*
3: Code for Timestepping with implicit backwards Euler.
4: */
5: #include src/ts/tsimpl.h
7: typedef struct {
8: Vec update; /* work vector where new solution is formed */
9: Vec func; /* work vector where F(t[i],u[i]) is stored */
10: Vec rhs; /* work vector for RHS; vec_sol/dt */
12: /* information used for Pseudo-timestepping */
14: int (*dt)(TS,PetscReal*,void*); /* compute next timestep, and related context */
15: void *dtctx;
16: int (*verify)(TS,Vec,void*,PetscReal*,int*); /* verify previous timestep and related context */
17: void *verifyctx;
19: PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */
20: PetscReal fnorm_previous;
22: PetscReal dt_increment; /* scaling that dt is incremented each time-step */
23: PetscTruth increment_dt_from_initial_dt;
24: } TS_Pseudo;
26: /* ------------------------------------------------------------------------------*/
28: /*@
29: TSPseudoComputeTimeStep - Computes the next timestep for a currently running
30: pseudo-timestepping process.
32: Collective on TS
34: Input Parameter:
35: . ts - timestep context
37: Output Parameter:
38: . dt - newly computed timestep
40: Level: advanced
42: Notes:
43: The routine to be called here to compute the timestep should be
44: set by calling TSPseudoSetTimeStep().
46: .keywords: timestep, pseudo, compute
48: .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
49: @*/
50: int TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
51: {
52: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53: int ierr;
56: PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
57: (*pseudo->dt)(ts,dt,pseudo->dtctx);
58: PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
59: return(0);
60: }
63: /* ------------------------------------------------------------------------------*/
64: /*@C
65: TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
67: Collective on TS
69: Input Parameters:
70: + ts - the timestep context
71: . dtctx - unused timestep context
72: - update - latest solution vector
74: Output Parameters:
75: + newdt - the timestep to use for the next step
76: - flag - flag indicating whether the last time step was acceptable
78: Level: advanced
80: Note:
81: This routine always returns a flag of 1, indicating an acceptable
82: timestep.
84: .keywords: timestep, pseudo, default, verify
86: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
87: @*/
88: int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,int *flag)
89: {
91: *flag = 1;
92: return(0);
93: }
96: /*@
97: TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
99: Collective on TS
101: Input Parameters:
102: + ts - timestep context
103: - update - latest solution vector
105: Output Parameters:
106: + dt - newly computed timestep (if it had to shrink)
107: - flag - indicates if current timestep was ok
109: Level: advanced
111: Notes:
112: The routine to be called here to compute the timestep should be
113: set by calling TSPseudoSetVerifyTimeStep().
115: .keywords: timestep, pseudo, verify
117: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
118: @*/
119: int TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,int *flag)
120: {
121: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
122: int ierr;
125: if (!pseudo->verify) {*flag = 1; return(0);}
127: (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);
129: return(0);
130: }
132: /* --------------------------------------------------------------------------------*/
134: static int TSStep_Pseudo(TS ts,int *steps,PetscReal *ptime)
135: {
136: Vec sol = ts->vec_sol;
137: int ierr,i,max_steps = ts->max_steps,its,ok,lits;
138: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
139: PetscReal current_time_step;
140:
142: *steps = -ts->steps;
144: VecCopy(sol,pseudo->update);
145: for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
146: TSPseudoComputeTimeStep(ts,&ts->time_step);
147: TSMonitor(ts,ts->steps,ts->ptime,sol);
148: current_time_step = ts->time_step;
149: while (PETSC_TRUE) {
150: ts->ptime += current_time_step;
151: SNESSolve(ts->snes,pseudo->update,&its);
152: SNESGetNumberLinearIterations(ts->snes,&lits);
153: ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits;
154: TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);
155: if (ok) break;
156: ts->ptime -= current_time_step;
157: current_time_step = ts->time_step;
158: }
159: VecCopy(pseudo->update,sol);
160: ts->steps++;
161: }
162: TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);
163: VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
164: TSMonitor(ts,ts->steps,ts->ptime,sol);
166: *steps += ts->steps;
167: *ptime = ts->ptime;
168: return(0);
169: }
171: /*------------------------------------------------------------*/
172: static int TSDestroy_Pseudo(TS ts)
173: {
174: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
175: int ierr;
178: if (pseudo->update) {VecDestroy(pseudo->update);}
179: if (pseudo->func) {VecDestroy(pseudo->func);}
180: if (pseudo->rhs) {VecDestroy(pseudo->rhs);}
181: PetscFree(pseudo);
182: return(0);
183: }
186: /*------------------------------------------------------------*/
188: /*
189: This defines the nonlinear equation that is to be solved with SNES
191: (U^{n+1} - U^{n})/dt - F(U^{n+1})
192: */
193: int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
194: {
195: TS ts = (TS) ctx;
196: PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
197: int ierr,i,n;
200: /* apply user provided function */
201: TSComputeRHSFunction(ts,ts->ptime,x,y);
202: /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
203: VecGetArray(ts->vec_sol,&un);
204: VecGetArray(x,&unp1);
205: VecGetArray(y,&Funp1);
206: VecGetLocalSize(x,&n);
207: for (i=0; i<n; i++) {
208: Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
209: }
210: VecRestoreArray(ts->vec_sol,&un);
211: VecRestoreArray(x,&unp1);
212: VecRestoreArray(y,&Funp1);
214: return(0);
215: }
217: /*
218: This constructs the Jacobian needed for SNES
220: J = I/dt - J_{F} where J_{F} is the given Jacobian of F.
221: */
222: int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
223: {
224: TS ts = (TS) ctx;
225: int ierr;
226: PetscScalar mone = -1.0,mdt = 1.0/ts->time_step;
229: /* construct users Jacobian */
230: TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);
232: /* shift and scale Jacobian */
233: MatScale(&mone,*AA);
234: MatShift(&mdt,*AA);
235: if (*BB != *AA && *str != SAME_PRECONDITIONER) {
236: MatScale(&mone,*BB);
237: MatShift(&mdt,*BB);
238: }
240: return(0);
241: }
244: static int TSSetUp_Pseudo(TS ts)
245: {
246: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
247: int ierr;
250: /* SNESSetFromOptions(ts->snes); */
251: VecDuplicate(ts->vec_sol,&pseudo->update);
252: VecDuplicate(ts->vec_sol,&pseudo->func);
253: SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);
254: SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);
255: return(0);
256: }
257: /*------------------------------------------------------------*/
259: int TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
260: {
261: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
262: int ierr;
265: (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %gn",step,ts->time_step,ptime,pseudo->fnorm);
266: return(0);
267: }
269: static int TSSetFromOptions_Pseudo(TS ts)
270: {
271: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
272: int ierr;
273: PetscTruth flg;
277: PetscOptionsHead("Pseudo-timestepping options");
278: PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);
279: if (flg) {
280: TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);
281: }
282: PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);
283: if (flg) {
284: TSPseudoIncrementDtFromInitialDt(ts);
285: }
286: PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);
287: PetscOptionsTail();
289: return(0);
290: }
292: static int TSView_Pseudo(TS ts,PetscViewer viewer)
293: {
295: return(0);
296: }
298: /* ----------------------------------------------------------------------------- */
299: /*@
300: TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
301: last timestep.
303: Collective on TS
305: Input Parameters:
306: + ts - timestep context
307: . dt - user-defined function to verify timestep
308: - ctx - [optional] user-defined context for private data
309: for the timestep verification routine (may be PETSC_NULL)
311: Level: advanced
313: Calling sequence of func:
314: . func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag);
316: . update - latest solution vector
317: . ctx - [optional] timestep context
318: . newdt - the timestep to use for the next step
319: . flag - flag indicating whether the last time step was acceptable
321: Notes:
322: The routine set here will be called by TSPseudoVerifyTimeStep()
323: during the timestepping process.
325: .keywords: timestep, pseudo, set, verify
327: .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
328: @*/
329: int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx)
330: {
331: int ierr,(*f)(TS,int (*)(TS,Vec,void*,PetscReal *,int *),void *);
336: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);
337: if (f) {
338: (*f)(ts,dt,ctx);
339: }
340: return(0);
341: }
343: /*@
344: TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
345: dt when using the TSPseudoDefaultTimeStep() routine.
347: Collective on TS
349: Input Parameters:
350: + ts - the timestep context
351: - inc - the scaling factor >= 1.0
353: Options Database Key:
354: $ -ts_pseudo_increment <increment>
356: Level: advanced
358: .keywords: timestep, pseudo, set, increment
360: .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
361: @*/
362: int TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
363: {
364: int ierr,(*f)(TS,PetscReal);
369: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);
370: if (f) {
371: (*f)(ts,inc);
372: }
373: return(0);
374: }
376: /*@
377: TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
378: is computed via the formula
379: $ dt = initial_dt*initial_fnorm/current_fnorm
380: rather than the default update,
381: $ dt = current_dt*previous_fnorm/current_fnorm.
383: Collective on TS
385: Input Parameter:
386: . ts - the timestep context
388: Options Database Key:
389: $ -ts_pseudo_increment_dt_from_initial_dt
391: Level: advanced
393: .keywords: timestep, pseudo, set, increment
395: .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
396: @*/
397: int TSPseudoIncrementDtFromInitialDt(TS ts)
398: {
399: int ierr,(*f)(TS);
404: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);
405: if (f) {
406: (*f)(ts);
407: }
408: return(0);
409: }
412: /*@
413: TSPseudoSetTimeStep - Sets the user-defined routine to be
414: called at each pseudo-timestep to update the timestep.
416: Collective on TS
418: Input Parameters:
419: + ts - timestep context
420: . dt - function to compute timestep
421: - ctx - [optional] user-defined context for private data
422: required by the function (may be PETSC_NULL)
424: Level: intermediate
426: Calling sequence of func:
427: . func (TS ts,PetscReal *newdt,void *ctx);
429: . newdt - the newly computed timestep
430: . ctx - [optional] timestep context
432: Notes:
433: The routine set here will be called by TSPseudoComputeTimeStep()
434: during the timestepping process.
436: .keywords: timestep, pseudo, set
438: .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
439: @*/
440: int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx)
441: {
442: int ierr,(*f)(TS,int (*)(TS,PetscReal *,void *),void *);
447: PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);
448: if (f) {
449: (*f)(ts,dt,ctx);
450: }
451: return(0);
452: }
454: /* ----------------------------------------------------------------------------- */
456: EXTERN_C_BEGIN
457: int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx)
458: {
459: TS_Pseudo *pseudo;
462: pseudo = (TS_Pseudo*)ts->data;
463: pseudo->verify = dt;
464: pseudo->verifyctx = ctx;
465: return(0);
466: }
467: EXTERN_C_END
469: EXTERN_C_BEGIN
470: int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
471: {
472: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
475: pseudo->dt_increment = inc;
476: return(0);
477: }
478: EXTERN_C_END
480: EXTERN_C_BEGIN
481: int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
482: {
483: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
486: pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
487: return(0);
488: }
489: EXTERN_C_END
491: EXTERN_C_BEGIN
492: int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx)
493: {
494: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
497: pseudo->dt = dt;
498: pseudo->dtctx = ctx;
499: return(0);
500: }
501: EXTERN_C_END
503: /* ----------------------------------------------------------------------------- */
505: EXTERN_C_BEGIN
506: int TSCreate_Pseudo(TS ts)
507: {
508: TS_Pseudo *pseudo;
509: int ierr;
512: ts->ops->destroy = TSDestroy_Pseudo;
513: ts->ops->view = TSView_Pseudo;
515: if (ts->problem_type == TS_LINEAR) {
516: SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
517: }
518: if (!ts->A) {
519: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
520: }
522: ts->ops->setup = TSSetUp_Pseudo;
523: ts->ops->step = TSStep_Pseudo;
524: ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
526: /* create the required nonlinear solver context */
527: SNESCreate(ts->comm,&ts->snes);
529: PetscNew(TS_Pseudo,&pseudo);
530: PetscLogObjectMemory(ts,sizeof(TS_Pseudo));
532: ierr = PetscMemzero(pseudo,sizeof(TS_Pseudo));
533: ts->data = (void*)pseudo;
535: pseudo->dt_increment = 1.1;
536: pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
537: pseudo->dt = TSPseudoDefaultTimeStep;
539: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
540: "TSPseudoSetVerifyTimeStep_Pseudo",
541: TSPseudoSetVerifyTimeStep_Pseudo);
542: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
543: "TSPseudoSetTimeStepIncrement_Pseudo",
544: TSPseudoSetTimeStepIncrement_Pseudo);
545: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
546: "TSPseudoIncrementDtFromInitialDt_Pseudo",
547: TSPseudoIncrementDtFromInitialDt_Pseudo);
548: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
549: "TSPseudoSetTimeStep_Pseudo",
550: TSPseudoSetTimeStep_Pseudo);
551: return(0);
552: }
553: EXTERN_C_END
555: /*@C
556: TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
557: Use with TSPseudoSetTimeStep().
559: Collective on TS
561: Input Parameters:
562: . ts - the timestep context
563: . dtctx - unused timestep context
565: Output Parameter:
566: . newdt - the timestep to use for the next step
568: Level: advanced
570: .keywords: timestep, pseudo, default
572: .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
573: @*/
574: int TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
575: {
576: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
577: PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
578: int ierr;
581: TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);
582: VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
583: if (pseudo->initial_fnorm == 0.0) {
584: /* first time through so compute initial function norm */
585: pseudo->initial_fnorm = pseudo->fnorm;
586: fnorm_previous = pseudo->fnorm;
587: }
588: if (pseudo->fnorm == 0.0) {
589: *newdt = 1.e12*inc*ts->time_step;
590: } else if (pseudo->increment_dt_from_initial_dt) {
591: *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
592: } else {
593: *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
594: }
595: pseudo->fnorm_previous = pseudo->fnorm;
596: return(0);
597: }