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