173f4d377SMatthew Knepley /*$Id: posindep.c,v 1.60 2001/09/11 16:34:22 bsmith Exp $*/ 22d3f70b5SBarry Smith /* 3fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 42d3f70b5SBarry Smith */ 5e090d566SSatish Balay #include "src/ts/tsimpl.h" /*I "petscts.h" I*/ 62d3f70b5SBarry Smith 72d3f70b5SBarry Smith typedef struct { 82d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 92d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 102d3f70b5SBarry Smith Vec rhs; /* work vector for RHS; vec_sol/dt */ 112d3f70b5SBarry Smith 122d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 132d3f70b5SBarry Smith 1487828ca2SBarry Smith int (*dt)(TS,PetscReal*,void*); /* compute next timestep, and related context */ 152d3f70b5SBarry Smith void *dtctx; 1687828ca2SBarry Smith int (*verify)(TS,Vec,void*,PetscReal*,int*); /* verify previous timestep and related context */ 177bf11e45SBarry Smith void *verifyctx; 182d3f70b5SBarry Smith 1987828ca2SBarry Smith PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */ 2087828ca2SBarry Smith PetscReal fnorm_previous; 2128aa8177SBarry Smith 2287828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 234bbc92c1SBarry Smith PetscTruth increment_dt_from_initial_dt; 247bf11e45SBarry Smith } TS_Pseudo; 252d3f70b5SBarry Smith 262d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 272d3f70b5SBarry Smith 284a2ae208SSatish Balay #undef __FUNCT__ 294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep" 307bf11e45SBarry Smith /*@ 317bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 32564e8f4eSLois Curfman McInnes pseudo-timestepping process. 332d3f70b5SBarry Smith 3415091d37SBarry Smith Collective on TS 3515091d37SBarry Smith 367bf11e45SBarry Smith Input Parameter: 377bf11e45SBarry Smith . ts - timestep context 387bf11e45SBarry Smith 397bf11e45SBarry Smith Output Parameter: 40fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 41fb4a63b6SLois Curfman McInnes 4215091d37SBarry Smith Level: advanced 43564e8f4eSLois Curfman McInnes 44564e8f4eSLois Curfman McInnes Notes: 45564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 46564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 47564e8f4eSLois Curfman McInnes 48fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 49564e8f4eSLois Curfman McInnes 50564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 517bf11e45SBarry Smith @*/ 5287828ca2SBarry Smith int TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 537bf11e45SBarry Smith { 547bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 557bf11e45SBarry Smith int ierr; 567bf11e45SBarry Smith 573a40ed3dSBarry Smith PetscFunctionBegin; 58d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 597bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 60d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 613a40ed3dSBarry Smith PetscFunctionReturn(0); 627bf11e45SBarry Smith } 637bf11e45SBarry Smith 647bf11e45SBarry Smith 657bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 664a2ae208SSatish Balay #undef __FUNCT__ 674a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 687bf11e45SBarry Smith /*@C 69639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 707bf11e45SBarry Smith 7115091d37SBarry Smith Collective on TS 7215091d37SBarry Smith 737bf11e45SBarry Smith Input Parameters: 7415091d37SBarry Smith + ts - the timestep context 757bf11e45SBarry Smith . dtctx - unused timestep context 7615091d37SBarry Smith - update - latest solution vector 777bf11e45SBarry Smith 78564e8f4eSLois Curfman McInnes Output Parameters: 7915091d37SBarry Smith + newdt - the timestep to use for the next step 8015091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 817bf11e45SBarry Smith 8215091d37SBarry Smith Level: advanced 83fee21e36SBarry Smith 84564e8f4eSLois Curfman McInnes Note: 85564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 86564e8f4eSLois Curfman McInnes timestep. 87564e8f4eSLois Curfman McInnes 88564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 89564e8f4eSLois Curfman McInnes 90564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 917bf11e45SBarry Smith @*/ 9287828ca2SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,int *flag) 937bf11e45SBarry Smith { 943a40ed3dSBarry Smith PetscFunctionBegin; 957bf11e45SBarry Smith *flag = 1; 963a40ed3dSBarry Smith PetscFunctionReturn(0); 977bf11e45SBarry Smith } 987bf11e45SBarry Smith 997bf11e45SBarry Smith 1004a2ae208SSatish Balay #undef __FUNCT__ 1014a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep" 1027bf11e45SBarry Smith /*@ 103564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1047bf11e45SBarry Smith 10515091d37SBarry Smith Collective on TS 10615091d37SBarry Smith 107fb4a63b6SLois Curfman McInnes Input Parameters: 10815091d37SBarry Smith + ts - timestep context 10915091d37SBarry Smith - update - latest solution vector 1107bf11e45SBarry Smith 111fb4a63b6SLois Curfman McInnes Output Parameters: 11215091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11315091d37SBarry Smith - flag - indicates if current timestep was ok 1147bf11e45SBarry Smith 11515091d37SBarry Smith Level: advanced 116fee21e36SBarry Smith 117564e8f4eSLois Curfman McInnes Notes: 118564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 119564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 120564e8f4eSLois Curfman McInnes 121564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 122564e8f4eSLois Curfman McInnes 123564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1247bf11e45SBarry Smith @*/ 12587828ca2SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,int *flag) 1267bf11e45SBarry Smith { 1277bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 1287bf11e45SBarry Smith int ierr; 1297bf11e45SBarry Smith 1303a40ed3dSBarry Smith PetscFunctionBegin; 1313a40ed3dSBarry Smith if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);} 1327bf11e45SBarry Smith 1337bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 1347bf11e45SBarry Smith 1353a40ed3dSBarry Smith PetscFunctionReturn(0); 1367bf11e45SBarry Smith } 1377bf11e45SBarry Smith 1387bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1397bf11e45SBarry Smith 1404a2ae208SSatish Balay #undef __FUNCT__ 1414a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo" 14287828ca2SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,PetscReal *ptime) 1432d3f70b5SBarry Smith { 1442d3f70b5SBarry Smith Vec sol = ts->vec_sol; 145f2267985SLois Curfman McInnes int ierr,i,max_steps = ts->max_steps,its,ok,lits; 1467bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 14787828ca2SBarry Smith PetscReal current_time_step; 1482d3f70b5SBarry Smith 1493a40ed3dSBarry Smith PetscFunctionBegin; 1502d3f70b5SBarry Smith *steps = -ts->steps; 1512d3f70b5SBarry Smith 1527bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr); 1532d3f70b5SBarry Smith for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) { 1547bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 155e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1567bf11e45SBarry Smith current_time_step = ts->time_step; 157e6e5fe25SBarry Smith while (PETSC_TRUE) { 1587bf11e45SBarry Smith ts->ptime += current_time_step; 1597bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its);CHKERRQ(ierr); 160f2267985SLois Curfman McInnes ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr); 161f2267985SLois Curfman McInnes ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits; 1627bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr); 1637bf11e45SBarry Smith if (ok) break; 1647bf11e45SBarry Smith ts->ptime -= current_time_step; 1657bf11e45SBarry Smith current_time_step = ts->time_step; 1667bf11e45SBarry Smith } 1677bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr); 1682d3f70b5SBarry Smith ts->steps++; 1692d3f70b5SBarry Smith } 170e6e5fe25SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 171e6e5fe25SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 172e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1732d3f70b5SBarry Smith 1742d3f70b5SBarry Smith *steps += ts->steps; 175142b95e3SSatish Balay *ptime = ts->ptime; 1763a40ed3dSBarry Smith PetscFunctionReturn(0); 1772d3f70b5SBarry Smith } 1782d3f70b5SBarry Smith 1792d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1804a2ae208SSatish Balay #undef __FUNCT__ 1814a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo" 182e1311b90SBarry Smith static int TSDestroy_Pseudo(TS ts) 1832d3f70b5SBarry Smith { 1847bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 1852d3f70b5SBarry Smith int ierr; 1862d3f70b5SBarry Smith 1873a40ed3dSBarry Smith PetscFunctionBegin; 188e4ef1970SSatish Balay if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 1897bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1907bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 191606d414cSSatish Balay ierr = PetscFree(pseudo);CHKERRQ(ierr); 1923a40ed3dSBarry Smith PetscFunctionReturn(0); 1932d3f70b5SBarry Smith } 1942d3f70b5SBarry Smith 1952d3f70b5SBarry Smith 1962d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1972d3f70b5SBarry Smith 1982d3f70b5SBarry Smith /* 1992d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2002d3f70b5SBarry Smith 2012d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2022d3f70b5SBarry Smith */ 2034a2ae208SSatish Balay #undef __FUNCT__ 2044a2ae208SSatish Balay #define __FUNCT__ "TSPseudoFunction" 2057bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2062d3f70b5SBarry Smith { 2072d3f70b5SBarry Smith TS ts = (TS) ctx; 208ea709b57SSatish Balay PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2092d3f70b5SBarry Smith int ierr,i,n; 2102d3f70b5SBarry Smith 2113a40ed3dSBarry Smith PetscFunctionBegin; 2122d3f70b5SBarry Smith /* apply user provided function */ 2132d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr); 2147bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2152d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr); 2162d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1);CHKERRQ(ierr); 2172d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr); 2182d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 2192d3f70b5SBarry Smith for (i=0; i<n; i++) { 2202d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2212d3f70b5SBarry Smith } 222888f2ed8SSatish Balay ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr); 223888f2ed8SSatish Balay ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr); 224888f2ed8SSatish Balay ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr); 2252d3f70b5SBarry Smith 2263a40ed3dSBarry Smith PetscFunctionReturn(0); 2272d3f70b5SBarry Smith } 2282d3f70b5SBarry Smith 2292d3f70b5SBarry Smith /* 2302d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2312d3f70b5SBarry Smith 2322d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2332d3f70b5SBarry Smith */ 2344a2ae208SSatish Balay #undef __FUNCT__ 2354a2ae208SSatish Balay #define __FUNCT__ "TSPseudoJacobian" 2367bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2372d3f70b5SBarry Smith { 2382d3f70b5SBarry Smith TS ts = (TS) ctx; 2392d3f70b5SBarry Smith int ierr; 240ea709b57SSatish Balay PetscScalar mone = -1.0,mdt = 1.0/ts->time_step; 2412d3f70b5SBarry Smith 2423a40ed3dSBarry Smith PetscFunctionBegin; 2432d3f70b5SBarry Smith /* construct users Jacobian */ 244a7a1495cSBarry Smith ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr); 2452d3f70b5SBarry Smith 246614f654bSBarry Smith /* shift and scale Jacobian */ 2472d3f70b5SBarry Smith ierr = MatScale(&mone,*AA);CHKERRQ(ierr); 2482d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA);CHKERRQ(ierr); 249614f654bSBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER) { 2502d3f70b5SBarry Smith ierr = MatScale(&mone,*BB);CHKERRQ(ierr); 2512d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB);CHKERRQ(ierr); 2522d3f70b5SBarry Smith } 2532d3f70b5SBarry Smith 2543a40ed3dSBarry Smith PetscFunctionReturn(0); 2552d3f70b5SBarry Smith } 2562d3f70b5SBarry Smith 2572d3f70b5SBarry Smith 2584a2ae208SSatish Balay #undef __FUNCT__ 2594a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 2607bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 2612d3f70b5SBarry Smith { 2627bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 263b1bcba4aSBarry Smith int ierr; 2642d3f70b5SBarry Smith 2653a40ed3dSBarry Smith PetscFunctionBegin; 266e6e5fe25SBarry Smith /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */ 2677bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 2687bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 2697bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 2707bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 2722d3f70b5SBarry Smith } 2732d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2742d3f70b5SBarry Smith 2754a2ae208SSatish Balay #undef __FUNCT__ 2764a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultMonitor" 27787828ca2SBarry Smith int TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx) 2782d3f70b5SBarry Smith { 2797bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 280d132466eSBarry Smith int ierr; 2812d3f70b5SBarry Smith 2823a40ed3dSBarry Smith PetscFunctionBegin; 283142b95e3SSatish Balay ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 2843a40ed3dSBarry Smith PetscFunctionReturn(0); 2852d3f70b5SBarry Smith } 2862d3f70b5SBarry Smith 2874a2ae208SSatish Balay #undef __FUNCT__ 2884a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 2897bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 2902d3f70b5SBarry Smith { 2914bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 292f1af5d2fSBarry Smith int ierr; 293f1af5d2fSBarry Smith PetscTruth flg; 2942d3f70b5SBarry Smith 2953a40ed3dSBarry Smith PetscFunctionBegin; 2962d3f70b5SBarry Smith 297b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 298b0a32e0cSBarry Smith ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr); 2992d3f70b5SBarry Smith if (flg) { 300329f5518SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 30128aa8177SBarry Smith } 302b0a32e0cSBarry Smith ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr); 303ca4b7087SBarry Smith if (flg) { 304ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 305ca4b7087SBarry Smith } 306419fbf4bSSatish Balay ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 307b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3082d3f70b5SBarry Smith 3093a40ed3dSBarry Smith PetscFunctionReturn(0); 3102d3f70b5SBarry Smith } 3112d3f70b5SBarry Smith 3124a2ae208SSatish Balay #undef __FUNCT__ 3134a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 314b0a32e0cSBarry Smith static int TSView_Pseudo(TS ts,PetscViewer viewer) 3152d3f70b5SBarry Smith { 3163a40ed3dSBarry Smith PetscFunctionBegin; 3173a40ed3dSBarry Smith PetscFunctionReturn(0); 3182d3f70b5SBarry Smith } 3192d3f70b5SBarry Smith 32082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3214a2ae208SSatish Balay #undef __FUNCT__ 3224a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 32382bf6240SBarry Smith /*@ 32482bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 32582bf6240SBarry Smith last timestep. 32682bf6240SBarry Smith 32715091d37SBarry Smith Collective on TS 32815091d37SBarry Smith 32982bf6240SBarry Smith Input Parameters: 33015091d37SBarry Smith + ts - timestep context 33182bf6240SBarry Smith . dt - user-defined function to verify timestep 33215091d37SBarry Smith - ctx - [optional] user-defined context for private data 33382bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 33482bf6240SBarry Smith 33515091d37SBarry Smith Level: advanced 336fee21e36SBarry Smith 33782bf6240SBarry Smith Calling sequence of func: 33887828ca2SBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag); 33982bf6240SBarry Smith 34082bf6240SBarry Smith . update - latest solution vector 34182bf6240SBarry Smith . ctx - [optional] timestep context 34282bf6240SBarry Smith . newdt - the timestep to use for the next step 34382bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 34482bf6240SBarry Smith 34582bf6240SBarry Smith Notes: 34682bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 34782bf6240SBarry Smith during the timestepping process. 34882bf6240SBarry Smith 34982bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 35082bf6240SBarry Smith 35182bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 35282bf6240SBarry Smith @*/ 35387828ca2SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx) 35482bf6240SBarry Smith { 35587828ca2SBarry Smith int ierr,(*f)(TS,int (*)(TS,Vec,void*,PetscReal *,int *),void *); 35682bf6240SBarry Smith 35782bf6240SBarry Smith PetscFunctionBegin; 35882bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 35982bf6240SBarry Smith 360c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 36182bf6240SBarry Smith if (f) { 36282bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 36382bf6240SBarry Smith } 36482bf6240SBarry Smith PetscFunctionReturn(0); 36582bf6240SBarry Smith } 36682bf6240SBarry Smith 3674a2ae208SSatish Balay #undef __FUNCT__ 3684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 36982bf6240SBarry Smith /*@ 37082bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 37182bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 37282bf6240SBarry Smith 373fee21e36SBarry Smith Collective on TS 374fee21e36SBarry Smith 37515091d37SBarry Smith Input Parameters: 37615091d37SBarry Smith + ts - the timestep context 37715091d37SBarry Smith - inc - the scaling factor >= 1.0 37815091d37SBarry Smith 37982bf6240SBarry Smith Options Database Key: 38082bf6240SBarry Smith $ -ts_pseudo_increment <increment> 38182bf6240SBarry Smith 38215091d37SBarry Smith Level: advanced 38315091d37SBarry Smith 38482bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 38582bf6240SBarry Smith 38682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 38782bf6240SBarry Smith @*/ 38887828ca2SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 38982bf6240SBarry Smith { 39087828ca2SBarry Smith int ierr,(*f)(TS,PetscReal); 39182bf6240SBarry Smith 39282bf6240SBarry Smith PetscFunctionBegin; 39382bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 39482bf6240SBarry Smith 395c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr); 39682bf6240SBarry Smith if (f) { 39782bf6240SBarry Smith ierr = (*f)(ts,inc);CHKERRQ(ierr); 39882bf6240SBarry Smith } 39982bf6240SBarry Smith PetscFunctionReturn(0); 40082bf6240SBarry Smith } 40182bf6240SBarry Smith 4024a2ae208SSatish Balay #undef __FUNCT__ 4034a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 40482bf6240SBarry Smith /*@ 40582bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 40682bf6240SBarry Smith is computed via the formula 40782bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 40882bf6240SBarry Smith rather than the default update, 40982bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 41082bf6240SBarry Smith 41115091d37SBarry Smith Collective on TS 41215091d37SBarry Smith 41382bf6240SBarry Smith Input Parameter: 41482bf6240SBarry Smith . ts - the timestep context 41582bf6240SBarry Smith 41682bf6240SBarry Smith Options Database Key: 41782bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 41882bf6240SBarry Smith 41915091d37SBarry Smith Level: advanced 42015091d37SBarry Smith 42182bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 42282bf6240SBarry Smith 42382bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 42482bf6240SBarry Smith @*/ 42582bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts) 42682bf6240SBarry Smith { 42782bf6240SBarry Smith int ierr,(*f)(TS); 42882bf6240SBarry Smith 42982bf6240SBarry Smith PetscFunctionBegin; 43082bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 43182bf6240SBarry Smith 432c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr); 43382bf6240SBarry Smith if (f) { 43482bf6240SBarry Smith ierr = (*f)(ts);CHKERRQ(ierr); 43582bf6240SBarry Smith } 43682bf6240SBarry Smith PetscFunctionReturn(0); 43782bf6240SBarry Smith } 43882bf6240SBarry Smith 43982bf6240SBarry Smith 4404a2ae208SSatish Balay #undef __FUNCT__ 4414a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 44282bf6240SBarry Smith /*@ 44382bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 44482bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 44582bf6240SBarry Smith 44615091d37SBarry Smith Collective on TS 44715091d37SBarry Smith 44882bf6240SBarry Smith Input Parameters: 44915091d37SBarry Smith + ts - timestep context 45082bf6240SBarry Smith . dt - function to compute timestep 45115091d37SBarry Smith - ctx - [optional] user-defined context for private data 45282bf6240SBarry Smith required by the function (may be PETSC_NULL) 45382bf6240SBarry Smith 45415091d37SBarry Smith Level: intermediate 455fee21e36SBarry Smith 45682bf6240SBarry Smith Calling sequence of func: 45787828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 45882bf6240SBarry Smith 45982bf6240SBarry Smith . newdt - the newly computed timestep 46082bf6240SBarry Smith . ctx - [optional] timestep context 46182bf6240SBarry Smith 46282bf6240SBarry Smith Notes: 46382bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 46482bf6240SBarry Smith during the timestepping process. 46582bf6240SBarry Smith 46682bf6240SBarry Smith .keywords: timestep, pseudo, set 46782bf6240SBarry Smith 46882bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 46982bf6240SBarry Smith @*/ 47087828ca2SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx) 47182bf6240SBarry Smith { 47287828ca2SBarry Smith int ierr,(*f)(TS,int (*)(TS,PetscReal *,void *),void *); 47382bf6240SBarry Smith 47482bf6240SBarry Smith PetscFunctionBegin; 47582bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 47682bf6240SBarry Smith 477c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 47882bf6240SBarry Smith if (f) { 47982bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 48082bf6240SBarry Smith } 48182bf6240SBarry Smith PetscFunctionReturn(0); 48282bf6240SBarry Smith } 48382bf6240SBarry Smith 48482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 48582bf6240SBarry Smith 486fb2e594dSBarry Smith EXTERN_C_BEGIN 4874a2ae208SSatish Balay #undef __FUNCT__ 4884a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 48987828ca2SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx) 49082bf6240SBarry Smith { 49182bf6240SBarry Smith TS_Pseudo *pseudo; 49282bf6240SBarry Smith 49382bf6240SBarry Smith PetscFunctionBegin; 49482bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 49582bf6240SBarry Smith pseudo->verify = dt; 49682bf6240SBarry Smith pseudo->verifyctx = ctx; 49782bf6240SBarry Smith PetscFunctionReturn(0); 49882bf6240SBarry Smith } 499fb2e594dSBarry Smith EXTERN_C_END 50082bf6240SBarry Smith 501fb2e594dSBarry Smith EXTERN_C_BEGIN 5024a2ae208SSatish Balay #undef __FUNCT__ 5034a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 50487828ca2SBarry Smith int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 50582bf6240SBarry Smith { 5064bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 50782bf6240SBarry Smith 50882bf6240SBarry Smith PetscFunctionBegin; 50982bf6240SBarry Smith pseudo->dt_increment = inc; 51082bf6240SBarry Smith PetscFunctionReturn(0); 51182bf6240SBarry Smith } 512fb2e594dSBarry Smith EXTERN_C_END 51382bf6240SBarry Smith 514fb2e594dSBarry Smith EXTERN_C_BEGIN 5154a2ae208SSatish Balay #undef __FUNCT__ 5164a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 51782bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 51882bf6240SBarry Smith { 5194bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 52082bf6240SBarry Smith 52182bf6240SBarry Smith PetscFunctionBegin; 5224bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 52382bf6240SBarry Smith PetscFunctionReturn(0); 52482bf6240SBarry Smith } 525fb2e594dSBarry Smith EXTERN_C_END 52682bf6240SBarry Smith 527fb2e594dSBarry Smith EXTERN_C_BEGIN 5284a2ae208SSatish Balay #undef __FUNCT__ 5294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 53087828ca2SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx) 53182bf6240SBarry Smith { 5324bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 53382bf6240SBarry Smith 53482bf6240SBarry Smith PetscFunctionBegin; 53582bf6240SBarry Smith pseudo->dt = dt; 53682bf6240SBarry Smith pseudo->dtctx = ctx; 53782bf6240SBarry Smith PetscFunctionReturn(0); 53882bf6240SBarry Smith } 539fb2e594dSBarry Smith EXTERN_C_END 54082bf6240SBarry Smith 54182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 54282bf6240SBarry Smith 543fb2e594dSBarry Smith EXTERN_C_BEGIN 5444a2ae208SSatish Balay #undef __FUNCT__ 5454a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 5467bf11e45SBarry Smith int TSCreate_Pseudo(TS ts) 5472d3f70b5SBarry Smith { 5487bf11e45SBarry Smith TS_Pseudo *pseudo; 5492d3f70b5SBarry Smith int ierr; 5502d3f70b5SBarry Smith 5513a40ed3dSBarry Smith PetscFunctionBegin; 552000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 553000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 5542d3f70b5SBarry Smith 5552d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 55629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 5572d3f70b5SBarry Smith } 5582d3f70b5SBarry Smith if (!ts->A) { 55929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian"); 5602d3f70b5SBarry Smith } 561273d9f13SBarry Smith 562000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 563000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 564000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 5657bf11e45SBarry Smith 5667bf11e45SBarry Smith /* create the required nonlinear solver context */ 567*4b27c08aSLois Curfman McInnes ierr = SNESCreate(ts->comm,&ts->snes);CHKERRQ(ierr); 5682d3f70b5SBarry Smith 569b0a32e0cSBarry Smith ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr); 570b0a32e0cSBarry Smith PetscLogObjectMemory(ts,sizeof(TS_Pseudo)); 571eed86810SBarry Smith 572549d3d68SSatish Balay ierr = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr); 5737bf11e45SBarry Smith ts->data = (void*)pseudo; 5742d3f70b5SBarry Smith 57528aa8177SBarry Smith pseudo->dt_increment = 1.1; 5764bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 57728aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 57882bf6240SBarry Smith 579f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 580e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 5810c97f09dSSatish Balay TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 582f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 583e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 5840c97f09dSSatish Balay TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 585f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 586e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 5870c97f09dSSatish Balay TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 5880c97f09dSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 5890c97f09dSSatish Balay "TSPseudoSetTimeStep_Pseudo", 5900c97f09dSSatish Balay TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 5913a40ed3dSBarry Smith PetscFunctionReturn(0); 5922d3f70b5SBarry Smith } 593fb2e594dSBarry Smith EXTERN_C_END 5942d3f70b5SBarry Smith 5954a2ae208SSatish Balay #undef __FUNCT__ 5964a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep" 59782bf6240SBarry Smith /*@C 59882bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 59982bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 60028aa8177SBarry Smith 60115091d37SBarry Smith Collective on TS 60215091d37SBarry Smith 60328aa8177SBarry Smith Input Parameters: 60428aa8177SBarry Smith . ts - the timestep context 60582bf6240SBarry Smith . dtctx - unused timestep context 60628aa8177SBarry Smith 60782bf6240SBarry Smith Output Parameter: 60882bf6240SBarry Smith . newdt - the timestep to use for the next step 60928aa8177SBarry Smith 61015091d37SBarry Smith Level: advanced 61115091d37SBarry Smith 61282bf6240SBarry Smith .keywords: timestep, pseudo, default 613564e8f4eSLois Curfman McInnes 61482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 61528aa8177SBarry Smith @*/ 61687828ca2SBarry Smith int TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 61728aa8177SBarry Smith { 61882bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 61987828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 62082bf6240SBarry Smith int ierr; 62128aa8177SBarry Smith 6223a40ed3dSBarry Smith PetscFunctionBegin; 62382bf6240SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 62482bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 62582bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 62682bf6240SBarry Smith /* first time through so compute initial function norm */ 62782bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 62882bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 62982bf6240SBarry Smith } 63082bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 63182bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 632004884caSBarry Smith } else if (pseudo->increment_dt_from_initial_dt) { 63382bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 63482bf6240SBarry Smith } else { 63582bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 63682bf6240SBarry Smith } 63782bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6383a40ed3dSBarry Smith PetscFunctionReturn(0); 63928aa8177SBarry Smith } 6402d3f70b5SBarry Smith 641004884caSBarry Smith 642004884caSBarry Smith 643004884caSBarry Smith 644004884caSBarry Smith 645004884caSBarry Smith 646004884caSBarry Smith 647004884caSBarry Smith 648004884caSBarry Smith 649004884caSBarry Smith 650004884caSBarry Smith 651004884caSBarry Smith 652004884caSBarry Smith 653004884caSBarry Smith 654004884caSBarry Smith 655004884caSBarry Smith 656004884caSBarry Smith 657