1*ea709b57SSatish Balay /*$Id: posindep.c,v 1.56 2001/08/06 21:18:14 bsmith Exp balay $*/ 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; 58b0a32e0cSBarry Smith ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 597bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 60b0a32e0cSBarry Smith 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); 1557bf11e45SBarry Smith current_time_step = ts->time_step; 1567bf11e45SBarry Smith while (1) { 1577bf11e45SBarry Smith ts->ptime += current_time_step; 1587bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its);CHKERRQ(ierr); 159f2267985SLois Curfman McInnes ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr); 160f2267985SLois Curfman McInnes ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits; 1617bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr); 1627bf11e45SBarry Smith if (ok) break; 1637bf11e45SBarry Smith ts->ptime -= current_time_step; 1647bf11e45SBarry Smith current_time_step = ts->time_step; 1657bf11e45SBarry Smith } 1667bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr); 1672d3f70b5SBarry Smith ts->steps++; 1682d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1692d3f70b5SBarry Smith } 1702d3f70b5SBarry Smith 1712d3f70b5SBarry Smith *steps += ts->steps; 172142b95e3SSatish Balay *ptime = ts->ptime; 1733a40ed3dSBarry Smith PetscFunctionReturn(0); 1742d3f70b5SBarry Smith } 1752d3f70b5SBarry Smith 1762d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1774a2ae208SSatish Balay #undef __FUNCT__ 1784a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo" 179e1311b90SBarry Smith static int TSDestroy_Pseudo(TS ts) 1802d3f70b5SBarry Smith { 1817bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 1822d3f70b5SBarry Smith int ierr; 1832d3f70b5SBarry Smith 1843a40ed3dSBarry Smith PetscFunctionBegin; 185e4ef1970SSatish Balay if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 1867bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1877bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 1882d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 189606d414cSSatish Balay ierr = PetscFree(pseudo);CHKERRQ(ierr); 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 1912d3f70b5SBarry Smith } 1922d3f70b5SBarry Smith 1932d3f70b5SBarry Smith 1942d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1952d3f70b5SBarry Smith /* 1962d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 1972d3f70b5SBarry Smith */ 1982d3f70b5SBarry Smith 1994a2ae208SSatish Balay #undef __FUNCT__ 2004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoMatMult" 2017bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 2022d3f70b5SBarry Smith { 2032d3f70b5SBarry Smith TS ts; 204*ea709b57SSatish Balay PetscScalar mdt,mone = -1.0; 2052d3f70b5SBarry Smith int ierr; 2062d3f70b5SBarry Smith 2073a40ed3dSBarry Smith PetscFunctionBegin; 2083a40ed3dSBarry Smith ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr); 2092d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2102d3f70b5SBarry Smith 2112d3f70b5SBarry Smith /* apply user provided function */ 2122d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y);CHKERRQ(ierr); 2132d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2142d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y);CHKERRQ(ierr); 2153a40ed3dSBarry Smith PetscFunctionReturn(0); 2162d3f70b5SBarry Smith } 2172d3f70b5SBarry Smith 2182d3f70b5SBarry Smith /* 2192d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2202d3f70b5SBarry Smith 2212d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2222d3f70b5SBarry Smith */ 2234a2ae208SSatish Balay #undef __FUNCT__ 2244a2ae208SSatish Balay #define __FUNCT__ "TSPseudoFunction" 2257bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2262d3f70b5SBarry Smith { 2272d3f70b5SBarry Smith TS ts = (TS) ctx; 228*ea709b57SSatish Balay PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2292d3f70b5SBarry Smith int ierr,i,n; 2302d3f70b5SBarry Smith 2313a40ed3dSBarry Smith PetscFunctionBegin; 2322d3f70b5SBarry Smith /* apply user provided function */ 2332d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr); 2347bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2352d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr); 2362d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1);CHKERRQ(ierr); 2372d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr); 2382d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 2392d3f70b5SBarry Smith for (i=0; i<n; i++) { 2402d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2412d3f70b5SBarry Smith } 242888f2ed8SSatish Balay ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr); 243888f2ed8SSatish Balay ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr); 244888f2ed8SSatish Balay ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr); 2452d3f70b5SBarry Smith 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2472d3f70b5SBarry Smith } 2482d3f70b5SBarry Smith 2492d3f70b5SBarry Smith /* 2502d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2512d3f70b5SBarry Smith 2522d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2532d3f70b5SBarry Smith */ 2544a2ae208SSatish Balay #undef __FUNCT__ 2554a2ae208SSatish Balay #define __FUNCT__ "TSPseudoJacobian" 2567bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2572d3f70b5SBarry Smith { 2582d3f70b5SBarry Smith TS ts = (TS) ctx; 2592d3f70b5SBarry Smith int ierr; 260*ea709b57SSatish Balay PetscScalar mone = -1.0,mdt = 1.0/ts->time_step; 261273d9f13SBarry Smith PetscTruth isshell; 2622d3f70b5SBarry Smith 2633a40ed3dSBarry Smith PetscFunctionBegin; 2642d3f70b5SBarry Smith /* construct users Jacobian */ 265a7a1495cSBarry Smith ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr); 2662d3f70b5SBarry Smith 2672d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 268273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)*AA,MATSHELL,&isshell);CHKERRQ(ierr); 269273d9f13SBarry Smith if (!isshell) { 2702d3f70b5SBarry Smith ierr = MatScale(&mone,*AA);CHKERRQ(ierr); 2712d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA);CHKERRQ(ierr); 2722d3f70b5SBarry Smith } 273273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)*BB,MATSHELL,&isshell);CHKERRQ(ierr); 274273d9f13SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && !isshell) { 2752d3f70b5SBarry Smith ierr = MatScale(&mone,*BB);CHKERRQ(ierr); 2762d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB);CHKERRQ(ierr); 2772d3f70b5SBarry Smith } 2782d3f70b5SBarry Smith 2793a40ed3dSBarry Smith PetscFunctionReturn(0); 2802d3f70b5SBarry Smith } 2812d3f70b5SBarry Smith 2822d3f70b5SBarry Smith 2834a2ae208SSatish Balay #undef __FUNCT__ 2844a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 2857bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 2862d3f70b5SBarry Smith { 2877bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 2882d3f70b5SBarry Smith int ierr,M,m; 2892d3f70b5SBarry Smith 2903a40ed3dSBarry Smith PetscFunctionBegin; 291b48991e4SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 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); 29937bd1cefSSatish Balay 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 3064a2ae208SSatish Balay #undef __FUNCT__ 3074a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultMonitor" 30887828ca2SBarry Smith int TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx) 3092d3f70b5SBarry Smith { 3107bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 311d132466eSBarry Smith int ierr; 3122d3f70b5SBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 314142b95e3SSatish Balay ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 3153a40ed3dSBarry Smith PetscFunctionReturn(0); 3162d3f70b5SBarry Smith } 3172d3f70b5SBarry Smith 3184a2ae208SSatish Balay #undef __FUNCT__ 3194a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 3207bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3212d3f70b5SBarry Smith { 3224bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 323f1af5d2fSBarry Smith int ierr; 324f1af5d2fSBarry Smith PetscTruth flg; 3252d3f70b5SBarry Smith 3263a40ed3dSBarry Smith PetscFunctionBegin; 3272d3f70b5SBarry Smith 328b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 329b0a32e0cSBarry Smith ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr); 3302d3f70b5SBarry Smith if (flg) { 331329f5518SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 33228aa8177SBarry Smith } 333b0a32e0cSBarry Smith ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr); 334ca4b7087SBarry Smith if (flg) { 335ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 336ca4b7087SBarry Smith } 33787828ca2SBarry Smith ierr = PetscOptionsPetscReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 338b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3392d3f70b5SBarry Smith 3403a40ed3dSBarry Smith PetscFunctionReturn(0); 3412d3f70b5SBarry Smith } 3422d3f70b5SBarry Smith 3434a2ae208SSatish Balay #undef __FUNCT__ 3444a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 345b0a32e0cSBarry Smith static int TSView_Pseudo(TS ts,PetscViewer viewer) 3462d3f70b5SBarry Smith { 3473a40ed3dSBarry Smith PetscFunctionBegin; 3483a40ed3dSBarry Smith PetscFunctionReturn(0); 3492d3f70b5SBarry Smith } 3502d3f70b5SBarry Smith 35182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3524a2ae208SSatish Balay #undef __FUNCT__ 3534a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 35482bf6240SBarry Smith /*@ 35582bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 35682bf6240SBarry Smith last timestep. 35782bf6240SBarry Smith 35815091d37SBarry Smith Collective on TS 35915091d37SBarry Smith 36082bf6240SBarry Smith Input Parameters: 36115091d37SBarry Smith + ts - timestep context 36282bf6240SBarry Smith . dt - user-defined function to verify timestep 36315091d37SBarry Smith - ctx - [optional] user-defined context for private data 36482bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 36582bf6240SBarry Smith 36615091d37SBarry Smith Level: advanced 367fee21e36SBarry Smith 36882bf6240SBarry Smith Calling sequence of func: 36987828ca2SBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag); 37082bf6240SBarry Smith 37182bf6240SBarry Smith . update - latest solution vector 37282bf6240SBarry Smith . ctx - [optional] timestep context 37382bf6240SBarry Smith . newdt - the timestep to use for the next step 37482bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 37582bf6240SBarry Smith 37682bf6240SBarry Smith Notes: 37782bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 37882bf6240SBarry Smith during the timestepping process. 37982bf6240SBarry Smith 38082bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 38182bf6240SBarry Smith 38282bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 38382bf6240SBarry Smith @*/ 38487828ca2SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx) 38582bf6240SBarry Smith { 38687828ca2SBarry Smith int ierr,(*f)(TS,int (*)(TS,Vec,void*,PetscReal *,int *),void *); 38782bf6240SBarry Smith 38882bf6240SBarry Smith PetscFunctionBegin; 38982bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 39082bf6240SBarry Smith 3911713a123SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)())&f);CHKERRQ(ierr); 39282bf6240SBarry Smith if (f) { 39382bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 39482bf6240SBarry Smith } 39582bf6240SBarry Smith PetscFunctionReturn(0); 39682bf6240SBarry Smith } 39782bf6240SBarry Smith 3984a2ae208SSatish Balay #undef __FUNCT__ 3994a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 40082bf6240SBarry Smith /*@ 40182bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 40282bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 40382bf6240SBarry Smith 404fee21e36SBarry Smith Collective on TS 405fee21e36SBarry Smith 40615091d37SBarry Smith Input Parameters: 40715091d37SBarry Smith + ts - the timestep context 40815091d37SBarry Smith - inc - the scaling factor >= 1.0 40915091d37SBarry Smith 41082bf6240SBarry Smith Options Database Key: 41182bf6240SBarry Smith $ -ts_pseudo_increment <increment> 41282bf6240SBarry Smith 41315091d37SBarry Smith Level: advanced 41415091d37SBarry Smith 41582bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 41682bf6240SBarry Smith 41782bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 41882bf6240SBarry Smith @*/ 41987828ca2SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 42082bf6240SBarry Smith { 42187828ca2SBarry Smith int ierr,(*f)(TS,PetscReal); 42282bf6240SBarry Smith 42382bf6240SBarry Smith PetscFunctionBegin; 42482bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 42582bf6240SBarry Smith 4261713a123SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)())&f);CHKERRQ(ierr); 42782bf6240SBarry Smith if (f) { 42882bf6240SBarry Smith ierr = (*f)(ts,inc);CHKERRQ(ierr); 42982bf6240SBarry Smith } 43082bf6240SBarry Smith PetscFunctionReturn(0); 43182bf6240SBarry Smith } 43282bf6240SBarry Smith 4334a2ae208SSatish Balay #undef __FUNCT__ 4344a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 43582bf6240SBarry Smith /*@ 43682bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 43782bf6240SBarry Smith is computed via the formula 43882bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 43982bf6240SBarry Smith rather than the default update, 44082bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 44182bf6240SBarry Smith 44215091d37SBarry Smith Collective on TS 44315091d37SBarry Smith 44482bf6240SBarry Smith Input Parameter: 44582bf6240SBarry Smith . ts - the timestep context 44682bf6240SBarry Smith 44782bf6240SBarry Smith Options Database Key: 44882bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 44982bf6240SBarry Smith 45015091d37SBarry Smith Level: advanced 45115091d37SBarry Smith 45282bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 45382bf6240SBarry Smith 45482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 45582bf6240SBarry Smith @*/ 45682bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts) 45782bf6240SBarry Smith { 45882bf6240SBarry Smith int ierr,(*f)(TS); 45982bf6240SBarry Smith 46082bf6240SBarry Smith PetscFunctionBegin; 46182bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 46282bf6240SBarry Smith 4631713a123SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)())&f);CHKERRQ(ierr); 46482bf6240SBarry Smith if (f) { 46582bf6240SBarry Smith ierr = (*f)(ts);CHKERRQ(ierr); 46682bf6240SBarry Smith } 46782bf6240SBarry Smith PetscFunctionReturn(0); 46882bf6240SBarry Smith } 46982bf6240SBarry Smith 47082bf6240SBarry Smith 4714a2ae208SSatish Balay #undef __FUNCT__ 4724a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 47382bf6240SBarry Smith /*@ 47482bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 47582bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 47682bf6240SBarry Smith 47715091d37SBarry Smith Collective on TS 47815091d37SBarry Smith 47982bf6240SBarry Smith Input Parameters: 48015091d37SBarry Smith + ts - timestep context 48182bf6240SBarry Smith . dt - function to compute timestep 48215091d37SBarry Smith - ctx - [optional] user-defined context for private data 48382bf6240SBarry Smith required by the function (may be PETSC_NULL) 48482bf6240SBarry Smith 48515091d37SBarry Smith Level: intermediate 486fee21e36SBarry Smith 48782bf6240SBarry Smith Calling sequence of func: 48887828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 48982bf6240SBarry Smith 49082bf6240SBarry Smith . newdt - the newly computed timestep 49182bf6240SBarry Smith . ctx - [optional] timestep context 49282bf6240SBarry Smith 49382bf6240SBarry Smith Notes: 49482bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 49582bf6240SBarry Smith during the timestepping process. 49682bf6240SBarry Smith 49782bf6240SBarry Smith .keywords: timestep, pseudo, set 49882bf6240SBarry Smith 49982bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 50082bf6240SBarry Smith @*/ 50187828ca2SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx) 50282bf6240SBarry Smith { 50387828ca2SBarry Smith int ierr,(*f)(TS,int (*)(TS,PetscReal *,void *),void *); 50482bf6240SBarry Smith 50582bf6240SBarry Smith PetscFunctionBegin; 50682bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 50782bf6240SBarry Smith 5081713a123SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)())&f);CHKERRQ(ierr); 50982bf6240SBarry Smith if (f) { 51082bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 51182bf6240SBarry Smith } 51282bf6240SBarry Smith PetscFunctionReturn(0); 51382bf6240SBarry Smith } 51482bf6240SBarry Smith 51582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 51682bf6240SBarry Smith 517fb2e594dSBarry Smith EXTERN_C_BEGIN 5184a2ae208SSatish Balay #undef __FUNCT__ 5194a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 52087828ca2SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx) 52182bf6240SBarry Smith { 52282bf6240SBarry Smith TS_Pseudo *pseudo; 52382bf6240SBarry Smith 52482bf6240SBarry Smith PetscFunctionBegin; 52582bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 52682bf6240SBarry Smith pseudo->verify = dt; 52782bf6240SBarry Smith pseudo->verifyctx = ctx; 52882bf6240SBarry Smith PetscFunctionReturn(0); 52982bf6240SBarry Smith } 530fb2e594dSBarry Smith EXTERN_C_END 53182bf6240SBarry Smith 532fb2e594dSBarry Smith EXTERN_C_BEGIN 5334a2ae208SSatish Balay #undef __FUNCT__ 5344a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 53587828ca2SBarry Smith int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 53682bf6240SBarry Smith { 5374bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 53882bf6240SBarry Smith 53982bf6240SBarry Smith PetscFunctionBegin; 54082bf6240SBarry Smith pseudo->dt_increment = inc; 54182bf6240SBarry Smith PetscFunctionReturn(0); 54282bf6240SBarry Smith } 543fb2e594dSBarry Smith EXTERN_C_END 54482bf6240SBarry Smith 545fb2e594dSBarry Smith EXTERN_C_BEGIN 5464a2ae208SSatish Balay #undef __FUNCT__ 5474a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 54882bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 54982bf6240SBarry Smith { 5504bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 55182bf6240SBarry Smith 55282bf6240SBarry Smith PetscFunctionBegin; 5534bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 55482bf6240SBarry Smith PetscFunctionReturn(0); 55582bf6240SBarry Smith } 556fb2e594dSBarry Smith EXTERN_C_END 55782bf6240SBarry Smith 558fb2e594dSBarry Smith EXTERN_C_BEGIN 5594a2ae208SSatish Balay #undef __FUNCT__ 5604a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 56187828ca2SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx) 56282bf6240SBarry Smith { 5634bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56482bf6240SBarry Smith 56582bf6240SBarry Smith PetscFunctionBegin; 56682bf6240SBarry Smith pseudo->dt = dt; 56782bf6240SBarry Smith pseudo->dtctx = ctx; 56882bf6240SBarry Smith PetscFunctionReturn(0); 56982bf6240SBarry Smith } 570fb2e594dSBarry Smith EXTERN_C_END 57182bf6240SBarry Smith 57282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 57382bf6240SBarry Smith 574fb2e594dSBarry Smith EXTERN_C_BEGIN 5754a2ae208SSatish Balay #undef __FUNCT__ 5764a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 5777bf11e45SBarry Smith int TSCreate_Pseudo(TS ts) 5782d3f70b5SBarry Smith { 5797bf11e45SBarry Smith TS_Pseudo *pseudo; 5802d3f70b5SBarry Smith int ierr; 581273d9f13SBarry Smith PetscTruth isshell; 5822d3f70b5SBarry Smith 5833a40ed3dSBarry Smith PetscFunctionBegin; 5847bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 5857bf11e45SBarry Smith ts->view = TSView_Pseudo; 5862d3f70b5SBarry Smith 5872d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 58829bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 5892d3f70b5SBarry Smith } 5902d3f70b5SBarry Smith if (!ts->A) { 59129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian"); 5922d3f70b5SBarry Smith } 593273d9f13SBarry Smith 594273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)ts->A,MATSHELL,&isshell);CHKERRQ(ierr); 595273d9f13SBarry Smith if (isshell) { 5962d3f70b5SBarry Smith ts->Ashell = ts->A; 5972d3f70b5SBarry Smith } 5987bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 5997bf11e45SBarry Smith ts->step = TSStep_Pseudo; 6007bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 6017bf11e45SBarry Smith 6027bf11e45SBarry Smith /* create the required nonlinear solver context */ 6032d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 6042d3f70b5SBarry Smith 605b0a32e0cSBarry Smith ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr); 606b0a32e0cSBarry Smith PetscLogObjectMemory(ts,sizeof(TS_Pseudo)); 607eed86810SBarry Smith 608549d3d68SSatish Balay ierr = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr); 6097bf11e45SBarry Smith ts->data = (void*)pseudo; 6102d3f70b5SBarry Smith 61128aa8177SBarry Smith pseudo->dt_increment = 1.1; 6124bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 61328aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 61482bf6240SBarry Smith 615f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 616e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 6170c97f09dSSatish Balay TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 618f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 619e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 6200c97f09dSSatish Balay TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 621f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 622e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 6230c97f09dSSatish Balay TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 6240c97f09dSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 6250c97f09dSSatish Balay "TSPseudoSetTimeStep_Pseudo", 6260c97f09dSSatish Balay TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6273a40ed3dSBarry Smith PetscFunctionReturn(0); 6282d3f70b5SBarry Smith } 629fb2e594dSBarry Smith EXTERN_C_END 6302d3f70b5SBarry Smith 6314a2ae208SSatish Balay #undef __FUNCT__ 6324a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep" 63382bf6240SBarry Smith /*@C 63482bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 63582bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 63628aa8177SBarry Smith 63715091d37SBarry Smith Collective on TS 63815091d37SBarry Smith 63928aa8177SBarry Smith Input Parameters: 64028aa8177SBarry Smith . ts - the timestep context 64182bf6240SBarry Smith . dtctx - unused timestep context 64228aa8177SBarry Smith 64382bf6240SBarry Smith Output Parameter: 64482bf6240SBarry Smith . newdt - the timestep to use for the next step 64528aa8177SBarry Smith 64615091d37SBarry Smith Level: advanced 64715091d37SBarry Smith 64882bf6240SBarry Smith .keywords: timestep, pseudo, default 649564e8f4eSLois Curfman McInnes 65082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 65128aa8177SBarry Smith @*/ 65287828ca2SBarry Smith int TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 65328aa8177SBarry Smith { 65482bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 65587828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 65682bf6240SBarry Smith int ierr; 65728aa8177SBarry Smith 6583a40ed3dSBarry Smith PetscFunctionBegin; 65982bf6240SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 66082bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 66182bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 66282bf6240SBarry Smith /* first time through so compute initial function norm */ 66382bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 66482bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 66582bf6240SBarry Smith } 66682bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 66782bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 668004884caSBarry Smith } else if (pseudo->increment_dt_from_initial_dt) { 66982bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 67082bf6240SBarry Smith } else { 67182bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 67282bf6240SBarry Smith } 67382bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6743a40ed3dSBarry Smith PetscFunctionReturn(0); 67528aa8177SBarry Smith } 6762d3f70b5SBarry Smith 677004884caSBarry Smith 678004884caSBarry Smith 679004884caSBarry Smith 680004884caSBarry Smith 681004884caSBarry Smith 682004884caSBarry Smith 683004884caSBarry Smith 684004884caSBarry Smith 685004884caSBarry Smith 686004884caSBarry Smith 687004884caSBarry Smith 688004884caSBarry Smith 689004884caSBarry Smith 690004884caSBarry Smith 691004884caSBarry Smith 692004884caSBarry Smith 693