163dd3a1aSKris Buschelman #define PETSCTS_DLL 263dd3a1aSKris Buschelman 32d3f70b5SBarry Smith /* 4fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 52d3f70b5SBarry Smith */ 6*7c4f633dSBarry Smith #include "private/tsimpl.h" /*I "petscts.h" I*/ 72d3f70b5SBarry Smith 82d3f70b5SBarry Smith typedef struct { 92d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 102d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 112d3f70b5SBarry Smith Vec rhs; /* work vector for RHS; vec_sol/dt */ 122d3f70b5SBarry Smith 132d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 142d3f70b5SBarry Smith 156849ba73SBarry Smith PetscErrorCode (*dt)(TS,PetscReal*,void*); /* compute next timestep, and related context */ 162d3f70b5SBarry Smith void *dtctx; 17a7cc72afSBarry Smith PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscTruth*); /* verify previous timestep and related context */ 187bf11e45SBarry Smith void *verifyctx; 192d3f70b5SBarry Smith 2087828ca2SBarry Smith PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */ 2187828ca2SBarry Smith PetscReal fnorm_previous; 2228aa8177SBarry Smith 2387828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 244bbc92c1SBarry Smith PetscTruth increment_dt_from_initial_dt; 257bf11e45SBarry Smith } TS_Pseudo; 262d3f70b5SBarry Smith 272d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 282d3f70b5SBarry Smith 294a2ae208SSatish Balay #undef __FUNCT__ 304a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep" 317bf11e45SBarry Smith /*@ 327bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 33564e8f4eSLois Curfman McInnes pseudo-timestepping process. 342d3f70b5SBarry Smith 3515091d37SBarry Smith Collective on TS 3615091d37SBarry Smith 377bf11e45SBarry Smith Input Parameter: 387bf11e45SBarry Smith . ts - timestep context 397bf11e45SBarry Smith 407bf11e45SBarry Smith Output Parameter: 41fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 42fb4a63b6SLois Curfman McInnes 4315091d37SBarry Smith Level: advanced 44564e8f4eSLois Curfman McInnes 45564e8f4eSLois Curfman McInnes Notes: 46564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 47564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 48564e8f4eSLois Curfman McInnes 49fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 50564e8f4eSLois Curfman McInnes 51564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 527bf11e45SBarry Smith @*/ 5363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 547bf11e45SBarry Smith { 557bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56dfbe8321SBarry Smith PetscErrorCode ierr; 577bf11e45SBarry Smith 583a40ed3dSBarry Smith PetscFunctionBegin; 59d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 607bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 61d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 623a40ed3dSBarry Smith PetscFunctionReturn(0); 637bf11e45SBarry Smith } 647bf11e45SBarry Smith 657bf11e45SBarry Smith 667bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 674a2ae208SSatish Balay #undef __FUNCT__ 684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 697bf11e45SBarry Smith /*@C 70639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 717bf11e45SBarry Smith 7215091d37SBarry Smith Collective on TS 7315091d37SBarry Smith 747bf11e45SBarry Smith Input Parameters: 7515091d37SBarry Smith + ts - the timestep context 767bf11e45SBarry Smith . dtctx - unused timestep context 7715091d37SBarry Smith - update - latest solution vector 787bf11e45SBarry Smith 79564e8f4eSLois Curfman McInnes Output Parameters: 8015091d37SBarry Smith + newdt - the timestep to use for the next step 8115091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 827bf11e45SBarry Smith 8315091d37SBarry Smith Level: advanced 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 @*/ 9363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscTruth *flag) 947bf11e45SBarry Smith { 953a40ed3dSBarry Smith PetscFunctionBegin; 96a7cc72afSBarry Smith *flag = PETSC_TRUE; 973a40ed3dSBarry Smith PetscFunctionReturn(0); 987bf11e45SBarry Smith } 997bf11e45SBarry Smith 1007bf11e45SBarry Smith 1014a2ae208SSatish Balay #undef __FUNCT__ 1024a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep" 1037bf11e45SBarry Smith /*@ 104564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1057bf11e45SBarry Smith 10615091d37SBarry Smith Collective on TS 10715091d37SBarry Smith 108fb4a63b6SLois Curfman McInnes Input Parameters: 10915091d37SBarry Smith + ts - timestep context 11015091d37SBarry Smith - update - latest solution vector 1117bf11e45SBarry Smith 112fb4a63b6SLois Curfman McInnes Output Parameters: 11315091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11415091d37SBarry Smith - flag - indicates if current timestep was ok 1157bf11e45SBarry Smith 11615091d37SBarry Smith Level: advanced 117fee21e36SBarry Smith 118564e8f4eSLois Curfman McInnes Notes: 119564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 120564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 121564e8f4eSLois Curfman McInnes 122564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 123564e8f4eSLois Curfman McInnes 124564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1257bf11e45SBarry Smith @*/ 12663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscTruth *flag) 1277bf11e45SBarry Smith { 1287bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 129dfbe8321SBarry Smith PetscErrorCode ierr; 1307bf11e45SBarry Smith 1313a40ed3dSBarry Smith PetscFunctionBegin; 132a7cc72afSBarry Smith if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);} 1337bf11e45SBarry Smith 1347bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 1357bf11e45SBarry Smith 1363a40ed3dSBarry Smith PetscFunctionReturn(0); 1377bf11e45SBarry Smith } 1387bf11e45SBarry Smith 1397bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1407bf11e45SBarry Smith 1414a2ae208SSatish Balay #undef __FUNCT__ 1424a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo" 143a7cc72afSBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime) 1442d3f70b5SBarry Smith { 1452d3f70b5SBarry Smith Vec sol = ts->vec_sol; 146dfbe8321SBarry Smith PetscErrorCode ierr; 147a7cc72afSBarry Smith PetscInt i,max_steps = ts->max_steps,its,lits; 148a7cc72afSBarry Smith PetscTruth ok; 1497bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 15087828ca2SBarry Smith PetscReal current_time_step; 1512d3f70b5SBarry Smith 1523a40ed3dSBarry Smith PetscFunctionBegin; 1532d3f70b5SBarry Smith *steps = -ts->steps; 1542d3f70b5SBarry Smith 1557bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr); 1562d3f70b5SBarry Smith for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) { 1577bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 158e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1597bf11e45SBarry Smith current_time_step = ts->time_step; 160e6e5fe25SBarry Smith while (PETSC_TRUE) { 1617bf11e45SBarry Smith ts->ptime += current_time_step; 162f69a0ea3SMatthew Knepley ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr); 163b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1649f954729SBarry Smith ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 165c9780b6fSBarry Smith ts->nonlinear_its += its; ts->linear_its += lits; 1667bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr); 1677bf11e45SBarry Smith if (ok) break; 1687bf11e45SBarry Smith ts->ptime -= current_time_step; 1697bf11e45SBarry Smith current_time_step = ts->time_step; 1707bf11e45SBarry Smith } 1717bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr); 1722d3f70b5SBarry Smith ts->steps++; 1732d3f70b5SBarry Smith } 174e6e5fe25SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 175e6e5fe25SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 176e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1772d3f70b5SBarry Smith 1782d3f70b5SBarry Smith *steps += ts->steps; 179142b95e3SSatish Balay *ptime = ts->ptime; 1803a40ed3dSBarry Smith PetscFunctionReturn(0); 1812d3f70b5SBarry Smith } 1822d3f70b5SBarry Smith 1832d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1844a2ae208SSatish Balay #undef __FUNCT__ 1854a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo" 1866849ba73SBarry Smith static PetscErrorCode TSDestroy_Pseudo(TS ts) 1872d3f70b5SBarry Smith { 1887bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 189dfbe8321SBarry Smith PetscErrorCode ierr; 1902d3f70b5SBarry Smith 1913a40ed3dSBarry Smith PetscFunctionBegin; 192e4ef1970SSatish Balay if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 1937bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1947bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 195606d414cSSatish Balay ierr = PetscFree(pseudo);CHKERRQ(ierr); 1963a40ed3dSBarry Smith PetscFunctionReturn(0); 1972d3f70b5SBarry Smith } 1982d3f70b5SBarry Smith 1992d3f70b5SBarry Smith 2002d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2012d3f70b5SBarry Smith 2022d3f70b5SBarry Smith /* 2032d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2042d3f70b5SBarry Smith 2052d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2062d3f70b5SBarry Smith */ 2074a2ae208SSatish Balay #undef __FUNCT__ 2084a2ae208SSatish Balay #define __FUNCT__ "TSPseudoFunction" 209dfbe8321SBarry Smith PetscErrorCode TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2102d3f70b5SBarry Smith { 2112d3f70b5SBarry Smith TS ts = (TS) ctx; 212ea709b57SSatish Balay PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 213dfbe8321SBarry Smith PetscErrorCode ierr; 214a7cc72afSBarry Smith PetscInt i,n; 2152d3f70b5SBarry Smith 2163a40ed3dSBarry Smith PetscFunctionBegin; 2172d3f70b5SBarry Smith /* apply user provided function */ 2182d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr); 2197bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2202d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr); 2212d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1);CHKERRQ(ierr); 2222d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr); 2232d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 2242d3f70b5SBarry Smith for (i=0; i<n; i++) { 2252d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2262d3f70b5SBarry Smith } 227888f2ed8SSatish Balay ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr); 228888f2ed8SSatish Balay ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr); 229888f2ed8SSatish Balay ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr); 2303a40ed3dSBarry Smith PetscFunctionReturn(0); 2312d3f70b5SBarry Smith } 2322d3f70b5SBarry Smith 2332d3f70b5SBarry Smith /* 2342d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2352d3f70b5SBarry Smith 2362d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2372d3f70b5SBarry Smith */ 2384a2ae208SSatish Balay #undef __FUNCT__ 2394a2ae208SSatish Balay #define __FUNCT__ "TSPseudoJacobian" 240dfbe8321SBarry Smith PetscErrorCode TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2412d3f70b5SBarry Smith { 2422d3f70b5SBarry Smith TS ts = (TS) ctx; 243dfbe8321SBarry Smith PetscErrorCode ierr; 2442d3f70b5SBarry Smith 2453a40ed3dSBarry Smith PetscFunctionBegin; 2462d3f70b5SBarry Smith /* construct users Jacobian */ 247a7a1495cSBarry Smith ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr); 2482d3f70b5SBarry Smith 249614f654bSBarry Smith /* shift and scale Jacobian */ 250b0f6e734SBarry Smith ierr = TSScaleShiftMatrices(ts,*AA,*BB,*str);CHKERRQ(ierr); 2513a40ed3dSBarry Smith PetscFunctionReturn(0); 2522d3f70b5SBarry Smith } 2532d3f70b5SBarry Smith 2542d3f70b5SBarry Smith 2554a2ae208SSatish Balay #undef __FUNCT__ 2564a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 2576849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 2582d3f70b5SBarry Smith { 2597bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 260dfbe8321SBarry Smith PetscErrorCode ierr; 2612d3f70b5SBarry Smith 2623a40ed3dSBarry Smith PetscFunctionBegin; 2637bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 2647bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 2657bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 2668beb423aSHong Zhang ierr = SNESSetJacobian(ts->snes,ts->Arhs,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 2682d3f70b5SBarry Smith } 2692d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2702d3f70b5SBarry Smith 2714a2ae208SSatish Balay #undef __FUNCT__ 272a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault" 273a6570f20SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 2742d3f70b5SBarry Smith { 2757bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 276dfbe8321SBarry Smith PetscErrorCode ierr; 277a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 2782d3f70b5SBarry Smith 2793a40ed3dSBarry Smith PetscFunctionBegin; 280a34d58ebSBarry Smith if (!ctx) { 2817adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 282a34d58ebSBarry Smith } 283a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 284a34d58ebSBarry Smith if (!ctx) { 285a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 286a34d58ebSBarry Smith } 2873a40ed3dSBarry Smith PetscFunctionReturn(0); 2882d3f70b5SBarry Smith } 2892d3f70b5SBarry Smith 2904a2ae208SSatish Balay #undef __FUNCT__ 2914a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 2926849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 2932d3f70b5SBarry Smith { 2944bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 295dfbe8321SBarry Smith PetscErrorCode ierr; 296f1af5d2fSBarry Smith PetscTruth flg; 297a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer; 2982d3f70b5SBarry Smith 2993a40ed3dSBarry Smith PetscFunctionBegin; 300b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 301a6570f20SBarry Smith ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",&flg);CHKERRQ(ierr); 3022d3f70b5SBarry Smith if (flg) { 3037adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 304a34d58ebSBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 30528aa8177SBarry Smith } 306b0a32e0cSBarry Smith ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr); 307ca4b7087SBarry Smith if (flg) { 308ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 309ca4b7087SBarry Smith } 310419fbf4bSSatish Balay ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 311b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3123a40ed3dSBarry Smith PetscFunctionReturn(0); 3132d3f70b5SBarry Smith } 3142d3f70b5SBarry Smith 3154a2ae208SSatish Balay #undef __FUNCT__ 3164a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 3176849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3182d3f70b5SBarry Smith { 3193a40ed3dSBarry Smith PetscFunctionBegin; 3203a40ed3dSBarry Smith PetscFunctionReturn(0); 3212d3f70b5SBarry Smith } 3222d3f70b5SBarry Smith 32382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3244a2ae208SSatish Balay #undef __FUNCT__ 3254a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 326ac226902SBarry Smith /*@C 32782bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 32882bf6240SBarry Smith last timestep. 32982bf6240SBarry Smith 33015091d37SBarry Smith Collective on TS 33115091d37SBarry Smith 33282bf6240SBarry Smith Input Parameters: 33315091d37SBarry Smith + ts - timestep context 33482bf6240SBarry Smith . dt - user-defined function to verify timestep 33515091d37SBarry Smith - ctx - [optional] user-defined context for private data 33682bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 33782bf6240SBarry Smith 33815091d37SBarry Smith Level: advanced 339fee21e36SBarry Smith 34082bf6240SBarry Smith Calling sequence of func: 341a7cc72afSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscTruth *flag); 34282bf6240SBarry Smith 34382bf6240SBarry Smith . update - latest solution vector 34482bf6240SBarry Smith . ctx - [optional] timestep context 34582bf6240SBarry Smith . newdt - the timestep to use for the next step 34682bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 34782bf6240SBarry Smith 34882bf6240SBarry Smith Notes: 34982bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 35082bf6240SBarry Smith during the timestepping process. 35182bf6240SBarry Smith 35282bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 35382bf6240SBarry Smith 35482bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 35582bf6240SBarry Smith @*/ 35663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscTruth*),void* ctx) 35782bf6240SBarry Smith { 358a7cc72afSBarry Smith PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscTruth *),void *); 35982bf6240SBarry Smith 36082bf6240SBarry Smith PetscFunctionBegin; 3614482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 36282bf6240SBarry Smith 363c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 36482bf6240SBarry Smith if (f) { 36582bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 36682bf6240SBarry Smith } 36782bf6240SBarry Smith PetscFunctionReturn(0); 36882bf6240SBarry Smith } 36982bf6240SBarry Smith 3704a2ae208SSatish Balay #undef __FUNCT__ 3714a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 37282bf6240SBarry Smith /*@ 37382bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 37482bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 37582bf6240SBarry Smith 376fee21e36SBarry Smith Collective on TS 377fee21e36SBarry Smith 37815091d37SBarry Smith Input Parameters: 37915091d37SBarry Smith + ts - the timestep context 38015091d37SBarry Smith - inc - the scaling factor >= 1.0 38115091d37SBarry Smith 38282bf6240SBarry Smith Options Database Key: 38382bf6240SBarry Smith $ -ts_pseudo_increment <increment> 38482bf6240SBarry Smith 38515091d37SBarry Smith Level: advanced 38615091d37SBarry Smith 38782bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 38882bf6240SBarry Smith 38982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 39082bf6240SBarry Smith @*/ 39163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 39282bf6240SBarry Smith { 393dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(TS,PetscReal); 39482bf6240SBarry Smith 39582bf6240SBarry Smith PetscFunctionBegin; 3964482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 39782bf6240SBarry Smith 398c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr); 39982bf6240SBarry Smith if (f) { 40082bf6240SBarry Smith ierr = (*f)(ts,inc);CHKERRQ(ierr); 40182bf6240SBarry Smith } 40282bf6240SBarry Smith PetscFunctionReturn(0); 40382bf6240SBarry Smith } 40482bf6240SBarry Smith 4054a2ae208SSatish Balay #undef __FUNCT__ 4064a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 40782bf6240SBarry Smith /*@ 40882bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 40982bf6240SBarry Smith is computed via the formula 41082bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 41182bf6240SBarry Smith rather than the default update, 41282bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 41382bf6240SBarry Smith 41415091d37SBarry Smith Collective on TS 41515091d37SBarry Smith 41682bf6240SBarry Smith Input Parameter: 41782bf6240SBarry Smith . ts - the timestep context 41882bf6240SBarry Smith 41982bf6240SBarry Smith Options Database Key: 42082bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 42182bf6240SBarry Smith 42215091d37SBarry Smith Level: advanced 42315091d37SBarry Smith 42482bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 42582bf6240SBarry Smith 42682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 42782bf6240SBarry Smith @*/ 42863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt(TS ts) 42982bf6240SBarry Smith { 430dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(TS); 43182bf6240SBarry Smith 43282bf6240SBarry Smith PetscFunctionBegin; 4334482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 43482bf6240SBarry Smith 435c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr); 43682bf6240SBarry Smith if (f) { 43782bf6240SBarry Smith ierr = (*f)(ts);CHKERRQ(ierr); 43882bf6240SBarry Smith } 43982bf6240SBarry Smith PetscFunctionReturn(0); 44082bf6240SBarry Smith } 44182bf6240SBarry Smith 44282bf6240SBarry Smith 4434a2ae208SSatish Balay #undef __FUNCT__ 4444a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 445ac226902SBarry Smith /*@C 44682bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 44782bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 44882bf6240SBarry Smith 44915091d37SBarry Smith Collective on TS 45015091d37SBarry Smith 45182bf6240SBarry Smith Input Parameters: 45215091d37SBarry Smith + ts - timestep context 45382bf6240SBarry Smith . dt - function to compute timestep 45415091d37SBarry Smith - ctx - [optional] user-defined context for private data 45582bf6240SBarry Smith required by the function (may be PETSC_NULL) 45682bf6240SBarry Smith 45715091d37SBarry Smith Level: intermediate 458fee21e36SBarry Smith 45982bf6240SBarry Smith Calling sequence of func: 46087828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 46182bf6240SBarry Smith 46282bf6240SBarry Smith . newdt - the newly computed timestep 46382bf6240SBarry Smith . ctx - [optional] timestep context 46482bf6240SBarry Smith 46582bf6240SBarry Smith Notes: 46682bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 46782bf6240SBarry Smith during the timestepping process. 46882bf6240SBarry Smith 46982bf6240SBarry Smith .keywords: timestep, pseudo, set 47082bf6240SBarry Smith 47182bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 47282bf6240SBarry Smith @*/ 47363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 47482bf6240SBarry Smith { 4756849ba73SBarry Smith PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *); 47682bf6240SBarry Smith 47782bf6240SBarry Smith PetscFunctionBegin; 4784482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 47982bf6240SBarry Smith 480c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 48182bf6240SBarry Smith if (f) { 48282bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 48382bf6240SBarry Smith } 48482bf6240SBarry Smith PetscFunctionReturn(0); 48582bf6240SBarry Smith } 48682bf6240SBarry Smith 48782bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 48882bf6240SBarry Smith 489a7cc72afSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 490fb2e594dSBarry Smith EXTERN_C_BEGIN 4914a2ae208SSatish Balay #undef __FUNCT__ 4924a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 49363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 49482bf6240SBarry Smith { 49582bf6240SBarry Smith TS_Pseudo *pseudo; 49682bf6240SBarry Smith 49782bf6240SBarry Smith PetscFunctionBegin; 49882bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 49982bf6240SBarry Smith pseudo->verify = dt; 50082bf6240SBarry Smith pseudo->verifyctx = ctx; 50182bf6240SBarry Smith PetscFunctionReturn(0); 50282bf6240SBarry Smith } 503fb2e594dSBarry Smith EXTERN_C_END 50482bf6240SBarry Smith 505fb2e594dSBarry Smith EXTERN_C_BEGIN 5064a2ae208SSatish Balay #undef __FUNCT__ 5074a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 50863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 50982bf6240SBarry Smith { 5104bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 51182bf6240SBarry Smith 51282bf6240SBarry Smith PetscFunctionBegin; 51382bf6240SBarry Smith pseudo->dt_increment = inc; 51482bf6240SBarry Smith PetscFunctionReturn(0); 51582bf6240SBarry Smith } 516fb2e594dSBarry Smith EXTERN_C_END 51782bf6240SBarry Smith 518fb2e594dSBarry Smith EXTERN_C_BEGIN 5194a2ae208SSatish Balay #undef __FUNCT__ 5204a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 52163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 52282bf6240SBarry Smith { 5234bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 52482bf6240SBarry Smith 52582bf6240SBarry Smith PetscFunctionBegin; 5264bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 52782bf6240SBarry Smith PetscFunctionReturn(0); 52882bf6240SBarry Smith } 529fb2e594dSBarry Smith EXTERN_C_END 53082bf6240SBarry Smith 5316849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 532fb2e594dSBarry Smith EXTERN_C_BEGIN 5334a2ae208SSatish Balay #undef __FUNCT__ 5344a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 53563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 53682bf6240SBarry Smith { 5374bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 53882bf6240SBarry Smith 53982bf6240SBarry Smith PetscFunctionBegin; 54082bf6240SBarry Smith pseudo->dt = dt; 54182bf6240SBarry Smith pseudo->dtctx = ctx; 54282bf6240SBarry Smith PetscFunctionReturn(0); 54382bf6240SBarry Smith } 544fb2e594dSBarry Smith EXTERN_C_END 54582bf6240SBarry Smith 54682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 54782bf6240SBarry Smith 548fb2e594dSBarry Smith EXTERN_C_BEGIN 5494a2ae208SSatish Balay #undef __FUNCT__ 5504a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 55163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Pseudo(TS ts) 5522d3f70b5SBarry Smith { 5537bf11e45SBarry Smith TS_Pseudo *pseudo; 554dfbe8321SBarry Smith PetscErrorCode ierr; 5552d3f70b5SBarry Smith 5563a40ed3dSBarry Smith PetscFunctionBegin; 557000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 558000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 5592d3f70b5SBarry Smith 5602d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 56129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 5622d3f70b5SBarry Smith } 5638beb423aSHong Zhang if (!ts->Arhs) { 56429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian"); 5652d3f70b5SBarry Smith } 566273d9f13SBarry Smith 567000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 568000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 569000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 5707bf11e45SBarry Smith 5717bf11e45SBarry Smith /* create the required nonlinear solver context */ 5727adad957SLisandro Dalcin ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 57317efb626SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 5742d3f70b5SBarry Smith 57538f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 5767bf11e45SBarry Smith ts->data = (void*)pseudo; 5772d3f70b5SBarry Smith 57828aa8177SBarry Smith pseudo->dt_increment = 1.1; 5794bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 58028aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 58182bf6240SBarry Smith 582f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 583e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 5840c97f09dSSatish Balay TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 585f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 586e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 5870c97f09dSSatish Balay TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 588f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 589e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 5900c97f09dSSatish Balay TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 5910c97f09dSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 5920c97f09dSSatish Balay "TSPseudoSetTimeStep_Pseudo", 5930c97f09dSSatish Balay TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 5943a40ed3dSBarry Smith PetscFunctionReturn(0); 5952d3f70b5SBarry Smith } 596fb2e594dSBarry Smith EXTERN_C_END 5972d3f70b5SBarry Smith 5984a2ae208SSatish Balay #undef __FUNCT__ 5994a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep" 60082bf6240SBarry Smith /*@C 60182bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 60282bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 60328aa8177SBarry Smith 60415091d37SBarry Smith Collective on TS 60515091d37SBarry Smith 60628aa8177SBarry Smith Input Parameters: 60728aa8177SBarry Smith . ts - the timestep context 60882bf6240SBarry Smith . dtctx - unused timestep context 60928aa8177SBarry Smith 61082bf6240SBarry Smith Output Parameter: 61182bf6240SBarry Smith . newdt - the timestep to use for the next step 61228aa8177SBarry Smith 61315091d37SBarry Smith Level: advanced 61415091d37SBarry Smith 61582bf6240SBarry Smith .keywords: timestep, pseudo, default 616564e8f4eSLois Curfman McInnes 61782bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 61828aa8177SBarry Smith @*/ 61963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 62028aa8177SBarry Smith { 62182bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 62287828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 623dfbe8321SBarry Smith PetscErrorCode ierr; 62428aa8177SBarry Smith 6253a40ed3dSBarry Smith PetscFunctionBegin; 62682bf6240SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 62782bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 62882bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 62982bf6240SBarry Smith /* first time through so compute initial function norm */ 63082bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 63182bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 63282bf6240SBarry Smith } 63382bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 63482bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 635004884caSBarry Smith } else if (pseudo->increment_dt_from_initial_dt) { 63682bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 63782bf6240SBarry Smith } else { 63882bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 63982bf6240SBarry Smith } 64082bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6413a40ed3dSBarry Smith PetscFunctionReturn(0); 64228aa8177SBarry Smith } 6432d3f70b5SBarry 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 657004884caSBarry Smith 658004884caSBarry Smith 659004884caSBarry Smith 660