xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 76be9ce4a233aaa47cda2bc7f5a27cd7faabecaa)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*76be9ce4SBarry Smith static char vcid[] = "$Id: posindep.c,v 1.22 1997/12/01 01:56:10 bsmith Exp bsmith $";
32d3f70b5SBarry Smith #endif
42d3f70b5SBarry Smith /*
5fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
62d3f70b5SBarry Smith */
72d3f70b5SBarry Smith #include <math.h>
82d3f70b5SBarry Smith #include "src/ts/tsimpl.h"                /*I   "ts.h"   I*/
92d3f70b5SBarry Smith #include "pinclude/pviewer.h"
102d3f70b5SBarry Smith 
112d3f70b5SBarry Smith 
122d3f70b5SBarry Smith typedef struct {
132d3f70b5SBarry Smith   Vec  update;      /* work vector where new solution is formed */
142d3f70b5SBarry Smith   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
152d3f70b5SBarry Smith   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
162d3f70b5SBarry Smith 
172d3f70b5SBarry Smith   /* information used for Pseudo-timestepping */
182d3f70b5SBarry Smith 
197bf11e45SBarry Smith   int    (*dt)(TS,double*,void*);              /* compute next timestep, and related context */
202d3f70b5SBarry Smith   void   *dtctx;
217bf11e45SBarry Smith   int    (*verify)(TS,Vec,void*,double*,int*); /* verify previous timestep and related context */
227bf11e45SBarry Smith   void   *verifyctx;
232d3f70b5SBarry Smith 
247bf11e45SBarry Smith   double initial_fnorm,fnorm;                  /* original and current norm of F(u) */
25ca4b7087SBarry Smith   double fnorm_previous;
2628aa8177SBarry Smith 
2728aa8177SBarry Smith   double dt_increment;                  /* scaling that dt is incremented each time-step */
28ca4b7087SBarry Smith   int    increment_dt_from_initial_dt;
297bf11e45SBarry Smith } TS_Pseudo;
302d3f70b5SBarry Smith 
312d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
3258562591SBarry Smith #undef __FUNC__
3358562591SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep"
342d3f70b5SBarry Smith /*@C
35fb4a63b6SLois Curfman McInnes    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
36fb4a63b6SLois Curfman McInnes    Use with TSPseudoSetTimeStep().
372d3f70b5SBarry Smith 
382d3f70b5SBarry Smith    Input Parameters:
392d3f70b5SBarry Smith .  ts - the timestep context
407bf11e45SBarry Smith .  dtctx - unused timestep context
412d3f70b5SBarry Smith 
422d3f70b5SBarry Smith    Output Parameter:
432d3f70b5SBarry Smith .  newdt - the timestep to use for the next step
442d3f70b5SBarry Smith 
45fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, default
46fb4a63b6SLois Curfman McInnes 
47564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
482d3f70b5SBarry Smith @*/
497bf11e45SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx)
502d3f70b5SBarry Smith {
517bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
52ca4b7087SBarry Smith   double    inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
532d3f70b5SBarry Smith   int       ierr;
542d3f70b5SBarry Smith 
553a40ed3dSBarry Smith   PetscFunctionBegin;
567bf11e45SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
577bf11e45SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr);
587bf11e45SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
592d3f70b5SBarry Smith     /* first time through so compute initial function norm */
607bf11e45SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
61ca4b7087SBarry Smith     fnorm_previous        = pseudo->fnorm;
622d3f70b5SBarry Smith   }
63ca4b7087SBarry Smith   if (pseudo->fnorm == 0.0) {
64ca4b7087SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
65ca4b7087SBarry Smith   }
66ca4b7087SBarry Smith   else if (pseudo->increment_dt_from_initial_dt) {
67ca4b7087SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
68ca4b7087SBarry Smith   } else {
69ca4b7087SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
70ca4b7087SBarry Smith   }
71ca4b7087SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
732d3f70b5SBarry Smith }
742d3f70b5SBarry Smith 
7558562591SBarry Smith #undef __FUNC__
7658562591SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep"
772d3f70b5SBarry Smith /*@
78fb4a63b6SLois Curfman McInnes    TSPseudoSetTimeStep - Sets the user-defined routine to be
79fb4a63b6SLois Curfman McInnes    called at each pseudo-timestep to update the timestep.
802d3f70b5SBarry Smith 
812d3f70b5SBarry Smith    Input Parameters:
822d3f70b5SBarry Smith .  ts - timestep context
832d3f70b5SBarry Smith .  dt - function to compute timestep
84564e8f4eSLois Curfman McInnes .  ctx - [optional] user-defined context for private data
85564e8f4eSLois Curfman McInnes          required by the function (may be PETSC_NULL)
86564e8f4eSLois Curfman McInnes 
87564e8f4eSLois Curfman McInnes    Calling sequence of func:
88564e8f4eSLois Curfman McInnes .  func (TS ts,double *newdt,void *ctx);
89564e8f4eSLois Curfman McInnes 
90564e8f4eSLois Curfman McInnes .  newdt - the newly computed timestep
91564e8f4eSLois Curfman McInnes .  ctx - [optional] timestep context
92564e8f4eSLois Curfman McInnes 
93564e8f4eSLois Curfman McInnes    Notes:
94564e8f4eSLois Curfman McInnes    The routine set here will be called by TSPseudoComputeTimeStep()
95564e8f4eSLois Curfman McInnes    during the timestepping process.
962d3f70b5SBarry Smith 
97fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, set
98fb4a63b6SLois Curfman McInnes 
99564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
1002d3f70b5SBarry Smith @*/
1017bf11e45SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx)
1022d3f70b5SBarry Smith {
1037bf11e45SBarry Smith   TS_Pseudo *pseudo;
1042d3f70b5SBarry Smith 
1053a40ed3dSBarry Smith   PetscFunctionBegin;
1062d3f70b5SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1073a40ed3dSBarry Smith   if (ts->type != TS_PSEUDO) PetscFunctionReturn(0);
1087bf11e45SBarry Smith 
1097bf11e45SBarry Smith   pseudo          = (TS_Pseudo*) ts->data;
1107bf11e45SBarry Smith   pseudo->dt      = dt;
1117bf11e45SBarry Smith   pseudo->dtctx   = ctx;
1123a40ed3dSBarry Smith   PetscFunctionReturn(0);
1132d3f70b5SBarry Smith }
1142d3f70b5SBarry Smith 
11558562591SBarry Smith #undef __FUNC__
11658562591SBarry Smith #define __FUNC__ "TSPseudoComputeTimeStep"
1177bf11e45SBarry Smith /*@
1187bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
119564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
1202d3f70b5SBarry Smith 
1217bf11e45SBarry Smith     Input Parameter:
1227bf11e45SBarry Smith .   ts - timestep context
1237bf11e45SBarry Smith 
1247bf11e45SBarry Smith     Output Parameter:
125fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
126fb4a63b6SLois Curfman McInnes 
127564e8f4eSLois Curfman McInnes 
128564e8f4eSLois Curfman McInnes     Notes:
129564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
130564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
131564e8f4eSLois Curfman McInnes 
132fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
133564e8f4eSLois Curfman McInnes 
134564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
1357bf11e45SBarry Smith @*/
1367bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt)
1377bf11e45SBarry Smith {
1387bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
1397bf11e45SBarry Smith   int       ierr;
1407bf11e45SBarry Smith 
1413a40ed3dSBarry Smith   PetscFunctionBegin;
1427bf11e45SBarry Smith   PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
1437bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr);
1447bf11e45SBarry Smith   PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
1453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1467bf11e45SBarry Smith }
1477bf11e45SBarry Smith 
1487bf11e45SBarry Smith 
1497bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
15058562591SBarry Smith #undef __FUNC__
15158562591SBarry Smith #define __FUNC__ "TSPseudoDefaultVerifyTimeStep"
1527bf11e45SBarry Smith /*@C
153639f9d9dSBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
1547bf11e45SBarry Smith 
1557bf11e45SBarry Smith    Input Parameters:
1567bf11e45SBarry Smith .  ts - the timestep context
1577bf11e45SBarry Smith .  dtctx - unused timestep context
158564e8f4eSLois Curfman McInnes .  update - latest solution vector
1597bf11e45SBarry Smith 
160564e8f4eSLois Curfman McInnes    Output Parameters:
1617bf11e45SBarry Smith .  newdt - the timestep to use for the next step
162564e8f4eSLois Curfman McInnes .  flag - flag indicating whether the last time step was acceptable
1637bf11e45SBarry Smith 
164564e8f4eSLois Curfman McInnes    Note:
165564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
166564e8f4eSLois Curfman McInnes    timestep.
167564e8f4eSLois Curfman McInnes 
168564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
169564e8f4eSLois Curfman McInnes 
170564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
1717bf11e45SBarry Smith @*/
1727bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag)
1737bf11e45SBarry Smith {
1743a40ed3dSBarry Smith   PetscFunctionBegin;
1757bf11e45SBarry Smith   *flag = 1;
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
1777bf11e45SBarry Smith }
1787bf11e45SBarry Smith 
1798304798fSSatish Balay #undef __FUNC__
1808304798fSSatish Balay #define __FUNC__ "TSPseudoSetVerifyTimeStep"
1817bf11e45SBarry Smith /*@
182564e8f4eSLois Curfman McInnes    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
183564e8f4eSLois Curfman McInnes    last timestep.
1847bf11e45SBarry Smith 
1857bf11e45SBarry Smith    Input Parameters:
1867bf11e45SBarry Smith .  ts - timestep context
187564e8f4eSLois Curfman McInnes .  dt - user-defined function to verify timestep
188564e8f4eSLois Curfman McInnes .  ctx - [optional] user-defined context for private data
189564e8f4eSLois Curfman McInnes          for the timestep verification routine (may be PETSC_NULL)
1907bf11e45SBarry Smith 
191564e8f4eSLois Curfman McInnes    Calling sequence of func:
192564e8f4eSLois Curfman McInnes .  func (TS ts,Vec update,void *ctx,double *newdt,int *flag);
193564e8f4eSLois Curfman McInnes 
194564e8f4eSLois Curfman McInnes .  update - latest solution vector
195564e8f4eSLois Curfman McInnes .  ctx - [optional] timestep context
196564e8f4eSLois Curfman McInnes .  newdt - the timestep to use for the next step
197564e8f4eSLois Curfman McInnes .  flag - flag indicating whether the last time step was acceptable
198564e8f4eSLois Curfman McInnes 
199564e8f4eSLois Curfman McInnes    Notes:
200564e8f4eSLois Curfman McInnes    The routine set here will be called by TSPseudoVerifyTimeStep()
201564e8f4eSLois Curfman McInnes    during the timestepping process.
202564e8f4eSLois Curfman McInnes 
203564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, verify
204564e8f4eSLois Curfman McInnes 
205564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
2067bf11e45SBarry Smith @*/
2077bf11e45SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx)
2087bf11e45SBarry Smith {
2097bf11e45SBarry Smith   TS_Pseudo *pseudo;
2107bf11e45SBarry Smith 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
2127bf11e45SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
2133a40ed3dSBarry Smith   if (ts->type != TS_PSEUDO) PetscFunctionReturn(0);
2147bf11e45SBarry Smith 
2157bf11e45SBarry Smith   pseudo              = (TS_Pseudo*) ts->data;
2167bf11e45SBarry Smith   pseudo->verify      = dt;
2177bf11e45SBarry Smith   pseudo->verifyctx   = ctx;
2183a40ed3dSBarry Smith   PetscFunctionReturn(0);
2197bf11e45SBarry Smith }
2207bf11e45SBarry Smith 
22158562591SBarry Smith #undef __FUNC__
22258562591SBarry Smith #define __FUNC__ "TSPseudoVerifyTimeStep"
2237bf11e45SBarry Smith /*@
224564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
2257bf11e45SBarry Smith 
226fb4a63b6SLois Curfman McInnes     Input Parameters:
2277bf11e45SBarry Smith .   ts - timestep context
228564e8f4eSLois Curfman McInnes .   update - latest solution vector
2297bf11e45SBarry Smith 
230fb4a63b6SLois Curfman McInnes     Output Parameters:
231fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep (if it had to shrink)
2327bf11e45SBarry Smith .   flag - indicates if current timestep was ok
2337bf11e45SBarry Smith 
234564e8f4eSLois Curfman McInnes     Notes:
235564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
236564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
237564e8f4eSLois Curfman McInnes 
238564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
239564e8f4eSLois Curfman McInnes 
240564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
2417bf11e45SBarry Smith @*/
2427bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag)
2437bf11e45SBarry Smith {
2447bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
2457bf11e45SBarry Smith   int       ierr;
2467bf11e45SBarry Smith 
2473a40ed3dSBarry Smith   PetscFunctionBegin;
2483a40ed3dSBarry Smith   if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);}
2497bf11e45SBarry Smith 
2507bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr);
2517bf11e45SBarry Smith 
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
2537bf11e45SBarry Smith }
2547bf11e45SBarry Smith 
2557bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
2567bf11e45SBarry Smith 
25758562591SBarry Smith #undef __FUNC__
25858562591SBarry Smith #define __FUNC__ "TSStep_Pseudo"
2597bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time)
2602d3f70b5SBarry Smith {
2612d3f70b5SBarry Smith   Vec       sol = ts->vec_sol;
262f2267985SLois Curfman McInnes   int       ierr,i,max_steps = ts->max_steps,its,ok,lits;
2637bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
2647bf11e45SBarry Smith   double    current_time_step;
2652d3f70b5SBarry Smith 
2663a40ed3dSBarry Smith   PetscFunctionBegin;
2672d3f70b5SBarry Smith   *steps = -ts->steps;
2682d3f70b5SBarry Smith 
2697bf11e45SBarry Smith   ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr);
2702d3f70b5SBarry Smith   for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) {
2717bf11e45SBarry Smith     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr);
2727bf11e45SBarry Smith     current_time_step = ts->time_step;
2737bf11e45SBarry Smith     while (1) {
2747bf11e45SBarry Smith       ts->ptime  += current_time_step;
2757bf11e45SBarry Smith       ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr);
276f2267985SLois Curfman McInnes       ierr = SNESGetNumberLinearIterations(ts->snes,&lits); CHKERRQ(ierr);
277f2267985SLois Curfman McInnes       ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits;
2787bf11e45SBarry Smith       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr);
2797bf11e45SBarry Smith       if (ok) break;
2807bf11e45SBarry Smith       ts->ptime        -= current_time_step;
2817bf11e45SBarry Smith       current_time_step = ts->time_step;
2827bf11e45SBarry Smith     }
2837bf11e45SBarry Smith     ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr);
2842d3f70b5SBarry Smith     ts->steps++;
2852d3f70b5SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
2862d3f70b5SBarry Smith   }
2872d3f70b5SBarry Smith 
2882d3f70b5SBarry Smith   *steps += ts->steps;
2892d3f70b5SBarry Smith   *time  = ts->ptime;
2903a40ed3dSBarry Smith   PetscFunctionReturn(0);
2912d3f70b5SBarry Smith }
2922d3f70b5SBarry Smith 
2932d3f70b5SBarry Smith /*------------------------------------------------------------*/
29458562591SBarry Smith #undef __FUNC__
295d4bb536fSBarry Smith #define __FUNC__ "TSDestroy_Pseudo"
2967bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj )
2972d3f70b5SBarry Smith {
2982d3f70b5SBarry Smith   TS        ts = (TS) obj;
2997bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
3002d3f70b5SBarry Smith   int       ierr;
3012d3f70b5SBarry Smith 
3023a40ed3dSBarry Smith   PetscFunctionBegin;
3037bf11e45SBarry Smith   ierr = VecDestroy(pseudo->update); CHKERRQ(ierr);
3047bf11e45SBarry Smith   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
3057bf11e45SBarry Smith   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
3062d3f70b5SBarry Smith   if (ts->Ashell)   {ierr = MatDestroy(ts->A); CHKERRQ(ierr);}
3077bf11e45SBarry Smith   PetscFree(pseudo);
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
3092d3f70b5SBarry Smith }
3102d3f70b5SBarry Smith 
3112d3f70b5SBarry Smith 
3122d3f70b5SBarry Smith /*------------------------------------------------------------*/
3132d3f70b5SBarry Smith /*
3142d3f70b5SBarry Smith     This matrix shell multiply where user provided Shell matrix
3152d3f70b5SBarry Smith */
3162d3f70b5SBarry Smith 
31758562591SBarry Smith #undef __FUNC__
31858562591SBarry Smith #define __FUNC__ "TSPseudoMatMult"
3197bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y)
3202d3f70b5SBarry Smith {
3212d3f70b5SBarry Smith   TS     ts;
3222d3f70b5SBarry Smith   Scalar mdt,mone = -1.0;
3232d3f70b5SBarry Smith   int    ierr;
3242d3f70b5SBarry Smith 
3253a40ed3dSBarry Smith   PetscFunctionBegin;
3263a40ed3dSBarry Smith   ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr);
3272d3f70b5SBarry Smith   mdt = 1.0/ts->time_step;
3282d3f70b5SBarry Smith 
3292d3f70b5SBarry Smith   /* apply user provided function */
3302d3f70b5SBarry Smith   ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr);
3312d3f70b5SBarry Smith   /* shift and scale by 1/dt - F */
3322d3f70b5SBarry Smith   ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr);
3333a40ed3dSBarry Smith   PetscFunctionReturn(0);
3342d3f70b5SBarry Smith }
3352d3f70b5SBarry Smith 
3362d3f70b5SBarry Smith /*
3372d3f70b5SBarry Smith     This defines the nonlinear equation that is to be solved with SNES
3382d3f70b5SBarry Smith 
3392d3f70b5SBarry Smith               (U^{n+1} - U^{n})/dt - F(U^{n+1})
3402d3f70b5SBarry Smith */
34158562591SBarry Smith #undef __FUNC__
34258562591SBarry Smith #define __FUNC__ "TSPseudoFunction"
3437bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
3442d3f70b5SBarry Smith {
3452d3f70b5SBarry Smith   TS     ts = (TS) ctx;
3462d3f70b5SBarry Smith   Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
3472d3f70b5SBarry Smith   int    ierr,i,n;
3482d3f70b5SBarry Smith 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
3502d3f70b5SBarry Smith   /* apply user provided function */
3512d3f70b5SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr);
3527bf11e45SBarry Smith   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
3532d3f70b5SBarry Smith   ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr);
3542d3f70b5SBarry Smith   ierr = VecGetArray(x,&unp1); CHKERRQ(ierr);
3552d3f70b5SBarry Smith   ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr);
3562d3f70b5SBarry Smith   ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
3572d3f70b5SBarry Smith   for ( i=0; i<n; i++ ) {
3582d3f70b5SBarry Smith     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
3592d3f70b5SBarry Smith   }
3602d3f70b5SBarry Smith   ierr = VecRestoreArray(ts->vec_sol,&un);
3612d3f70b5SBarry Smith   ierr = VecRestoreArray(x,&unp1);
3622d3f70b5SBarry Smith   ierr = VecRestoreArray(y,&Funp1);
3632d3f70b5SBarry Smith 
3643a40ed3dSBarry Smith   PetscFunctionReturn(0);
3652d3f70b5SBarry Smith }
3662d3f70b5SBarry Smith 
3672d3f70b5SBarry Smith /*
3682d3f70b5SBarry Smith    This constructs the Jacobian needed for SNES
3692d3f70b5SBarry Smith 
3702d3f70b5SBarry Smith              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
3712d3f70b5SBarry Smith */
37258562591SBarry Smith #undef __FUNC__
37358562591SBarry Smith #define __FUNC__ "TSPseudoJacobian"
3747bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
3752d3f70b5SBarry Smith {
3762d3f70b5SBarry Smith   TS      ts = (TS) ctx;
3772d3f70b5SBarry Smith   int     ierr;
3782d3f70b5SBarry Smith   Scalar  mone = -1.0, mdt = 1.0/ts->time_step;
3792d3f70b5SBarry Smith   MatType mtype;
3802d3f70b5SBarry Smith 
3813a40ed3dSBarry Smith   PetscFunctionBegin;
3822d3f70b5SBarry Smith   /* construct users Jacobian */
3832d3f70b5SBarry Smith   if (ts->rhsjacobian) {
3842d3f70b5SBarry Smith     ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr);
3852d3f70b5SBarry Smith   }
3862d3f70b5SBarry Smith 
3872d3f70b5SBarry Smith   /* shift and scale Jacobian, if not a shell matrix */
3882d3f70b5SBarry Smith   ierr = MatGetType(*AA,&mtype,PETSC_NULL);
3892d3f70b5SBarry Smith   if (mtype != MATSHELL) {
3902d3f70b5SBarry Smith     ierr = MatScale(&mone,*AA); CHKERRQ(ierr);
3912d3f70b5SBarry Smith     ierr = MatShift(&mdt,*AA); CHKERRQ(ierr);
3922d3f70b5SBarry Smith   }
3932d3f70b5SBarry Smith   ierr = MatGetType(*BB,&mtype,PETSC_NULL);
3942d3f70b5SBarry Smith   if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) {
3952d3f70b5SBarry Smith     ierr = MatScale(&mone,*BB); CHKERRQ(ierr);
3962d3f70b5SBarry Smith     ierr = MatShift(&mdt,*BB); CHKERRQ(ierr);
3972d3f70b5SBarry Smith   }
3982d3f70b5SBarry Smith 
3993a40ed3dSBarry Smith   PetscFunctionReturn(0);
4002d3f70b5SBarry Smith }
4012d3f70b5SBarry Smith 
4022d3f70b5SBarry Smith 
40358562591SBarry Smith #undef __FUNC__
40458562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo"
4057bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts)
4062d3f70b5SBarry Smith {
4077bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
4082d3f70b5SBarry Smith   int       ierr, M, m;
4092d3f70b5SBarry Smith 
4103a40ed3dSBarry Smith   PetscFunctionBegin;
4117bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr);
4127bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr);
4137bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
4142d3f70b5SBarry Smith   if (ts->Ashell) { /* construct new shell matrix */
4152d3f70b5SBarry Smith     ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr);
4162d3f70b5SBarry Smith     ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr);
4172d3f70b5SBarry Smith     ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr);
4187bf11e45SBarry Smith     ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr);
4192d3f70b5SBarry Smith   }
4207bf11e45SBarry Smith   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
4213a40ed3dSBarry Smith   PetscFunctionReturn(0);
4222d3f70b5SBarry Smith }
4232d3f70b5SBarry Smith /*------------------------------------------------------------*/
4242d3f70b5SBarry Smith 
42558562591SBarry Smith #undef __FUNC__
426d4bb536fSBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor"
4272d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
4282d3f70b5SBarry Smith {
4297bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
4302d3f70b5SBarry Smith 
4313a40ed3dSBarry Smith   PetscFunctionBegin;
432*76be9ce4SBarry Smith   (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm);
4333a40ed3dSBarry Smith   PetscFunctionReturn(0);
4342d3f70b5SBarry Smith }
4352d3f70b5SBarry Smith 
43658562591SBarry Smith #undef __FUNC__
43758562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo"
4387bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts)
4392d3f70b5SBarry Smith {
4402d3f70b5SBarry Smith   int    ierr,flg;
44128aa8177SBarry Smith   double inc;
4422d3f70b5SBarry Smith 
4433a40ed3dSBarry Smith   PetscFunctionBegin;
4442d3f70b5SBarry Smith   ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr);
4452d3f70b5SBarry Smith 
4462d3f70b5SBarry Smith   ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr);
4472d3f70b5SBarry Smith   if (flg) {
44828aa8177SBarry Smith     ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr);
44928aa8177SBarry Smith   }
45028aa8177SBarry Smith   ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg);  CHKERRQ(ierr);
45128aa8177SBarry Smith   if (flg) {
45228aa8177SBarry Smith     ierr = TSPseudoSetTimeStepIncrement(ts,inc);  CHKERRQ(ierr);
4532d3f70b5SBarry Smith   }
454ca4b7087SBarry Smith   ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr);
455ca4b7087SBarry Smith   if (flg) {
456ca4b7087SBarry Smith     ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr);
457ca4b7087SBarry Smith   }
4583a40ed3dSBarry Smith   PetscFunctionReturn(0);
4592d3f70b5SBarry Smith }
4602d3f70b5SBarry Smith 
46158562591SBarry Smith #undef __FUNC__
462d4bb536fSBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo"
463ca4b7087SBarry Smith static int TSPrintHelp_Pseudo(TS ts,char *p)
4642d3f70b5SBarry Smith {
4653a40ed3dSBarry Smith   PetscFunctionBegin;
466*76be9ce4SBarry Smith   (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n");
467*76be9ce4SBarry Smith   (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p);
468*76be9ce4SBarry Smith   (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p);
469*76be9ce4SBarry Smith   (*PetscHelpPrintf)(ts->comm,"     initial fnorm/current fnorm to determine new timestep\n");
470*76be9ce4SBarry Smith   (*PetscHelpPrintf)(ts->comm,"     default is current_dt * previous fnorm/current fnorm\n");
4713a40ed3dSBarry Smith   PetscFunctionReturn(0);
4722d3f70b5SBarry Smith }
4732d3f70b5SBarry Smith 
47458562591SBarry Smith #undef __FUNC__
475d4bb536fSBarry Smith #define __FUNC__ "TSView_Pseudo"
4767bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer)
4772d3f70b5SBarry Smith {
4783a40ed3dSBarry Smith   PetscFunctionBegin;
4793a40ed3dSBarry Smith   PetscFunctionReturn(0);
4802d3f70b5SBarry Smith }
4812d3f70b5SBarry Smith 
4822d3f70b5SBarry Smith /* ------------------------------------------------------------ */
48358562591SBarry Smith #undef __FUNC__
48458562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo"
4857bf11e45SBarry Smith int TSCreate_Pseudo(TS ts )
4862d3f70b5SBarry Smith {
4877bf11e45SBarry Smith   TS_Pseudo *pseudo;
4882d3f70b5SBarry Smith   int       ierr;
4892d3f70b5SBarry Smith   MatType   mtype;
4902d3f70b5SBarry Smith 
4913a40ed3dSBarry Smith   PetscFunctionBegin;
4927bf11e45SBarry Smith   ts->type 	      = TS_PSEUDO;
4937bf11e45SBarry Smith   ts->destroy         = TSDestroy_Pseudo;
4947bf11e45SBarry Smith   ts->printhelp       = TSPrintHelp_Pseudo;
4957bf11e45SBarry Smith   ts->view            = TSView_Pseudo;
4962d3f70b5SBarry Smith 
4972d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
498a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems");
4992d3f70b5SBarry Smith   }
5002d3f70b5SBarry Smith   if (!ts->A) {
501a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian");
5022d3f70b5SBarry Smith   }
5032d3f70b5SBarry Smith   ierr = MatGetType(ts->A,&mtype,PETSC_NULL);
5042d3f70b5SBarry Smith   if (mtype == MATSHELL) {
5052d3f70b5SBarry Smith     ts->Ashell = ts->A;
5062d3f70b5SBarry Smith   }
5077bf11e45SBarry Smith   ts->setup           = TSSetUp_Pseudo;
5087bf11e45SBarry Smith   ts->step            = TSStep_Pseudo;
5097bf11e45SBarry Smith   ts->setfromoptions  = TSSetFromOptions_Pseudo;
5107bf11e45SBarry Smith 
5117bf11e45SBarry Smith   /* create the required nonlinear solver context */
5122d3f70b5SBarry Smith   ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr);
5132d3f70b5SBarry Smith 
5147bf11e45SBarry Smith   pseudo   = PetscNew(TS_Pseudo); CHKPTRQ(pseudo);
515eed86810SBarry Smith   PLogObjectMemory(ts,sizeof(TS_Pseudo));
516eed86810SBarry Smith 
5177bf11e45SBarry Smith   PetscMemzero(pseudo,sizeof(TS_Pseudo));
5187bf11e45SBarry Smith   ts->data = (void *) pseudo;
5192d3f70b5SBarry Smith 
52028aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
521ca4b7087SBarry Smith   pseudo->increment_dt_from_initial_dt = 0;
52228aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
5233a40ed3dSBarry Smith   PetscFunctionReturn(0);
5242d3f70b5SBarry Smith }
5252d3f70b5SBarry Smith 
5262d3f70b5SBarry Smith 
52758562591SBarry Smith #undef __FUNC__
52858562591SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement"
52928aa8177SBarry Smith /*@
53028aa8177SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
53128aa8177SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
53228aa8177SBarry Smith 
53328aa8177SBarry Smith     Input Parameters:
53428aa8177SBarry Smith .   ts - the timestep context
53528aa8177SBarry Smith .   inc - the scaling factor >= 1.0
53628aa8177SBarry Smith 
53728aa8177SBarry Smith     Options Database Key:
53828aa8177SBarry Smith $    -ts_pseudo_increment <increment>
53928aa8177SBarry Smith 
540564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, increment
541564e8f4eSLois Curfman McInnes 
54228aa8177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
54328aa8177SBarry Smith @*/
54428aa8177SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc)
54528aa8177SBarry Smith {
54628aa8177SBarry Smith   TS_Pseudo *pseudo;
54728aa8177SBarry Smith 
5483a40ed3dSBarry Smith   PetscFunctionBegin;
54928aa8177SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
5503a40ed3dSBarry Smith   if (ts->type != TS_PSEUDO) PetscFunctionReturn(0);
55128aa8177SBarry Smith 
55228aa8177SBarry Smith   pseudo               = (TS_Pseudo*) ts->data;
55328aa8177SBarry Smith   pseudo->dt_increment = inc;
5543a40ed3dSBarry Smith   PetscFunctionReturn(0);
55528aa8177SBarry Smith }
5562d3f70b5SBarry Smith 
557ca4b7087SBarry Smith #undef __FUNC__
558ca4b7087SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt"
559ca4b7087SBarry Smith /*@
560b5ed506cSLois Curfman McInnes     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
561b5ed506cSLois Curfman McInnes     is computed via the formula
562b5ed506cSLois Curfman McInnes $         dt = initial_dt*initial_fnorm/current_fnorm
563b5ed506cSLois Curfman McInnes       rather than the default update,
564b5ed506cSLois Curfman McInnes $         dt = current_dt*previous_fnorm/current_fnorm.
565ca4b7087SBarry Smith 
5660704ad96SLois Curfman McInnes     Input Parameter:
567ca4b7087SBarry Smith .   ts - the timestep context
568ca4b7087SBarry Smith 
569ca4b7087SBarry Smith     Options Database Key:
570ca4b7087SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
571ca4b7087SBarry Smith 
572ca4b7087SBarry Smith .keywords: timestep, pseudo, set, increment
573ca4b7087SBarry Smith 
574ca4b7087SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
575ca4b7087SBarry Smith @*/
576ca4b7087SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts)
577ca4b7087SBarry Smith {
578ca4b7087SBarry Smith   TS_Pseudo *pseudo;
579ca4b7087SBarry Smith 
5803a40ed3dSBarry Smith   PetscFunctionBegin;
581ca4b7087SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
5823a40ed3dSBarry Smith   if (ts->type != TS_PSEUDO) PetscFunctionReturn(0);
583ca4b7087SBarry Smith 
584ca4b7087SBarry Smith   pseudo                               = (TS_Pseudo*) ts->data;
585ca4b7087SBarry Smith   pseudo->increment_dt_from_initial_dt = 1;
5863a40ed3dSBarry Smith   PetscFunctionReturn(0);
587ca4b7087SBarry Smith }
5882d3f70b5SBarry Smith 
5892d3f70b5SBarry Smith 
590