1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*e1311b90SBarry Smith static char vcid[] = "$Id: posindep.c,v 1.25 1998/03/20 22:51:32 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" 173*e1311b90SBarry Smith static int TSDestroy_Pseudo(TS ts ) 1742d3f70b5SBarry Smith { 1757bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1762d3f70b5SBarry Smith int ierr; 1772d3f70b5SBarry Smith 1783a40ed3dSBarry Smith PetscFunctionBegin; 1797bf11e45SBarry Smith ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 1807bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1817bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 1822d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 1837bf11e45SBarry Smith PetscFree(pseudo); 1843a40ed3dSBarry Smith PetscFunctionReturn(0); 1852d3f70b5SBarry Smith } 1862d3f70b5SBarry Smith 1872d3f70b5SBarry Smith 1882d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1892d3f70b5SBarry Smith /* 1902d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 1912d3f70b5SBarry Smith */ 1922d3f70b5SBarry Smith 19358562591SBarry Smith #undef __FUNC__ 19458562591SBarry Smith #define __FUNC__ "TSPseudoMatMult" 1957bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 1962d3f70b5SBarry Smith { 1972d3f70b5SBarry Smith TS ts; 1982d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 1992d3f70b5SBarry Smith int ierr; 2002d3f70b5SBarry Smith 2013a40ed3dSBarry Smith PetscFunctionBegin; 2023a40ed3dSBarry Smith ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr); 2032d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2042d3f70b5SBarry Smith 2052d3f70b5SBarry Smith /* apply user provided function */ 2062d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 2072d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2082d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 2093a40ed3dSBarry Smith PetscFunctionReturn(0); 2102d3f70b5SBarry Smith } 2112d3f70b5SBarry Smith 2122d3f70b5SBarry Smith /* 2132d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2142d3f70b5SBarry Smith 2152d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2162d3f70b5SBarry Smith */ 21758562591SBarry Smith #undef __FUNC__ 21858562591SBarry Smith #define __FUNC__ "TSPseudoFunction" 2197bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2202d3f70b5SBarry Smith { 2212d3f70b5SBarry Smith TS ts = (TS) ctx; 2222d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2232d3f70b5SBarry Smith int ierr,i,n; 2242d3f70b5SBarry Smith 2253a40ed3dSBarry Smith PetscFunctionBegin; 2262d3f70b5SBarry Smith /* apply user provided function */ 2272d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 2287bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2292d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 2302d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 2312d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 2322d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 2332d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 2342d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2352d3f70b5SBarry Smith } 2362d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 2372d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 2382d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 2392d3f70b5SBarry Smith 2403a40ed3dSBarry Smith PetscFunctionReturn(0); 2412d3f70b5SBarry Smith } 2422d3f70b5SBarry Smith 2432d3f70b5SBarry Smith /* 2442d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2452d3f70b5SBarry Smith 2462d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2472d3f70b5SBarry Smith */ 24858562591SBarry Smith #undef __FUNC__ 24958562591SBarry Smith #define __FUNC__ "TSPseudoJacobian" 2507bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2512d3f70b5SBarry Smith { 2522d3f70b5SBarry Smith TS ts = (TS) ctx; 2532d3f70b5SBarry Smith int ierr; 2542d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 2552d3f70b5SBarry Smith MatType mtype; 2562d3f70b5SBarry Smith 2573a40ed3dSBarry Smith PetscFunctionBegin; 2582d3f70b5SBarry Smith /* construct users Jacobian */ 2592d3f70b5SBarry Smith if (ts->rhsjacobian) { 2602d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 2612d3f70b5SBarry Smith } 2622d3f70b5SBarry Smith 2632d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 2642d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 2652d3f70b5SBarry Smith if (mtype != MATSHELL) { 2662d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 2672d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 2682d3f70b5SBarry Smith } 2692d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 2702d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 2712d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 2722d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 2732d3f70b5SBarry Smith } 2742d3f70b5SBarry Smith 2753a40ed3dSBarry Smith PetscFunctionReturn(0); 2762d3f70b5SBarry Smith } 2772d3f70b5SBarry Smith 2782d3f70b5SBarry Smith 27958562591SBarry Smith #undef __FUNC__ 28058562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo" 2817bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 2822d3f70b5SBarry Smith { 2837bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2842d3f70b5SBarry Smith int ierr, M, m; 2852d3f70b5SBarry Smith 2863a40ed3dSBarry Smith PetscFunctionBegin; 2877bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 2887bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 2897bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 2902d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 2912d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 2922d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 2932d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 2947bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 2952d3f70b5SBarry Smith } 2967bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 2973a40ed3dSBarry Smith PetscFunctionReturn(0); 2982d3f70b5SBarry Smith } 2992d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3002d3f70b5SBarry Smith 30158562591SBarry Smith #undef __FUNC__ 302d4bb536fSBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor" 3032d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 3042d3f70b5SBarry Smith { 3057bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3062d3f70b5SBarry Smith 3073a40ed3dSBarry Smith PetscFunctionBegin; 30876be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 3093a40ed3dSBarry Smith PetscFunctionReturn(0); 3102d3f70b5SBarry Smith } 3112d3f70b5SBarry Smith 31258562591SBarry Smith #undef __FUNC__ 31358562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo" 3147bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3152d3f70b5SBarry Smith { 3162d3f70b5SBarry Smith int ierr,flg; 31728aa8177SBarry Smith double inc; 3182d3f70b5SBarry Smith 3193a40ed3dSBarry Smith PetscFunctionBegin; 3202d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 3212d3f70b5SBarry Smith 3222d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 3232d3f70b5SBarry Smith if (flg) { 32428aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 32528aa8177SBarry Smith } 32628aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 32728aa8177SBarry Smith if (flg) { 32828aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 3292d3f70b5SBarry Smith } 330ca4b7087SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr); 331ca4b7087SBarry Smith if (flg) { 332ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr); 333ca4b7087SBarry Smith } 3343a40ed3dSBarry Smith PetscFunctionReturn(0); 3352d3f70b5SBarry Smith } 3362d3f70b5SBarry Smith 33758562591SBarry Smith #undef __FUNC__ 338d4bb536fSBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo" 339ca4b7087SBarry Smith static int TSPrintHelp_Pseudo(TS ts,char *p) 3402d3f70b5SBarry Smith { 3413a40ed3dSBarry Smith PetscFunctionBegin; 34276be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n"); 34376be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p); 34476be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p); 34576be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," initial fnorm/current fnorm to determine new timestep\n"); 34676be9ce4SBarry Smith (*PetscHelpPrintf)(ts->comm," default is current_dt * previous fnorm/current fnorm\n"); 3473a40ed3dSBarry Smith PetscFunctionReturn(0); 3482d3f70b5SBarry Smith } 3492d3f70b5SBarry Smith 35058562591SBarry Smith #undef __FUNC__ 351d4bb536fSBarry Smith #define __FUNC__ "TSView_Pseudo" 352*e1311b90SBarry Smith static int TSView_Pseudo(TS ts,Viewer viewer) 3532d3f70b5SBarry Smith { 3543a40ed3dSBarry Smith PetscFunctionBegin; 3553a40ed3dSBarry Smith PetscFunctionReturn(0); 3562d3f70b5SBarry Smith } 3572d3f70b5SBarry Smith 35882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 35982bf6240SBarry Smith #undef __FUNC__ 36082bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep" 36182bf6240SBarry Smith /*@ 36282bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 36382bf6240SBarry Smith last timestep. 36482bf6240SBarry Smith 36582bf6240SBarry Smith Input Parameters: 36682bf6240SBarry Smith . ts - timestep context 36782bf6240SBarry Smith . dt - user-defined function to verify timestep 36882bf6240SBarry Smith . ctx - [optional] user-defined context for private data 36982bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 37082bf6240SBarry Smith 37182bf6240SBarry Smith Calling sequence of func: 37282bf6240SBarry Smith . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 37382bf6240SBarry Smith 37482bf6240SBarry Smith . update - latest solution vector 37582bf6240SBarry Smith . ctx - [optional] timestep context 37682bf6240SBarry Smith . newdt - the timestep to use for the next step 37782bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 37882bf6240SBarry Smith 37982bf6240SBarry Smith Notes: 38082bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 38182bf6240SBarry Smith during the timestepping process. 38282bf6240SBarry Smith 38382bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 38482bf6240SBarry Smith 38582bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 38682bf6240SBarry Smith @*/ 38782bf6240SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 38882bf6240SBarry Smith { 38982bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,Vec,void*,double *,int *),void *); 39082bf6240SBarry Smith 39182bf6240SBarry Smith PetscFunctionBegin; 39282bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 39382bf6240SBarry Smith 394*e1311b90SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep",(void **)&f);CHKERRQ(ierr); 39582bf6240SBarry Smith if (f) { 39682bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 39782bf6240SBarry Smith } 39882bf6240SBarry Smith PetscFunctionReturn(0); 39982bf6240SBarry Smith } 40082bf6240SBarry Smith 40182bf6240SBarry Smith #undef __FUNC__ 40282bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement" 40382bf6240SBarry Smith /*@ 40482bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 40582bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 40682bf6240SBarry Smith 40782bf6240SBarry Smith Input Parameters: 40882bf6240SBarry Smith . ts - the timestep context 40982bf6240SBarry Smith . inc - the scaling factor >= 1.0 41082bf6240SBarry Smith 41182bf6240SBarry Smith Options Database Key: 41282bf6240SBarry Smith $ -ts_pseudo_increment <increment> 41382bf6240SBarry Smith 41482bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 41582bf6240SBarry Smith 41682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 41782bf6240SBarry Smith @*/ 41882bf6240SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 41982bf6240SBarry Smith { 42082bf6240SBarry Smith int ierr, (*f)(TS,double); 42182bf6240SBarry Smith 42282bf6240SBarry Smith PetscFunctionBegin; 42382bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 42482bf6240SBarry Smith 425*e1311b90SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncremen",(void **)&f);CHKERRQ(ierr); 42682bf6240SBarry Smith if (f) { 42782bf6240SBarry Smith ierr = (*f)(ts,inc);CHKERRQ(ierr); 42882bf6240SBarry Smith } 42982bf6240SBarry Smith PetscFunctionReturn(0); 43082bf6240SBarry Smith } 43182bf6240SBarry Smith 43282bf6240SBarry Smith #undef __FUNC__ 43382bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt" 43482bf6240SBarry Smith /*@ 43582bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 43682bf6240SBarry Smith is computed via the formula 43782bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 43882bf6240SBarry Smith rather than the default update, 43982bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 44082bf6240SBarry Smith 44182bf6240SBarry Smith Input Parameter: 44282bf6240SBarry Smith . ts - the timestep context 44382bf6240SBarry Smith 44482bf6240SBarry Smith Options Database Key: 44582bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 44682bf6240SBarry Smith 44782bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 44882bf6240SBarry Smith 44982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 45082bf6240SBarry Smith @*/ 45182bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts) 45282bf6240SBarry Smith { 45382bf6240SBarry Smith int ierr, (*f)(TS); 45482bf6240SBarry Smith 45582bf6240SBarry Smith PetscFunctionBegin; 45682bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 45782bf6240SBarry Smith 458*e1311b90SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt",(void **)&f);CHKERRQ(ierr); 45982bf6240SBarry Smith if (f) { 46082bf6240SBarry Smith ierr = (*f)(ts);CHKERRQ(ierr); 46182bf6240SBarry Smith } 46282bf6240SBarry Smith PetscFunctionReturn(0); 46382bf6240SBarry Smith } 46482bf6240SBarry Smith 46582bf6240SBarry Smith 46682bf6240SBarry Smith #undef __FUNC__ 46782bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep" 46882bf6240SBarry Smith /*@ 46982bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 47082bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 47182bf6240SBarry Smith 47282bf6240SBarry Smith Input Parameters: 47382bf6240SBarry Smith . ts - timestep context 47482bf6240SBarry Smith . dt - function to compute timestep 47582bf6240SBarry Smith . ctx - [optional] user-defined context for private data 47682bf6240SBarry Smith required by the function (may be PETSC_NULL) 47782bf6240SBarry Smith 47882bf6240SBarry Smith Calling sequence of func: 47982bf6240SBarry Smith . func (TS ts,double *newdt,void *ctx); 48082bf6240SBarry Smith 48182bf6240SBarry Smith . newdt - the newly computed timestep 48282bf6240SBarry Smith . ctx - [optional] timestep context 48382bf6240SBarry Smith 48482bf6240SBarry Smith Notes: 48582bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 48682bf6240SBarry Smith during the timestepping process. 48782bf6240SBarry Smith 48882bf6240SBarry Smith .keywords: timestep, pseudo, set 48982bf6240SBarry Smith 49082bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 49182bf6240SBarry Smith @*/ 49282bf6240SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 49382bf6240SBarry Smith { 49482bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,double *,void *),void *); 49582bf6240SBarry Smith 49682bf6240SBarry Smith PetscFunctionBegin; 49782bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 49882bf6240SBarry Smith 499*e1311b90SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep",(void **)&f);CHKERRQ(ierr); 50082bf6240SBarry Smith if (f) { 50182bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 50282bf6240SBarry Smith } 50382bf6240SBarry Smith PetscFunctionReturn(0); 50482bf6240SBarry Smith } 50582bf6240SBarry Smith 50682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 50782bf6240SBarry Smith 50882bf6240SBarry Smith #undef __FUNC__ 50982bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep_Pseudo" 51082bf6240SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 51182bf6240SBarry Smith { 51282bf6240SBarry Smith TS_Pseudo *pseudo; 51382bf6240SBarry Smith 51482bf6240SBarry Smith PetscFunctionBegin; 51582bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 51682bf6240SBarry Smith pseudo->verify = dt; 51782bf6240SBarry Smith pseudo->verifyctx = ctx; 51882bf6240SBarry Smith PetscFunctionReturn(0); 51982bf6240SBarry Smith } 52082bf6240SBarry Smith 52182bf6240SBarry Smith #undef __FUNC__ 52282bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement_Pseudo" 52382bf6240SBarry Smith int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,double inc) 52482bf6240SBarry Smith { 52582bf6240SBarry Smith TS_Pseudo *pseudo; 52682bf6240SBarry Smith 52782bf6240SBarry Smith PetscFunctionBegin; 52882bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 52982bf6240SBarry Smith pseudo->dt_increment = inc; 53082bf6240SBarry Smith PetscFunctionReturn(0); 53182bf6240SBarry Smith } 53282bf6240SBarry Smith 53382bf6240SBarry Smith #undef __FUNC__ 53482bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 53582bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 53682bf6240SBarry Smith { 53782bf6240SBarry Smith TS_Pseudo *pseudo; 53882bf6240SBarry Smith 53982bf6240SBarry Smith PetscFunctionBegin; 54082bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 54182bf6240SBarry Smith pseudo->increment_dt_from_initial_dt = 1; 54282bf6240SBarry Smith PetscFunctionReturn(0); 54382bf6240SBarry Smith } 54482bf6240SBarry Smith 54582bf6240SBarry Smith #undef __FUNC__ 54682bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep_Pseudo" 54782bf6240SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,double*,void*),void* ctx) 54882bf6240SBarry Smith { 54982bf6240SBarry Smith TS_Pseudo *pseudo; 55082bf6240SBarry Smith 55182bf6240SBarry Smith PetscFunctionBegin; 55282bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 55382bf6240SBarry Smith pseudo->dt = dt; 55482bf6240SBarry Smith pseudo->dtctx = ctx; 55582bf6240SBarry Smith PetscFunctionReturn(0); 55682bf6240SBarry Smith } 55782bf6240SBarry Smith 55882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 55982bf6240SBarry Smith 56058562591SBarry Smith #undef __FUNC__ 56158562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo" 5627bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 5632d3f70b5SBarry Smith { 5647bf11e45SBarry Smith TS_Pseudo *pseudo; 5652d3f70b5SBarry Smith int ierr; 5662d3f70b5SBarry Smith MatType mtype; 5672d3f70b5SBarry Smith 5683a40ed3dSBarry Smith PetscFunctionBegin; 5697bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 5707bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 5717bf11e45SBarry Smith ts->view = TSView_Pseudo; 5722d3f70b5SBarry Smith 5732d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 574a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems"); 5752d3f70b5SBarry Smith } 5762d3f70b5SBarry Smith if (!ts->A) { 577a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian"); 5782d3f70b5SBarry Smith } 5792d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 5802d3f70b5SBarry Smith if (mtype == MATSHELL) { 5812d3f70b5SBarry Smith ts->Ashell = ts->A; 5822d3f70b5SBarry Smith } 5837bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 5847bf11e45SBarry Smith ts->step = TSStep_Pseudo; 5857bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 5867bf11e45SBarry Smith 5877bf11e45SBarry Smith /* create the required nonlinear solver context */ 5882d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 5892d3f70b5SBarry Smith 5907bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 591eed86810SBarry Smith PLogObjectMemory(ts,sizeof(TS_Pseudo)); 592eed86810SBarry Smith 5937bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 5947bf11e45SBarry Smith ts->data = (void *) pseudo; 5952d3f70b5SBarry Smith 59628aa8177SBarry Smith pseudo->dt_increment = 1.1; 597ca4b7087SBarry Smith pseudo->increment_dt_from_initial_dt = 0; 59828aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 59982bf6240SBarry Smith 600*e1311b90SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep", 601*e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 602*e1311b90SBarry Smith (void*)TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 603*e1311b90SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement", 604*e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 605*e1311b90SBarry Smith (void*)TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 606*e1311b90SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt", 607*e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 608*e1311b90SBarry Smith (void*)TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 609*e1311b90SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep","TSPseudoSetTimeStep_Pseudo", 610*e1311b90SBarry Smith (void*)TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 6122d3f70b5SBarry Smith } 6132d3f70b5SBarry Smith 61458562591SBarry Smith #undef __FUNC__ 61582bf6240SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep" 61682bf6240SBarry Smith /*@C 61782bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 61882bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 61928aa8177SBarry Smith 62028aa8177SBarry Smith Input Parameters: 62128aa8177SBarry Smith . ts - the timestep context 62282bf6240SBarry Smith . dtctx - unused timestep context 62328aa8177SBarry Smith 62482bf6240SBarry Smith Output Parameter: 62582bf6240SBarry Smith . newdt - the timestep to use for the next step 62628aa8177SBarry Smith 62782bf6240SBarry Smith .keywords: timestep, pseudo, default 628564e8f4eSLois Curfman McInnes 62982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 63028aa8177SBarry Smith @*/ 63182bf6240SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 63228aa8177SBarry Smith { 63382bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 63482bf6240SBarry Smith double inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 63582bf6240SBarry Smith int ierr; 63628aa8177SBarry Smith 6373a40ed3dSBarry Smith PetscFunctionBegin; 63882bf6240SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 63982bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 64082bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 64182bf6240SBarry Smith /* first time through so compute initial function norm */ 64282bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 64382bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 64482bf6240SBarry Smith } 64582bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 64682bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 64782bf6240SBarry Smith } 64882bf6240SBarry Smith else if (pseudo->increment_dt_from_initial_dt) { 64982bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 65082bf6240SBarry Smith } else { 65182bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 65282bf6240SBarry Smith } 65382bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6543a40ed3dSBarry Smith PetscFunctionReturn(0); 65528aa8177SBarry Smith } 6562d3f70b5SBarry Smith 657