1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*e4ef1970SSatish Balay static char vcid[] = "$Id: posindep.c,v 1.28 1998/05/29 20:38:26 bsmith Exp balay $"; 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 /* ------------------------------------------------------------------------------*/ 322d3f70b5SBarry Smith 3358562591SBarry Smith #undef __FUNC__ 3458562591SBarry Smith #define __FUNC__ "TSPseudoComputeTimeStep" 357bf11e45SBarry Smith /*@ 367bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 37564e8f4eSLois Curfman McInnes pseudo-timestepping process. 382d3f70b5SBarry Smith 397bf11e45SBarry Smith Input Parameter: 407bf11e45SBarry Smith . ts - timestep context 417bf11e45SBarry Smith 427bf11e45SBarry Smith Output Parameter: 43fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 44fb4a63b6SLois Curfman McInnes 45fee21e36SBarry Smith Collective on TS 46564e8f4eSLois Curfman McInnes 47564e8f4eSLois Curfman McInnes Notes: 48564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 49564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 50564e8f4eSLois Curfman McInnes 51fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 52564e8f4eSLois Curfman McInnes 53564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 547bf11e45SBarry Smith @*/ 557bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 567bf11e45SBarry Smith { 577bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 587bf11e45SBarry Smith int ierr; 597bf11e45SBarry Smith 603a40ed3dSBarry Smith PetscFunctionBegin; 617bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 627bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr); 637bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 643a40ed3dSBarry Smith PetscFunctionReturn(0); 657bf11e45SBarry Smith } 667bf11e45SBarry Smith 677bf11e45SBarry Smith 687bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 6958562591SBarry Smith #undef __FUNC__ 7058562591SBarry Smith #define __FUNC__ "TSPseudoDefaultVerifyTimeStep" 717bf11e45SBarry Smith /*@C 72639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 737bf11e45SBarry Smith 747bf11e45SBarry Smith Input Parameters: 757bf11e45SBarry Smith . ts - the timestep context 767bf11e45SBarry Smith . dtctx - unused timestep context 77564e8f4eSLois Curfman McInnes . update - latest solution vector 787bf11e45SBarry Smith 79564e8f4eSLois Curfman McInnes Output Parameters: 807bf11e45SBarry Smith . newdt - the timestep to use for the next step 81564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 827bf11e45SBarry Smith 83fee21e36SBarry Smith Collective on TS 84fee21e36SBarry Smith 85564e8f4eSLois Curfman McInnes Note: 86564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 87564e8f4eSLois Curfman McInnes timestep. 88564e8f4eSLois Curfman McInnes 89564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 90564e8f4eSLois Curfman McInnes 91564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 927bf11e45SBarry Smith @*/ 937bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag) 947bf11e45SBarry Smith { 953a40ed3dSBarry Smith PetscFunctionBegin; 967bf11e45SBarry Smith *flag = 1; 973a40ed3dSBarry Smith PetscFunctionReturn(0); 987bf11e45SBarry Smith } 997bf11e45SBarry Smith 1007bf11e45SBarry Smith 10158562591SBarry Smith #undef __FUNC__ 10258562591SBarry Smith #define __FUNC__ "TSPseudoVerifyTimeStep" 1037bf11e45SBarry Smith /*@ 104564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1057bf11e45SBarry Smith 106fb4a63b6SLois Curfman McInnes Input Parameters: 1077bf11e45SBarry Smith . ts - timestep context 108564e8f4eSLois Curfman McInnes . update - latest solution vector 1097bf11e45SBarry Smith 110fb4a63b6SLois Curfman McInnes Output Parameters: 111fb4a63b6SLois Curfman McInnes . dt - newly computed timestep (if it had to shrink) 1127bf11e45SBarry Smith . flag - indicates if current timestep was ok 1137bf11e45SBarry Smith 114fee21e36SBarry Smith Collective on TS 115fee21e36SBarry Smith 116564e8f4eSLois Curfman McInnes Notes: 117564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 118564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 119564e8f4eSLois Curfman McInnes 120564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 121564e8f4eSLois Curfman McInnes 122564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1237bf11e45SBarry Smith @*/ 1247bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 1257bf11e45SBarry Smith { 1267bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1277bf11e45SBarry Smith int ierr; 1287bf11e45SBarry Smith 1293a40ed3dSBarry Smith PetscFunctionBegin; 1303a40ed3dSBarry Smith if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);} 1317bf11e45SBarry Smith 1327bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 1337bf11e45SBarry Smith 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 1357bf11e45SBarry Smith } 1367bf11e45SBarry Smith 1377bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1387bf11e45SBarry Smith 13958562591SBarry Smith #undef __FUNC__ 14058562591SBarry Smith #define __FUNC__ "TSStep_Pseudo" 1417bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 1422d3f70b5SBarry Smith { 1432d3f70b5SBarry Smith Vec sol = ts->vec_sol; 144f2267985SLois Curfman McInnes int ierr,i,max_steps = ts->max_steps,its,ok,lits; 1457bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1467bf11e45SBarry Smith double current_time_step; 1472d3f70b5SBarry Smith 1483a40ed3dSBarry Smith PetscFunctionBegin; 1492d3f70b5SBarry Smith *steps = -ts->steps; 1502d3f70b5SBarry Smith 1517bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 1522d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 1537bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 1547bf11e45SBarry Smith current_time_step = ts->time_step; 1557bf11e45SBarry Smith while (1) { 1567bf11e45SBarry Smith ts->ptime += current_time_step; 1577bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 158f2267985SLois Curfman McInnes ierr = SNESGetNumberLinearIterations(ts->snes,&lits); CHKERRQ(ierr); 159f2267985SLois Curfman McInnes ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits; 1607bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 1617bf11e45SBarry Smith if (ok) break; 1627bf11e45SBarry Smith ts->ptime -= current_time_step; 1637bf11e45SBarry Smith current_time_step = ts->time_step; 1647bf11e45SBarry Smith } 1657bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 1662d3f70b5SBarry Smith ts->steps++; 1672d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1682d3f70b5SBarry Smith } 1692d3f70b5SBarry Smith 1702d3f70b5SBarry Smith *steps += ts->steps; 1712d3f70b5SBarry Smith *time = ts->ptime; 1723a40ed3dSBarry Smith PetscFunctionReturn(0); 1732d3f70b5SBarry Smith } 1742d3f70b5SBarry Smith 1752d3f70b5SBarry Smith /*------------------------------------------------------------*/ 17658562591SBarry Smith #undef __FUNC__ 177d4bb536fSBarry Smith #define __FUNC__ "TSDestroy_Pseudo" 178e1311b90SBarry Smith static int TSDestroy_Pseudo(TS ts ) 1792d3f70b5SBarry Smith { 1807bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1812d3f70b5SBarry Smith int ierr; 1822d3f70b5SBarry Smith 1833a40ed3dSBarry Smith PetscFunctionBegin; 184*e4ef1970SSatish Balay if (pseudo->update) {ierr = VecDestroy(pseudo->update); CHKERRQ(ierr);} 1857bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1867bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 1872d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 1887bf11e45SBarry Smith PetscFree(pseudo); 1893a40ed3dSBarry Smith PetscFunctionReturn(0); 1902d3f70b5SBarry Smith } 1912d3f70b5SBarry Smith 1922d3f70b5SBarry Smith 1932d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1942d3f70b5SBarry Smith /* 1952d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 1962d3f70b5SBarry Smith */ 1972d3f70b5SBarry Smith 19858562591SBarry Smith #undef __FUNC__ 19958562591SBarry Smith #define __FUNC__ "TSPseudoMatMult" 2007bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 2012d3f70b5SBarry Smith { 2022d3f70b5SBarry Smith TS ts; 2032d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 2042d3f70b5SBarry Smith int ierr; 2052d3f70b5SBarry Smith 2063a40ed3dSBarry Smith PetscFunctionBegin; 2073a40ed3dSBarry Smith ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr); 2082d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2092d3f70b5SBarry Smith 2102d3f70b5SBarry Smith /* apply user provided function */ 2112d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 2122d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2132d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 2143a40ed3dSBarry Smith PetscFunctionReturn(0); 2152d3f70b5SBarry Smith } 2162d3f70b5SBarry Smith 2172d3f70b5SBarry Smith /* 2182d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2192d3f70b5SBarry Smith 2202d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2212d3f70b5SBarry Smith */ 22258562591SBarry Smith #undef __FUNC__ 22358562591SBarry Smith #define __FUNC__ "TSPseudoFunction" 2247bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2252d3f70b5SBarry Smith { 2262d3f70b5SBarry Smith TS ts = (TS) ctx; 2272d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2282d3f70b5SBarry Smith int ierr,i,n; 2292d3f70b5SBarry Smith 2303a40ed3dSBarry Smith PetscFunctionBegin; 2312d3f70b5SBarry Smith /* apply user provided function */ 2322d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 2337bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2342d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 2352d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 2362d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 2372d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 2382d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 2392d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2402d3f70b5SBarry Smith } 2412d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 2422d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 2432d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 2442d3f70b5SBarry Smith 2453a40ed3dSBarry Smith PetscFunctionReturn(0); 2462d3f70b5SBarry Smith } 2472d3f70b5SBarry Smith 2482d3f70b5SBarry Smith /* 2492d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2502d3f70b5SBarry Smith 2512d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2522d3f70b5SBarry Smith */ 25358562591SBarry Smith #undef __FUNC__ 25458562591SBarry Smith #define __FUNC__ "TSPseudoJacobian" 2557bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2562d3f70b5SBarry Smith { 2572d3f70b5SBarry Smith TS ts = (TS) ctx; 2582d3f70b5SBarry Smith int ierr; 2592d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 2602d3f70b5SBarry Smith MatType mtype; 2612d3f70b5SBarry Smith 2623a40ed3dSBarry Smith PetscFunctionBegin; 2632d3f70b5SBarry Smith /* construct users Jacobian */ 2642d3f70b5SBarry Smith if (ts->rhsjacobian) { 2652d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 2662d3f70b5SBarry Smith } 2672d3f70b5SBarry Smith 2682d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 2692d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 2702d3f70b5SBarry Smith if (mtype != MATSHELL) { 2712d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 2722d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 2732d3f70b5SBarry Smith } 2742d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 2752d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 2762d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 2772d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 2782d3f70b5SBarry Smith } 2792d3f70b5SBarry Smith 2803a40ed3dSBarry Smith PetscFunctionReturn(0); 2812d3f70b5SBarry Smith } 2822d3f70b5SBarry Smith 2832d3f70b5SBarry Smith 28458562591SBarry Smith #undef __FUNC__ 28558562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo" 2867bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 2872d3f70b5SBarry Smith { 2887bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2892d3f70b5SBarry Smith int ierr, M, m; 2902d3f70b5SBarry Smith 2913a40ed3dSBarry Smith PetscFunctionBegin; 2927bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 2937bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 2947bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 2952d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 2962d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 2972d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 2982d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 2997bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 3002d3f70b5SBarry Smith } 3017bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 3023a40ed3dSBarry Smith PetscFunctionReturn(0); 3032d3f70b5SBarry Smith } 3042d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3052d3f70b5SBarry Smith 30658562591SBarry Smith #undef __FUNC__ 307d4bb536fSBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor" 3082d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 3092d3f70b5SBarry Smith { 3107bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3112d3f70b5SBarry Smith 3123a40ed3dSBarry Smith PetscFunctionBegin; 31376be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 3143a40ed3dSBarry Smith PetscFunctionReturn(0); 3152d3f70b5SBarry Smith } 3162d3f70b5SBarry Smith 31758562591SBarry Smith #undef __FUNC__ 31858562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo" 3197bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3202d3f70b5SBarry Smith { 3212d3f70b5SBarry Smith int ierr,flg; 32228aa8177SBarry Smith double inc; 3232d3f70b5SBarry Smith 3243a40ed3dSBarry Smith PetscFunctionBegin; 3252d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 3262d3f70b5SBarry Smith 3272d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 3282d3f70b5SBarry Smith if (flg) { 32928aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 33028aa8177SBarry Smith } 33128aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 33228aa8177SBarry Smith if (flg) { 33328aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 3342d3f70b5SBarry Smith } 335ca4b7087SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr); 336ca4b7087SBarry Smith if (flg) { 337ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr); 338ca4b7087SBarry Smith } 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 3402d3f70b5SBarry Smith } 3412d3f70b5SBarry Smith 34258562591SBarry Smith #undef __FUNC__ 343d4bb536fSBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo" 344ca4b7087SBarry Smith static int TSPrintHelp_Pseudo(TS ts,char *p) 3452d3f70b5SBarry Smith { 3463a40ed3dSBarry Smith PetscFunctionBegin; 34776be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n"); 34876be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p); 34976be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p); 35076be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," initial fnorm/current fnorm to determine new timestep\n"); 35176be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," default is current_dt * previous fnorm/current fnorm\n"); 3523a40ed3dSBarry Smith PetscFunctionReturn(0); 3532d3f70b5SBarry Smith } 3542d3f70b5SBarry Smith 35558562591SBarry Smith #undef __FUNC__ 356d4bb536fSBarry Smith #define __FUNC__ "TSView_Pseudo" 357e1311b90SBarry Smith static int TSView_Pseudo(TS ts,Viewer viewer) 3582d3f70b5SBarry Smith { 3593a40ed3dSBarry Smith PetscFunctionBegin; 3603a40ed3dSBarry Smith PetscFunctionReturn(0); 3612d3f70b5SBarry Smith } 3622d3f70b5SBarry Smith 36382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 36482bf6240SBarry Smith #undef __FUNC__ 36582bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep" 36682bf6240SBarry Smith /*@ 36782bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 36882bf6240SBarry Smith last timestep. 36982bf6240SBarry Smith 37082bf6240SBarry Smith Input Parameters: 37182bf6240SBarry Smith . ts - timestep context 37282bf6240SBarry Smith . dt - user-defined function to verify timestep 37382bf6240SBarry Smith . ctx - [optional] user-defined context for private data 37482bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 37582bf6240SBarry Smith 376fee21e36SBarry Smith Collective on TS 377fee21e36SBarry Smith 37882bf6240SBarry Smith Calling sequence of func: 37982bf6240SBarry Smith . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 38082bf6240SBarry Smith 38182bf6240SBarry Smith . update - latest solution vector 38282bf6240SBarry Smith . ctx - [optional] timestep context 38382bf6240SBarry Smith . newdt - the timestep to use for the next step 38482bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 38582bf6240SBarry Smith 38682bf6240SBarry Smith Notes: 38782bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 38882bf6240SBarry Smith during the timestepping process. 38982bf6240SBarry Smith 39082bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 39182bf6240SBarry Smith 39282bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 39382bf6240SBarry Smith @*/ 39482bf6240SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 39582bf6240SBarry Smith { 39682bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,Vec,void*,double *,int *),void *); 39782bf6240SBarry Smith 39882bf6240SBarry Smith PetscFunctionBegin; 39982bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 40082bf6240SBarry Smith 40194d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void **)&f);CHKERRQ(ierr); 40282bf6240SBarry Smith if (f) { 40382bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 40482bf6240SBarry Smith } 40582bf6240SBarry Smith PetscFunctionReturn(0); 40682bf6240SBarry Smith } 40782bf6240SBarry Smith 40882bf6240SBarry Smith #undef __FUNC__ 40982bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement" 41082bf6240SBarry Smith /*@ 41182bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 41282bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 41382bf6240SBarry Smith 41482bf6240SBarry Smith Input Parameters: 41582bf6240SBarry Smith . ts - the timestep context 41682bf6240SBarry Smith . inc - the scaling factor >= 1.0 41782bf6240SBarry Smith 418fee21e36SBarry Smith Collective on TS 419fee21e36SBarry Smith 42082bf6240SBarry Smith Options Database Key: 42182bf6240SBarry Smith $ -ts_pseudo_increment <increment> 42282bf6240SBarry Smith 42382bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 42482bf6240SBarry Smith 42582bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 42682bf6240SBarry Smith @*/ 42782bf6240SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 42882bf6240SBarry Smith { 42982bf6240SBarry Smith int ierr, (*f)(TS,double); 43082bf6240SBarry Smith 43182bf6240SBarry Smith PetscFunctionBegin; 43282bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 43382bf6240SBarry Smith 43494d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void **)&f);CHKERRQ(ierr); 43582bf6240SBarry Smith if (f) { 43682bf6240SBarry Smith ierr = (*f)(ts,inc);CHKERRQ(ierr); 43782bf6240SBarry Smith } 43882bf6240SBarry Smith PetscFunctionReturn(0); 43982bf6240SBarry Smith } 44082bf6240SBarry Smith 44182bf6240SBarry Smith #undef __FUNC__ 44282bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt" 44382bf6240SBarry Smith /*@ 44482bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 44582bf6240SBarry Smith is computed via the formula 44682bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 44782bf6240SBarry Smith rather than the default update, 44882bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 44982bf6240SBarry Smith 45082bf6240SBarry Smith Input Parameter: 45182bf6240SBarry Smith . ts - the timestep context 45282bf6240SBarry Smith 453fee21e36SBarry Smith Collective on TS 454fee21e36SBarry Smith 45582bf6240SBarry Smith Options Database Key: 45682bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 45782bf6240SBarry Smith 45882bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 45982bf6240SBarry Smith 46082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 46182bf6240SBarry Smith @*/ 46282bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts) 46382bf6240SBarry Smith { 46482bf6240SBarry Smith int ierr, (*f)(TS); 46582bf6240SBarry Smith 46682bf6240SBarry Smith PetscFunctionBegin; 46782bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 46882bf6240SBarry Smith 46994d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void **)&f);CHKERRQ(ierr); 47082bf6240SBarry Smith if (f) { 47182bf6240SBarry Smith ierr = (*f)(ts);CHKERRQ(ierr); 47282bf6240SBarry Smith } 47382bf6240SBarry Smith PetscFunctionReturn(0); 47482bf6240SBarry Smith } 47582bf6240SBarry Smith 47682bf6240SBarry Smith 47782bf6240SBarry Smith #undef __FUNC__ 47882bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep" 47982bf6240SBarry Smith /*@ 48082bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 48182bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 48282bf6240SBarry Smith 48382bf6240SBarry Smith Input Parameters: 48482bf6240SBarry Smith . ts - timestep context 48582bf6240SBarry Smith . dt - function to compute timestep 48682bf6240SBarry Smith . ctx - [optional] user-defined context for private data 48782bf6240SBarry Smith required by the function (may be PETSC_NULL) 48882bf6240SBarry Smith 489fee21e36SBarry Smith Collective on TS 490fee21e36SBarry Smith 49182bf6240SBarry Smith Calling sequence of func: 49282bf6240SBarry Smith . func (TS ts,double *newdt,void *ctx); 49382bf6240SBarry Smith 49482bf6240SBarry Smith . newdt - the newly computed timestep 49582bf6240SBarry Smith . ctx - [optional] timestep context 49682bf6240SBarry Smith 49782bf6240SBarry Smith Notes: 49882bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 49982bf6240SBarry Smith during the timestepping process. 50082bf6240SBarry Smith 50182bf6240SBarry Smith .keywords: timestep, pseudo, set 50282bf6240SBarry Smith 50382bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 50482bf6240SBarry Smith @*/ 50582bf6240SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 50682bf6240SBarry Smith { 50782bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,double *,void *),void *); 50882bf6240SBarry Smith 50982bf6240SBarry Smith PetscFunctionBegin; 51082bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 51182bf6240SBarry Smith 51294d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void **)&f);CHKERRQ(ierr); 51382bf6240SBarry Smith if (f) { 51482bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 51582bf6240SBarry Smith } 51682bf6240SBarry Smith PetscFunctionReturn(0); 51782bf6240SBarry Smith } 51882bf6240SBarry Smith 51982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 52082bf6240SBarry Smith 52182bf6240SBarry Smith #undef __FUNC__ 52282bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep_Pseudo" 52382bf6240SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 52482bf6240SBarry Smith { 52582bf6240SBarry Smith TS_Pseudo *pseudo; 52682bf6240SBarry Smith 52782bf6240SBarry Smith PetscFunctionBegin; 52882bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 52982bf6240SBarry Smith pseudo->verify = dt; 53082bf6240SBarry Smith pseudo->verifyctx = ctx; 53182bf6240SBarry Smith PetscFunctionReturn(0); 53282bf6240SBarry Smith } 53382bf6240SBarry Smith 53482bf6240SBarry Smith #undef __FUNC__ 53582bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement_Pseudo" 53682bf6240SBarry Smith int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,double inc) 53782bf6240SBarry Smith { 53882bf6240SBarry Smith TS_Pseudo *pseudo; 53982bf6240SBarry Smith 54082bf6240SBarry Smith PetscFunctionBegin; 54182bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 54282bf6240SBarry Smith pseudo->dt_increment = inc; 54382bf6240SBarry Smith PetscFunctionReturn(0); 54482bf6240SBarry Smith } 54582bf6240SBarry Smith 54682bf6240SBarry Smith #undef __FUNC__ 54782bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 54882bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 54982bf6240SBarry Smith { 55082bf6240SBarry Smith TS_Pseudo *pseudo; 55182bf6240SBarry Smith 55282bf6240SBarry Smith PetscFunctionBegin; 55382bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 55482bf6240SBarry Smith pseudo->increment_dt_from_initial_dt = 1; 55582bf6240SBarry Smith PetscFunctionReturn(0); 55682bf6240SBarry Smith } 55782bf6240SBarry Smith 55882bf6240SBarry Smith #undef __FUNC__ 55982bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep_Pseudo" 56082bf6240SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,double*,void*),void* ctx) 56182bf6240SBarry Smith { 56282bf6240SBarry Smith TS_Pseudo *pseudo; 56382bf6240SBarry Smith 56482bf6240SBarry Smith PetscFunctionBegin; 56582bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 56682bf6240SBarry Smith pseudo->dt = dt; 56782bf6240SBarry Smith pseudo->dtctx = ctx; 56882bf6240SBarry Smith PetscFunctionReturn(0); 56982bf6240SBarry Smith } 57082bf6240SBarry Smith 57182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 57282bf6240SBarry Smith 57358562591SBarry Smith #undef __FUNC__ 57458562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo" 5757bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 5762d3f70b5SBarry Smith { 5777bf11e45SBarry Smith TS_Pseudo *pseudo; 5782d3f70b5SBarry Smith int ierr; 5792d3f70b5SBarry Smith MatType mtype; 5802d3f70b5SBarry Smith 5813a40ed3dSBarry Smith PetscFunctionBegin; 5827bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 5837bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 5847bf11e45SBarry Smith ts->view = TSView_Pseudo; 5852d3f70b5SBarry Smith 5862d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 587a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems"); 5882d3f70b5SBarry Smith } 5892d3f70b5SBarry Smith if (!ts->A) { 590a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian"); 5912d3f70b5SBarry Smith } 5922d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 5932d3f70b5SBarry Smith if (mtype == MATSHELL) { 5942d3f70b5SBarry Smith ts->Ashell = ts->A; 5952d3f70b5SBarry Smith } 5967bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 5977bf11e45SBarry Smith ts->step = TSStep_Pseudo; 5987bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 5997bf11e45SBarry Smith 6007bf11e45SBarry Smith /* create the required nonlinear solver context */ 6012d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 6022d3f70b5SBarry Smith 6037bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 604eed86810SBarry Smith PLogObjectMemory(ts,sizeof(TS_Pseudo)); 605eed86810SBarry Smith 6067bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 6077bf11e45SBarry Smith ts->data = (void *) pseudo; 6082d3f70b5SBarry Smith 60928aa8177SBarry Smith pseudo->dt_increment = 1.1; 610ca4b7087SBarry Smith pseudo->increment_dt_from_initial_dt = 0; 61128aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 61282bf6240SBarry Smith 61394d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 614e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 615e1311b90SBarry Smith (void*)TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 61694d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 617e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 618e1311b90SBarry Smith (void*)TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 61994d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 620e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 621e1311b90SBarry Smith (void*)TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 62294d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C","TSPseudoSetTimeStep_Pseudo", 623e1311b90SBarry Smith (void*)TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6243a40ed3dSBarry Smith PetscFunctionReturn(0); 6252d3f70b5SBarry Smith } 6262d3f70b5SBarry Smith 62758562591SBarry Smith #undef __FUNC__ 62882bf6240SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep" 62982bf6240SBarry Smith /*@C 63082bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 63182bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 63228aa8177SBarry Smith 63328aa8177SBarry Smith Input Parameters: 63428aa8177SBarry Smith . ts - the timestep context 63582bf6240SBarry Smith . dtctx - unused timestep context 63628aa8177SBarry Smith 637fee21e36SBarry Smith Collective on TS 638fee21e36SBarry Smith 63982bf6240SBarry Smith Output Parameter: 64082bf6240SBarry Smith . newdt - the timestep to use for the next step 64128aa8177SBarry Smith 64282bf6240SBarry Smith .keywords: timestep, pseudo, default 643564e8f4eSLois Curfman McInnes 64482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 64528aa8177SBarry Smith @*/ 64682bf6240SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 64728aa8177SBarry Smith { 64882bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 64982bf6240SBarry Smith double inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 65082bf6240SBarry Smith int ierr; 65128aa8177SBarry Smith 6523a40ed3dSBarry Smith PetscFunctionBegin; 65382bf6240SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 65482bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 65582bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 65682bf6240SBarry Smith /* first time through so compute initial function norm */ 65782bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 65882bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 65982bf6240SBarry Smith } 66082bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 66182bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 66282bf6240SBarry Smith } 66382bf6240SBarry Smith else if (pseudo->increment_dt_from_initial_dt) { 66482bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 66582bf6240SBarry Smith } else { 66682bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 66782bf6240SBarry Smith } 66882bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6693a40ed3dSBarry Smith PetscFunctionReturn(0); 67028aa8177SBarry Smith } 6712d3f70b5SBarry Smith 672