xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 8ba3a721c9105e9240dbd90dceddbf8833cabcd1)
1 #ifndef lint
2 static char vcid[] = "$Id: posindep.c,v 1.14 1997/01/24 04:30:12 curfman Exp curfman $";
3 #endif
4 /*
5        Code for Timestepping with implicit backwards Euler.
6 */
7 #include <math.h>
8 #include "src/ts/tsimpl.h"                /*I   "ts.h"   I*/
9 #include "pinclude/pviewer.h"
10 
11 
12 typedef struct {
13   Vec  update;      /* work vector where new solution is formed */
14   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
15   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
16 
17   /* information used for Pseudo-timestepping */
18 
19   int    (*dt)(TS,double*,void*);              /* compute next timestep, and related context */
20   void   *dtctx;
21   int    (*verify)(TS,Vec,void*,double*,int*); /* verify previous timestep and related context */
22   void   *verifyctx;
23 
24   double initial_fnorm,fnorm;                  /* original and current norm of F(u) */
25   double fnorm_previous;
26 
27   double dt_increment;                  /* scaling that dt is incremented each time-step */
28   int    increment_dt_from_initial_dt;
29 } TS_Pseudo;
30 
31 /* ------------------------------------------------------------------------------*/
32 #undef __FUNC__
33 #define __FUNC__ "TSPseudoDefaultTimeStep"
34 /*@C
35    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
36    Use with TSPseudoSetTimeStep().
37 
38    Input Parameters:
39 .  ts - the timestep context
40 .  dtctx - unused timestep context
41 
42    Output Parameter:
43 .  newdt - the timestep to use for the next step
44 
45 .keywords: timestep, pseudo, default
46 
47 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
48 @*/
49 int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx)
50 {
51   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
52   double    inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
53   int       ierr;
54 
55   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
56   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr);
57   if (pseudo->initial_fnorm == 0.0) {
58     /* first time through so compute initial function norm */
59     pseudo->initial_fnorm = pseudo->fnorm;
60     fnorm_previous        = pseudo->fnorm;
61   }
62   if (pseudo->fnorm == 0.0) {
63     *newdt = 1.e12*inc*ts->time_step;
64   }
65   else if (pseudo->increment_dt_from_initial_dt) {
66     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
67   } else {
68     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
69   }
70   pseudo->fnorm_previous = pseudo->fnorm;
71   return 0;
72 }
73 
74 #undef __FUNC__
75 #define __FUNC__ "TSPseudoSetTimeStep"
76 /*@
77    TSPseudoSetTimeStep - Sets the user-defined routine to be
78    called at each pseudo-timestep to update the timestep.
79 
80    Input Parameters:
81 .  ts - timestep context
82 .  dt - function to compute timestep
83 .  ctx - [optional] user-defined context for private data
84          required by the function (may be PETSC_NULL)
85 
86    Calling sequence of func:
87 .  func (TS ts,double *newdt,void *ctx);
88 
89 .  newdt - the newly computed timestep
90 .  ctx - [optional] timestep context
91 
92    Notes:
93    The routine set here will be called by TSPseudoComputeTimeStep()
94    during the timestepping process.
95 
96 .keywords: timestep, pseudo, set
97 
98 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
99 @*/
100 int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx)
101 {
102   TS_Pseudo *pseudo;
103 
104   PetscValidHeaderSpecific(ts,TS_COOKIE);
105   if (ts->type != TS_PSEUDO) return 0;
106 
107   pseudo          = (TS_Pseudo*) ts->data;
108   pseudo->dt      = dt;
109   pseudo->dtctx   = ctx;
110   return 0;
111 }
112 
113 #undef __FUNC__
114 #define __FUNC__ "TSPseudoComputeTimeStep"
115 /*@
116     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
117     pseudo-timestepping process.
118 
119     Input Parameter:
120 .   ts - timestep context
121 
122     Output Parameter:
123 .   dt - newly computed timestep
124 
125 
126     Notes:
127     The routine to be called here to compute the timestep should be
128     set by calling TSPseudoSetTimeStep().
129 
130 .keywords: timestep, pseudo, compute
131 
132 .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
133 @*/
134 int TSPseudoComputeTimeStep(TS ts,double *dt)
135 {
136   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
137   int       ierr;
138 
139   PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
140   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr);
141   PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
142   return 0;
143 }
144 
145 
146 /* ------------------------------------------------------------------------------*/
147 #undef __FUNC__
148 #define __FUNC__ "TSPseudoDefaultVerifyTimeStep"
149 /*@C
150    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
151 
152    Input Parameters:
153 .  ts - the timestep context
154 .  dtctx - unused timestep context
155 .  update - latest solution vector
156 
157    Output Parameters:
158 .  newdt - the timestep to use for the next step
159 .  flag - flag indicating whether the last time step was acceptable
160 
161    Note:
162    This routine always returns a flag of 1, indicating an acceptable
163    timestep.
164 
165 .keywords: timestep, pseudo, default, verify
166 
167 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
168 @*/
169 int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag)
170 {
171   *flag = 1;
172   return 0;
173 }
174 
175 /*@
176    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
177    last timestep.
178 
179    Input Parameters:
180 .  ts - timestep context
181 .  dt - user-defined function to verify timestep
182 .  ctx - [optional] user-defined context for private data
183          for the timestep verification routine (may be PETSC_NULL)
184 
185    Calling sequence of func:
186 .  func (TS ts,Vec update,void *ctx,double *newdt,int *flag);
187 
188 .  update - latest solution vector
189 .  ctx - [optional] timestep context
190 .  newdt - the timestep to use for the next step
191 .  flag - flag indicating whether the last time step was acceptable
192 
193    Notes:
194    The routine set here will be called by TSPseudoVerifyTimeStep()
195    during the timestepping process.
196 
197 .keywords: timestep, pseudo, set, verify
198 
199 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
200 @*/
201 int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx)
202 {
203   TS_Pseudo *pseudo;
204 
205   PetscValidHeaderSpecific(ts,TS_COOKIE);
206   if (ts->type != TS_PSEUDO) return 0;
207 
208   pseudo              = (TS_Pseudo*) ts->data;
209   pseudo->verify      = dt;
210   pseudo->verifyctx   = ctx;
211   return 0;
212 }
213 
214 #undef __FUNC__
215 #define __FUNC__ "TSPseudoVerifyTimeStep"
216 /*@
217     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
218 
219     Input Parameters:
220 .   ts - timestep context
221 .   update - latest solution vector
222 
223     Output Parameters:
224 .   dt - newly computed timestep (if it had to shrink)
225 .   flag - indicates if current timestep was ok
226 
227     Notes:
228     The routine to be called here to compute the timestep should be
229     set by calling TSPseudoSetVerifyTimeStep().
230 
231 .keywords: timestep, pseudo, verify
232 
233 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
234 @*/
235 int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag)
236 {
237   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
238   int       ierr;
239 
240   if (!pseudo->verify) {*flag = 1; return 0;}
241 
242   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr);
243 
244   return 0;
245 }
246 
247 /* --------------------------------------------------------------------------------*/
248 
249 #undef __FUNC__
250 #define __FUNC__ "TSStep_Pseudo"
251 static int TSStep_Pseudo(TS ts,int *steps,double *time)
252 {
253   Vec       sol = ts->vec_sol;
254   int       ierr,i,max_steps = ts->max_steps,its,ok,lits;
255   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
256   double    current_time_step;
257 
258   *steps = -ts->steps;
259 
260   ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr);
261   for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) {
262     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr);
263     current_time_step = ts->time_step;
264     while (1) {
265       ts->ptime  += current_time_step;
266       ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr);
267       ierr = SNESGetNumberLinearIterations(ts->snes,&lits); CHKERRQ(ierr);
268       ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits;
269       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr);
270       if (ok) break;
271       ts->ptime        -= current_time_step;
272       current_time_step = ts->time_step;
273     }
274     ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr);
275     ts->steps++;
276     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
277   }
278 
279   *steps += ts->steps;
280   *time  = ts->ptime;
281   return 0;
282 }
283 
284 /*------------------------------------------------------------*/
285 #undef __FUNC__
286 #define __FUNC__ "TSDestroy_Pseudo"
287 static int TSDestroy_Pseudo(PetscObject obj )
288 {
289   TS        ts = (TS) obj;
290   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
291   int       ierr;
292 
293   ierr = VecDestroy(pseudo->update); CHKERRQ(ierr);
294   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
295   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
296   if (ts->Ashell)   {ierr = MatDestroy(ts->A); CHKERRQ(ierr);}
297   PetscFree(pseudo);
298   return 0;
299 }
300 
301 
302 /*------------------------------------------------------------*/
303 /*
304     This matrix shell multiply where user provided Shell matrix
305 */
306 
307 #undef __FUNC__
308 #define __FUNC__ "TSPseudoMatMult"
309 int TSPseudoMatMult(Mat mat,Vec x,Vec y)
310 {
311   TS     ts;
312   Scalar mdt,mone = -1.0;
313   int    ierr;
314 
315   MatShellGetContext(mat,(void **)&ts);
316   mdt = 1.0/ts->time_step;
317 
318   /* apply user provided function */
319   ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr);
320   /* shift and scale by 1/dt - F */
321   ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr);
322   return 0;
323 }
324 
325 /*
326     This defines the nonlinear equation that is to be solved with SNES
327 
328               (U^{n+1} - U^{n})/dt - F(U^{n+1})
329 */
330 #undef __FUNC__
331 #define __FUNC__ "TSPseudoFunction"
332 int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
333 {
334   TS     ts = (TS) ctx;
335   Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
336   int    ierr,i,n;
337 
338   /* apply user provided function */
339   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr);
340   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
341   ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr);
342   ierr = VecGetArray(x,&unp1); CHKERRQ(ierr);
343   ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr);
344   ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
345   for ( i=0; i<n; i++ ) {
346     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
347   }
348   ierr = VecRestoreArray(ts->vec_sol,&un);
349   ierr = VecRestoreArray(x,&unp1);
350   ierr = VecRestoreArray(y,&Funp1);
351 
352   return 0;
353 }
354 
355 /*
356    This constructs the Jacobian needed for SNES
357 
358              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
359 */
360 #undef __FUNC__
361 #define __FUNC__ "TSPseudoJacobian"
362 int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
363 {
364   TS      ts = (TS) ctx;
365   int     ierr;
366   Scalar  mone = -1.0, mdt = 1.0/ts->time_step;
367   MatType mtype;
368 
369   /* construct users Jacobian */
370   if (ts->rhsjacobian) {
371     ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr);
372   }
373 
374   /* shift and scale Jacobian, if not a shell matrix */
375   ierr = MatGetType(*AA,&mtype,PETSC_NULL);
376   if (mtype != MATSHELL) {
377     ierr = MatScale(&mone,*AA); CHKERRQ(ierr);
378     ierr = MatShift(&mdt,*AA); CHKERRQ(ierr);
379   }
380   ierr = MatGetType(*BB,&mtype,PETSC_NULL);
381   if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) {
382     ierr = MatScale(&mone,*BB); CHKERRQ(ierr);
383     ierr = MatShift(&mdt,*BB); CHKERRQ(ierr);
384   }
385 
386   return 0;
387 }
388 
389 
390 #undef __FUNC__
391 #define __FUNC__ "TSSetUp_Pseudo"
392 static int TSSetUp_Pseudo(TS ts)
393 {
394   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
395   int         ierr, M, m;
396 
397   ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr);
398   ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr);
399   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
400   if (ts->Ashell) { /* construct new shell matrix */
401     ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr);
402     ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr);
403     ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr);
404     ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr);
405   }
406   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
407   return 0;
408 }
409 /*------------------------------------------------------------*/
410 
411 #undef __FUNC__
412 #define __FUNC__ "TSPseudoDefaultMonitor"
413 int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
414 {
415   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
416 
417   PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm);
418   return 0;
419 }
420 
421 #undef __FUNC__
422 #define __FUNC__ "TSSetFromOptions_Pseudo"
423 static int TSSetFromOptions_Pseudo(TS ts)
424 {
425   int    ierr,flg;
426   double inc;
427 
428   ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr);
429 
430   ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr);
431   if (flg) {
432     ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr);
433   }
434   ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg);  CHKERRQ(ierr);
435   if (flg) {
436     ierr = TSPseudoSetTimeStepIncrement(ts,inc);  CHKERRQ(ierr);
437   }
438   ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr);
439   if (flg) {
440     ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr);
441   }
442   return 0;
443 }
444 
445 #undef __FUNC__
446 #define __FUNC__ "TSPrintHelp_Pseudo"
447 static int TSPrintHelp_Pseudo(TS ts,char *p)
448 {
449   PetscPrintf(ts->comm," Options for TS Pseudo timestepper:\n");
450   PetscPrintf(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p);
451   PetscPrintf(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p);
452   PetscPrintf(ts->comm,"     initial fnorm/current fnorm to determine new timestep\n");
453   PetscPrintf(ts->comm,"     default is current_dt * previous fnorm/current fnorm\n");
454   return 0;
455 }
456 
457 #undef __FUNC__
458 #define __FUNC__ "TSView_Pseudo"
459 static int TSView_Pseudo(PetscObject obj,Viewer viewer)
460 {
461   return 0;
462 }
463 
464 /* ------------------------------------------------------------ */
465 #undef __FUNC__
466 #define __FUNC__ "TSCreate_Pseudo"
467 int TSCreate_Pseudo(TS ts )
468 {
469   TS_Pseudo *pseudo;
470   int       ierr;
471   MatType   mtype;
472 
473   ts->type 	      = TS_PSEUDO;
474   ts->destroy         = TSDestroy_Pseudo;
475   ts->printhelp       = TSPrintHelp_Pseudo;
476   ts->view            = TSView_Pseudo;
477 
478   if (ts->problem_type == TS_LINEAR) {
479     SETERRQ(1,0,"Only for nonlinear problems");
480   }
481   if (!ts->A) {
482     SETERRQ(1,0,"Must set Jacobian");
483   }
484   ierr = MatGetType(ts->A,&mtype,PETSC_NULL);
485   if (mtype == MATSHELL) {
486     ts->Ashell = ts->A;
487   }
488   ts->setup           = TSSetUp_Pseudo;
489   ts->step            = TSStep_Pseudo;
490   ts->setfromoptions  = TSSetFromOptions_Pseudo;
491 
492   /* create the required nonlinear solver context */
493   ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr);
494 
495   pseudo   = PetscNew(TS_Pseudo); CHKPTRQ(pseudo);
496   PetscMemzero(pseudo,sizeof(TS_Pseudo));
497   ts->data = (void *) pseudo;
498 
499   pseudo->dt_increment                 = 1.1;
500   pseudo->increment_dt_from_initial_dt = 0;
501   pseudo->dt                           = TSPseudoDefaultTimeStep;
502   return 0;
503 }
504 
505 
506 #undef __FUNC__
507 #define __FUNC__ "TSPseudoSetTimeStepIncrement"
508 /*@
509     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
510     dt when using the TSPseudoDefaultTimeStep() routine.
511 
512     Input Parameters:
513 .   ts - the timestep context
514 .   inc - the scaling factor >= 1.0
515 
516     Options Database Key:
517 $    -ts_pseudo_increment <increment>
518 
519 .keywords: timestep, pseudo, set, increment
520 
521 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
522 @*/
523 int TSPseudoSetTimeStepIncrement(TS ts,double inc)
524 {
525   TS_Pseudo *pseudo;
526 
527   PetscValidHeaderSpecific(ts,TS_COOKIE);
528   if (ts->type != TS_PSEUDO) return 0;
529 
530   pseudo               = (TS_Pseudo*) ts->data;
531   pseudo->dt_increment = inc;
532   return 0;
533 }
534 
535 #undef __FUNC__
536 #define __FUNC__ "TSPseudoIncrementDtFromInitialDt"
537 /*@
538     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
539     is computed via the formula
540 $         dt = initial_dt*initial_fnorm/current_fnorm
541       rather than the default update,
542 $         dt = current_dt*previous_fnorm/current_fnorm.
543 
544     Input Parameter:
545 .   ts - the timestep context
546 
547     Options Database Key:
548 $    -ts_pseudo_increment_dt_from_initial_dt
549 
550 .keywords: timestep, pseudo, set, increment
551 
552 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
553 @*/
554 int TSPseudoIncrementDtFromInitialDt(TS ts)
555 {
556   TS_Pseudo *pseudo;
557 
558   PetscValidHeaderSpecific(ts,TS_COOKIE);
559   if (ts->type != TS_PSEUDO) return 0;
560 
561   pseudo                               = (TS_Pseudo*) ts->data;
562   pseudo->increment_dt_from_initial_dt = 1;
563   return 0;
564 }
565 
566 
567