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