1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*82bf6240SBarry Smith static char vcid[] = "$Id: posindep.c,v 1.23 1998/01/14 02:43:54 bsmith Exp bsmith $"; 32d3f70b5SBarry Smith #endif 42d3f70b5SBarry Smith /* 5fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 62d3f70b5SBarry Smith */ 72d3f70b5SBarry Smith #include <math.h> 82d3f70b5SBarry Smith #include "src/ts/tsimpl.h" /*I "ts.h" I*/ 92d3f70b5SBarry Smith #include "pinclude/pviewer.h" 102d3f70b5SBarry Smith 112d3f70b5SBarry Smith 122d3f70b5SBarry Smith typedef struct { 132d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 142d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 152d3f70b5SBarry Smith Vec rhs; /* work vector for RHS; vec_sol/dt */ 162d3f70b5SBarry Smith 172d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 182d3f70b5SBarry Smith 197bf11e45SBarry Smith int (*dt)(TS,double*,void*); /* compute next timestep, and related context */ 202d3f70b5SBarry Smith void *dtctx; 217bf11e45SBarry Smith int (*verify)(TS,Vec,void*,double*,int*); /* verify previous timestep and related context */ 227bf11e45SBarry Smith void *verifyctx; 232d3f70b5SBarry Smith 247bf11e45SBarry Smith double initial_fnorm,fnorm; /* original and current norm of F(u) */ 25ca4b7087SBarry Smith double fnorm_previous; 2628aa8177SBarry Smith 2728aa8177SBarry Smith double dt_increment; /* scaling that dt is incremented each time-step */ 28ca4b7087SBarry Smith int increment_dt_from_initial_dt; 297bf11e45SBarry Smith } TS_Pseudo; 302d3f70b5SBarry Smith 312d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 322d3f70b5SBarry Smith 3358562591SBarry Smith #undef __FUNC__ 3458562591SBarry Smith #define __FUNC__ "TSPseudoComputeTimeStep" 357bf11e45SBarry Smith /*@ 367bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 37564e8f4eSLois Curfman McInnes pseudo-timestepping process. 382d3f70b5SBarry Smith 397bf11e45SBarry Smith Input Parameter: 407bf11e45SBarry Smith . ts - timestep context 417bf11e45SBarry Smith 427bf11e45SBarry Smith Output Parameter: 43fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 44fb4a63b6SLois Curfman McInnes 45564e8f4eSLois Curfman McInnes 46564e8f4eSLois Curfman McInnes Notes: 47564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 48564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 49564e8f4eSLois Curfman McInnes 50fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 51564e8f4eSLois Curfman McInnes 52564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 537bf11e45SBarry Smith @*/ 547bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 557bf11e45SBarry Smith { 567bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 577bf11e45SBarry Smith int ierr; 587bf11e45SBarry Smith 593a40ed3dSBarry Smith PetscFunctionBegin; 607bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 617bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr); 627bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 633a40ed3dSBarry Smith PetscFunctionReturn(0); 647bf11e45SBarry Smith } 657bf11e45SBarry Smith 667bf11e45SBarry Smith 677bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 6858562591SBarry Smith #undef __FUNC__ 6958562591SBarry Smith #define __FUNC__ "TSPseudoDefaultVerifyTimeStep" 707bf11e45SBarry Smith /*@C 71639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 727bf11e45SBarry Smith 737bf11e45SBarry Smith Input Parameters: 747bf11e45SBarry Smith . ts - the timestep context 757bf11e45SBarry Smith . dtctx - unused timestep context 76564e8f4eSLois Curfman McInnes . update - latest solution vector 777bf11e45SBarry Smith 78564e8f4eSLois Curfman McInnes Output Parameters: 797bf11e45SBarry Smith . newdt - the timestep to use for the next step 80564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 817bf11e45SBarry Smith 82564e8f4eSLois Curfman McInnes Note: 83564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 84564e8f4eSLois Curfman McInnes timestep. 85564e8f4eSLois Curfman McInnes 86564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 87564e8f4eSLois Curfman McInnes 88564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 897bf11e45SBarry Smith @*/ 907bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag) 917bf11e45SBarry Smith { 923a40ed3dSBarry Smith PetscFunctionBegin; 937bf11e45SBarry Smith *flag = 1; 943a40ed3dSBarry Smith PetscFunctionReturn(0); 957bf11e45SBarry Smith } 967bf11e45SBarry Smith 977bf11e45SBarry Smith 9858562591SBarry Smith #undef __FUNC__ 9958562591SBarry Smith #define __FUNC__ "TSPseudoVerifyTimeStep" 1007bf11e45SBarry Smith /*@ 101564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1027bf11e45SBarry Smith 103fb4a63b6SLois Curfman McInnes Input Parameters: 1047bf11e45SBarry Smith . ts - timestep context 105564e8f4eSLois Curfman McInnes . update - latest solution vector 1067bf11e45SBarry Smith 107fb4a63b6SLois Curfman McInnes Output Parameters: 108fb4a63b6SLois Curfman McInnes . dt - newly computed timestep (if it had to shrink) 1097bf11e45SBarry Smith . flag - indicates if current timestep was ok 1107bf11e45SBarry Smith 111564e8f4eSLois Curfman McInnes Notes: 112564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 113564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 114564e8f4eSLois Curfman McInnes 115564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 116564e8f4eSLois Curfman McInnes 117564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1187bf11e45SBarry Smith @*/ 1197bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 1207bf11e45SBarry Smith { 1217bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1227bf11e45SBarry Smith int ierr; 1237bf11e45SBarry Smith 1243a40ed3dSBarry Smith PetscFunctionBegin; 1253a40ed3dSBarry Smith if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);} 1267bf11e45SBarry Smith 1277bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 1287bf11e45SBarry Smith 1293a40ed3dSBarry Smith PetscFunctionReturn(0); 1307bf11e45SBarry Smith } 1317bf11e45SBarry Smith 1327bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1337bf11e45SBarry Smith 13458562591SBarry Smith #undef __FUNC__ 13558562591SBarry Smith #define __FUNC__ "TSStep_Pseudo" 1367bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 1372d3f70b5SBarry Smith { 1382d3f70b5SBarry Smith Vec sol = ts->vec_sol; 139f2267985SLois Curfman McInnes int ierr,i,max_steps = ts->max_steps,its,ok,lits; 1407bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1417bf11e45SBarry Smith double current_time_step; 1422d3f70b5SBarry Smith 1433a40ed3dSBarry Smith PetscFunctionBegin; 1442d3f70b5SBarry Smith *steps = -ts->steps; 1452d3f70b5SBarry Smith 1467bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 1472d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 1487bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 1497bf11e45SBarry Smith current_time_step = ts->time_step; 1507bf11e45SBarry Smith while (1) { 1517bf11e45SBarry Smith ts->ptime += current_time_step; 1527bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 153f2267985SLois Curfman McInnes ierr = SNESGetNumberLinearIterations(ts->snes,&lits); CHKERRQ(ierr); 154f2267985SLois Curfman McInnes ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits; 1557bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 1567bf11e45SBarry Smith if (ok) break; 1577bf11e45SBarry Smith ts->ptime -= current_time_step; 1587bf11e45SBarry Smith current_time_step = ts->time_step; 1597bf11e45SBarry Smith } 1607bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 1612d3f70b5SBarry Smith ts->steps++; 1622d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1632d3f70b5SBarry Smith } 1642d3f70b5SBarry Smith 1652d3f70b5SBarry Smith *steps += ts->steps; 1662d3f70b5SBarry Smith *time = ts->ptime; 1673a40ed3dSBarry Smith PetscFunctionReturn(0); 1682d3f70b5SBarry Smith } 1692d3f70b5SBarry Smith 1702d3f70b5SBarry Smith /*------------------------------------------------------------*/ 17158562591SBarry Smith #undef __FUNC__ 172d4bb536fSBarry Smith #define __FUNC__ "TSDestroy_Pseudo" 1737bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj ) 1742d3f70b5SBarry Smith { 1752d3f70b5SBarry Smith TS ts = (TS) obj; 1767bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1772d3f70b5SBarry Smith int ierr; 1782d3f70b5SBarry Smith 1793a40ed3dSBarry Smith PetscFunctionBegin; 1807bf11e45SBarry Smith ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 1817bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1827bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 1832d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 1847bf11e45SBarry Smith PetscFree(pseudo); 1853a40ed3dSBarry Smith PetscFunctionReturn(0); 1862d3f70b5SBarry Smith } 1872d3f70b5SBarry Smith 1882d3f70b5SBarry Smith 1892d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1902d3f70b5SBarry Smith /* 1912d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 1922d3f70b5SBarry Smith */ 1932d3f70b5SBarry Smith 19458562591SBarry Smith #undef __FUNC__ 19558562591SBarry Smith #define __FUNC__ "TSPseudoMatMult" 1967bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 1972d3f70b5SBarry Smith { 1982d3f70b5SBarry Smith TS ts; 1992d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 2002d3f70b5SBarry Smith int ierr; 2012d3f70b5SBarry Smith 2023a40ed3dSBarry Smith PetscFunctionBegin; 2033a40ed3dSBarry Smith ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr); 2042d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2052d3f70b5SBarry Smith 2062d3f70b5SBarry Smith /* apply user provided function */ 2072d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 2082d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2092d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 2103a40ed3dSBarry Smith PetscFunctionReturn(0); 2112d3f70b5SBarry Smith } 2122d3f70b5SBarry Smith 2132d3f70b5SBarry Smith /* 2142d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2152d3f70b5SBarry Smith 2162d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2172d3f70b5SBarry Smith */ 21858562591SBarry Smith #undef __FUNC__ 21958562591SBarry Smith #define __FUNC__ "TSPseudoFunction" 2207bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2212d3f70b5SBarry Smith { 2222d3f70b5SBarry Smith TS ts = (TS) ctx; 2232d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2242d3f70b5SBarry Smith int ierr,i,n; 2252d3f70b5SBarry Smith 2263a40ed3dSBarry Smith PetscFunctionBegin; 2272d3f70b5SBarry Smith /* apply user provided function */ 2282d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 2297bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2302d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 2312d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 2322d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 2332d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 2342d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 2352d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2362d3f70b5SBarry Smith } 2372d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 2382d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 2392d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 2402d3f70b5SBarry Smith 2413a40ed3dSBarry Smith PetscFunctionReturn(0); 2422d3f70b5SBarry Smith } 2432d3f70b5SBarry Smith 2442d3f70b5SBarry Smith /* 2452d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2462d3f70b5SBarry Smith 2472d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2482d3f70b5SBarry Smith */ 24958562591SBarry Smith #undef __FUNC__ 25058562591SBarry Smith #define __FUNC__ "TSPseudoJacobian" 2517bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2522d3f70b5SBarry Smith { 2532d3f70b5SBarry Smith TS ts = (TS) ctx; 2542d3f70b5SBarry Smith int ierr; 2552d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 2562d3f70b5SBarry Smith MatType mtype; 2572d3f70b5SBarry Smith 2583a40ed3dSBarry Smith PetscFunctionBegin; 2592d3f70b5SBarry Smith /* construct users Jacobian */ 2602d3f70b5SBarry Smith if (ts->rhsjacobian) { 2612d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 2622d3f70b5SBarry Smith } 2632d3f70b5SBarry Smith 2642d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 2652d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 2662d3f70b5SBarry Smith if (mtype != MATSHELL) { 2672d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 2682d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 2692d3f70b5SBarry Smith } 2702d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 2712d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 2722d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 2732d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 2742d3f70b5SBarry Smith } 2752d3f70b5SBarry Smith 2763a40ed3dSBarry Smith PetscFunctionReturn(0); 2772d3f70b5SBarry Smith } 2782d3f70b5SBarry Smith 2792d3f70b5SBarry Smith 28058562591SBarry Smith #undef __FUNC__ 28158562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo" 2827bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 2832d3f70b5SBarry Smith { 2847bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2852d3f70b5SBarry Smith int ierr, M, m; 2862d3f70b5SBarry Smith 2873a40ed3dSBarry Smith PetscFunctionBegin; 2887bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 2897bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 2907bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 2912d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 2922d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 2932d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 2942d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 2957bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 2962d3f70b5SBarry Smith } 2977bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 2983a40ed3dSBarry Smith PetscFunctionReturn(0); 2992d3f70b5SBarry Smith } 3002d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3012d3f70b5SBarry Smith 30258562591SBarry Smith #undef __FUNC__ 303d4bb536fSBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor" 3042d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 3052d3f70b5SBarry Smith { 3067bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3072d3f70b5SBarry Smith 3083a40ed3dSBarry Smith PetscFunctionBegin; 30976be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 3103a40ed3dSBarry Smith PetscFunctionReturn(0); 3112d3f70b5SBarry Smith } 3122d3f70b5SBarry Smith 31358562591SBarry Smith #undef __FUNC__ 31458562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo" 3157bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3162d3f70b5SBarry Smith { 3172d3f70b5SBarry Smith int ierr,flg; 31828aa8177SBarry Smith double inc; 3192d3f70b5SBarry Smith 3203a40ed3dSBarry Smith PetscFunctionBegin; 3212d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 3222d3f70b5SBarry Smith 3232d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 3242d3f70b5SBarry Smith if (flg) { 32528aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 32628aa8177SBarry Smith } 32728aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 32828aa8177SBarry Smith if (flg) { 32928aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 3302d3f70b5SBarry Smith } 331ca4b7087SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr); 332ca4b7087SBarry Smith if (flg) { 333ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr); 334ca4b7087SBarry Smith } 3353a40ed3dSBarry Smith PetscFunctionReturn(0); 3362d3f70b5SBarry Smith } 3372d3f70b5SBarry Smith 33858562591SBarry Smith #undef __FUNC__ 339d4bb536fSBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo" 340ca4b7087SBarry Smith static int TSPrintHelp_Pseudo(TS ts,char *p) 3412d3f70b5SBarry Smith { 3423a40ed3dSBarry Smith PetscFunctionBegin; 34376be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n"); 34476be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p); 34576be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p); 34676be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," initial fnorm/current fnorm to determine new timestep\n"); 34776be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," default is current_dt * previous fnorm/current fnorm\n"); 3483a40ed3dSBarry Smith PetscFunctionReturn(0); 3492d3f70b5SBarry Smith } 3502d3f70b5SBarry Smith 35158562591SBarry Smith #undef __FUNC__ 352d4bb536fSBarry Smith #define __FUNC__ "TSView_Pseudo" 3537bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer) 3542d3f70b5SBarry Smith { 3553a40ed3dSBarry Smith PetscFunctionBegin; 3563a40ed3dSBarry Smith PetscFunctionReturn(0); 3572d3f70b5SBarry Smith } 3582d3f70b5SBarry Smith 359*82bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 360*82bf6240SBarry Smith #undef __FUNC__ 361*82bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep" 362*82bf6240SBarry Smith /*@ 363*82bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 364*82bf6240SBarry Smith last timestep. 365*82bf6240SBarry Smith 366*82bf6240SBarry Smith Input Parameters: 367*82bf6240SBarry Smith . ts - timestep context 368*82bf6240SBarry Smith . dt - user-defined function to verify timestep 369*82bf6240SBarry Smith . ctx - [optional] user-defined context for private data 370*82bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 371*82bf6240SBarry Smith 372*82bf6240SBarry Smith Calling sequence of func: 373*82bf6240SBarry Smith . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 374*82bf6240SBarry Smith 375*82bf6240SBarry Smith . update - latest solution vector 376*82bf6240SBarry Smith . ctx - [optional] timestep context 377*82bf6240SBarry Smith . newdt - the timestep to use for the next step 378*82bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 379*82bf6240SBarry Smith 380*82bf6240SBarry Smith Notes: 381*82bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 382*82bf6240SBarry Smith during the timestepping process. 383*82bf6240SBarry Smith 384*82bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 385*82bf6240SBarry Smith 386*82bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 387*82bf6240SBarry Smith @*/ 388*82bf6240SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 389*82bf6240SBarry Smith { 390*82bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,Vec,void*,double *,int *),void *); 391*82bf6240SBarry Smith 392*82bf6240SBarry Smith PetscFunctionBegin; 393*82bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 394*82bf6240SBarry Smith 395*82bf6240SBarry Smith ierr = DLRegisterFind(ts->qlist,"TSPseudoSetVerifyTimeStep",(int (**)(void *))&f);CHKERRQ(ierr); 396*82bf6240SBarry Smith if (f) { 397*82bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 398*82bf6240SBarry Smith } 399*82bf6240SBarry Smith PetscFunctionReturn(0); 400*82bf6240SBarry Smith } 401*82bf6240SBarry Smith 402*82bf6240SBarry Smith #undef __FUNC__ 403*82bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement" 404*82bf6240SBarry Smith /*@ 405*82bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 406*82bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 407*82bf6240SBarry Smith 408*82bf6240SBarry Smith Input Parameters: 409*82bf6240SBarry Smith . ts - the timestep context 410*82bf6240SBarry Smith . inc - the scaling factor >= 1.0 411*82bf6240SBarry Smith 412*82bf6240SBarry Smith Options Database Key: 413*82bf6240SBarry Smith $ -ts_pseudo_increment <increment> 414*82bf6240SBarry Smith 415*82bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 416*82bf6240SBarry Smith 417*82bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 418*82bf6240SBarry Smith @*/ 419*82bf6240SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 420*82bf6240SBarry Smith { 421*82bf6240SBarry Smith int ierr, (*f)(TS,double); 422*82bf6240SBarry Smith 423*82bf6240SBarry Smith PetscFunctionBegin; 424*82bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 425*82bf6240SBarry Smith 426*82bf6240SBarry Smith ierr = DLRegisterFind(ts->qlist,"TSPseudoSetTimeStepIncrement",(int (**)(void *))&f);CHKERRQ(ierr); 427*82bf6240SBarry Smith if (f) { 428*82bf6240SBarry Smith ierr = (*f)(ts,inc);CHKERRQ(ierr); 429*82bf6240SBarry Smith } 430*82bf6240SBarry Smith PetscFunctionReturn(0); 431*82bf6240SBarry Smith } 432*82bf6240SBarry Smith 433*82bf6240SBarry Smith #undef __FUNC__ 434*82bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt" 435*82bf6240SBarry Smith /*@ 436*82bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 437*82bf6240SBarry Smith is computed via the formula 438*82bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 439*82bf6240SBarry Smith rather than the default update, 440*82bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 441*82bf6240SBarry Smith 442*82bf6240SBarry Smith Input Parameter: 443*82bf6240SBarry Smith . ts - the timestep context 444*82bf6240SBarry Smith 445*82bf6240SBarry Smith Options Database Key: 446*82bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 447*82bf6240SBarry Smith 448*82bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 449*82bf6240SBarry Smith 450*82bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 451*82bf6240SBarry Smith @*/ 452*82bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts) 453*82bf6240SBarry Smith { 454*82bf6240SBarry Smith int ierr, (*f)(TS); 455*82bf6240SBarry Smith 456*82bf6240SBarry Smith PetscFunctionBegin; 457*82bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 458*82bf6240SBarry Smith 459*82bf6240SBarry Smith ierr = DLRegisterFind(ts->qlist,"TSPseudoIncrementDtFromInitialDt",(int (**)(void *))&f);CHKERRQ(ierr); 460*82bf6240SBarry Smith if (f) { 461*82bf6240SBarry Smith ierr = (*f)(ts);CHKERRQ(ierr); 462*82bf6240SBarry Smith } 463*82bf6240SBarry Smith PetscFunctionReturn(0); 464*82bf6240SBarry Smith } 465*82bf6240SBarry Smith 466*82bf6240SBarry Smith 467*82bf6240SBarry Smith #undef __FUNC__ 468*82bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep" 469*82bf6240SBarry Smith /*@ 470*82bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 471*82bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 472*82bf6240SBarry Smith 473*82bf6240SBarry Smith Input Parameters: 474*82bf6240SBarry Smith . ts - timestep context 475*82bf6240SBarry Smith . dt - function to compute timestep 476*82bf6240SBarry Smith . ctx - [optional] user-defined context for private data 477*82bf6240SBarry Smith required by the function (may be PETSC_NULL) 478*82bf6240SBarry Smith 479*82bf6240SBarry Smith Calling sequence of func: 480*82bf6240SBarry Smith . func (TS ts,double *newdt,void *ctx); 481*82bf6240SBarry Smith 482*82bf6240SBarry Smith . newdt - the newly computed timestep 483*82bf6240SBarry Smith . ctx - [optional] timestep context 484*82bf6240SBarry Smith 485*82bf6240SBarry Smith Notes: 486*82bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 487*82bf6240SBarry Smith during the timestepping process. 488*82bf6240SBarry Smith 489*82bf6240SBarry Smith .keywords: timestep, pseudo, set 490*82bf6240SBarry Smith 491*82bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 492*82bf6240SBarry Smith @*/ 493*82bf6240SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 494*82bf6240SBarry Smith { 495*82bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,double *,void *),void *); 496*82bf6240SBarry Smith 497*82bf6240SBarry Smith PetscFunctionBegin; 498*82bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 499*82bf6240SBarry Smith 500*82bf6240SBarry Smith ierr = DLRegisterFind(ts->qlist,"TSPseudoSetTimeStep",(int (**)(void *))&f);CHKERRQ(ierr); 501*82bf6240SBarry Smith if (f) { 502*82bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 503*82bf6240SBarry Smith } 504*82bf6240SBarry Smith PetscFunctionReturn(0); 505*82bf6240SBarry Smith } 506*82bf6240SBarry Smith 507*82bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 508*82bf6240SBarry Smith 509*82bf6240SBarry Smith #undef __FUNC__ 510*82bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep_Pseudo" 511*82bf6240SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 512*82bf6240SBarry Smith { 513*82bf6240SBarry Smith TS_Pseudo *pseudo; 514*82bf6240SBarry Smith 515*82bf6240SBarry Smith PetscFunctionBegin; 516*82bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 517*82bf6240SBarry Smith pseudo->verify = dt; 518*82bf6240SBarry Smith pseudo->verifyctx = ctx; 519*82bf6240SBarry Smith PetscFunctionReturn(0); 520*82bf6240SBarry Smith } 521*82bf6240SBarry Smith 522*82bf6240SBarry Smith #undef __FUNC__ 523*82bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement_Pseudo" 524*82bf6240SBarry Smith int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,double inc) 525*82bf6240SBarry Smith { 526*82bf6240SBarry Smith TS_Pseudo *pseudo; 527*82bf6240SBarry Smith 528*82bf6240SBarry Smith PetscFunctionBegin; 529*82bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 530*82bf6240SBarry Smith pseudo->dt_increment = inc; 531*82bf6240SBarry Smith PetscFunctionReturn(0); 532*82bf6240SBarry Smith } 533*82bf6240SBarry Smith 534*82bf6240SBarry Smith #undef __FUNC__ 535*82bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 536*82bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 537*82bf6240SBarry Smith { 538*82bf6240SBarry Smith TS_Pseudo *pseudo; 539*82bf6240SBarry Smith 540*82bf6240SBarry Smith PetscFunctionBegin; 541*82bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 542*82bf6240SBarry Smith pseudo->increment_dt_from_initial_dt = 1; 543*82bf6240SBarry Smith PetscFunctionReturn(0); 544*82bf6240SBarry Smith } 545*82bf6240SBarry Smith 546*82bf6240SBarry Smith #undef __FUNC__ 547*82bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep_Pseudo" 548*82bf6240SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,double*,void*),void* ctx) 549*82bf6240SBarry Smith { 550*82bf6240SBarry Smith TS_Pseudo *pseudo; 551*82bf6240SBarry Smith 552*82bf6240SBarry Smith PetscFunctionBegin; 553*82bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 554*82bf6240SBarry Smith pseudo->dt = dt; 555*82bf6240SBarry Smith pseudo->dtctx = ctx; 556*82bf6240SBarry Smith PetscFunctionReturn(0); 557*82bf6240SBarry Smith } 558*82bf6240SBarry Smith 559*82bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 560*82bf6240SBarry Smith 56158562591SBarry Smith #undef __FUNC__ 56258562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo" 5637bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 5642d3f70b5SBarry Smith { 5657bf11e45SBarry Smith TS_Pseudo *pseudo; 5662d3f70b5SBarry Smith int ierr; 5672d3f70b5SBarry Smith MatType mtype; 5682d3f70b5SBarry Smith 5693a40ed3dSBarry Smith PetscFunctionBegin; 5707bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 5717bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 5727bf11e45SBarry Smith ts->view = TSView_Pseudo; 5732d3f70b5SBarry Smith 5742d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 575a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems"); 5762d3f70b5SBarry Smith } 5772d3f70b5SBarry Smith if (!ts->A) { 578a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian"); 5792d3f70b5SBarry Smith } 5802d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 5812d3f70b5SBarry Smith if (mtype == MATSHELL) { 5822d3f70b5SBarry Smith ts->Ashell = ts->A; 5832d3f70b5SBarry Smith } 5847bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 5857bf11e45SBarry Smith ts->step = TSStep_Pseudo; 5867bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 5877bf11e45SBarry Smith 5887bf11e45SBarry Smith /* create the required nonlinear solver context */ 5892d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 5902d3f70b5SBarry Smith 5917bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 592eed86810SBarry Smith PLogObjectMemory(ts,sizeof(TS_Pseudo)); 593eed86810SBarry Smith 5947bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 5957bf11e45SBarry Smith ts->data = (void *) pseudo; 5962d3f70b5SBarry Smith 59728aa8177SBarry Smith pseudo->dt_increment = 1.1; 598ca4b7087SBarry Smith pseudo->increment_dt_from_initial_dt = 0; 59928aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 600*82bf6240SBarry Smith 601*82bf6240SBarry Smith ierr = DLRegister(&ts->qlist,"TSPseudoSetVerifyTimeStep","TSPseudoSetVerifyTimeStep_Pseudo", 602*82bf6240SBarry Smith TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 603*82bf6240SBarry Smith ierr = DLRegister(&ts->qlist,"TSPseudoSetTimeStepIncrement","TSPseudoSetTimeStepIncrement_Pseudo", 604*82bf6240SBarry Smith TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 605*82bf6240SBarry Smith ierr = DLRegister(&ts->qlist,"TSPseudoIncrementDtFromInitialDt","TSPseudoIncrementDtFromInitialDt_Pseudo", 606*82bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 607*82bf6240SBarry Smith ierr = DLRegister(&ts->qlist,"TSPseudoSetTimeStep","TSPseudoSetTimeStep_Pseudo",TSPseudoSetTimeStep_Pseudo); 608*82bf6240SBarry Smith CHKERRQ(ierr); 6093a40ed3dSBarry Smith PetscFunctionReturn(0); 6102d3f70b5SBarry Smith } 6112d3f70b5SBarry Smith 61258562591SBarry Smith #undef __FUNC__ 613*82bf6240SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep" 614*82bf6240SBarry Smith /*@C 615*82bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 616*82bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 61728aa8177SBarry Smith 61828aa8177SBarry Smith Input Parameters: 61928aa8177SBarry Smith . ts - the timestep context 620*82bf6240SBarry Smith . dtctx - unused timestep context 62128aa8177SBarry Smith 622*82bf6240SBarry Smith Output Parameter: 623*82bf6240SBarry Smith . newdt - the timestep to use for the next step 62428aa8177SBarry Smith 625*82bf6240SBarry Smith .keywords: timestep, pseudo, default 626564e8f4eSLois Curfman McInnes 627*82bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 62828aa8177SBarry Smith @*/ 629*82bf6240SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 63028aa8177SBarry Smith { 631*82bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 632*82bf6240SBarry Smith double inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 633*82bf6240SBarry Smith int ierr; 63428aa8177SBarry Smith 6353a40ed3dSBarry Smith PetscFunctionBegin; 636*82bf6240SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 637*82bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 638*82bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 639*82bf6240SBarry Smith /* first time through so compute initial function norm */ 640*82bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 641*82bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 642*82bf6240SBarry Smith } 643*82bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 644*82bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 645*82bf6240SBarry Smith } 646*82bf6240SBarry Smith else if (pseudo->increment_dt_from_initial_dt) { 647*82bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 648*82bf6240SBarry Smith } else { 649*82bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 650*82bf6240SBarry Smith } 651*82bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6523a40ed3dSBarry Smith PetscFunctionReturn(0); 65328aa8177SBarry Smith } 6542d3f70b5SBarry Smith 655