12d3f70b5SBarry Smith #ifndef lint 2*28aa8177SBarry Smith static char vcid[] = "$Id: posindep.c,v 1.3 1996/09/30 20:23:41 curfman 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) */ 25*28aa8177SBarry Smith 26*28aa8177SBarry Smith double dt_increment; /* scaling that dt is incremented each time-step */ 277bf11e45SBarry Smith } TS_Pseudo; 282d3f70b5SBarry Smith 292d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 302d3f70b5SBarry Smith /*@C 31fb4a63b6SLois Curfman McInnes TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 32fb4a63b6SLois Curfman McInnes Use with TSPseudoSetTimeStep(). 332d3f70b5SBarry Smith 342d3f70b5SBarry Smith Input Parameters: 352d3f70b5SBarry Smith . ts - the timestep context 367bf11e45SBarry Smith . dtctx - unused timestep context 372d3f70b5SBarry Smith 382d3f70b5SBarry Smith Output Parameter: 392d3f70b5SBarry Smith . newdt - the timestep to use for the next step 402d3f70b5SBarry Smith 41fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, default 42fb4a63b6SLois Curfman McInnes 43fb4a63b6SLois Curfman McInnes .seealso: TSPseudoSetTimeStep() 442d3f70b5SBarry Smith @*/ 457bf11e45SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 462d3f70b5SBarry Smith { 477bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 48*28aa8177SBarry Smith double inc = pseudo->dt_increment; 492d3f70b5SBarry Smith int ierr; 502d3f70b5SBarry Smith 517bf11e45SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 527bf11e45SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 537bf11e45SBarry Smith if (pseudo->initial_fnorm == 0.0) { 542d3f70b5SBarry Smith /* first time through so compute initial function norm */ 557bf11e45SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 562d3f70b5SBarry Smith } 57*28aa8177SBarry Smith if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 58*28aa8177SBarry Smith else *newdt = inc*ts->time_step*pseudo->initial_fnorm/pseudo->fnorm; 592d3f70b5SBarry Smith return 0; 602d3f70b5SBarry Smith } 612d3f70b5SBarry Smith 622d3f70b5SBarry Smith /*@ 63fb4a63b6SLois Curfman McInnes TSPseudoSetTimeStep - Sets the user-defined routine to be 64fb4a63b6SLois Curfman McInnes called at each pseudo-timestep to update the timestep. 652d3f70b5SBarry Smith 662d3f70b5SBarry Smith Input Parameters: 672d3f70b5SBarry Smith . ts - timestep context 682d3f70b5SBarry Smith . dt - function to compute timestep 692d3f70b5SBarry Smith . ctx - [optional] context required by function 702d3f70b5SBarry Smith 71fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, set 72fb4a63b6SLois Curfman McInnes 73fb4a63b6SLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep() 742d3f70b5SBarry Smith @*/ 757bf11e45SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 762d3f70b5SBarry Smith { 777bf11e45SBarry Smith TS_Pseudo *pseudo; 782d3f70b5SBarry Smith 792d3f70b5SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 807bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 817bf11e45SBarry Smith 827bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 837bf11e45SBarry Smith pseudo->dt = dt; 847bf11e45SBarry Smith pseudo->dtctx = ctx; 852d3f70b5SBarry Smith return 0; 862d3f70b5SBarry Smith } 872d3f70b5SBarry Smith 887bf11e45SBarry Smith /*@ 897bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 907bf11e45SBarry Smith pseudo-timestepping. 912d3f70b5SBarry Smith 927bf11e45SBarry Smith Input Parameter: 937bf11e45SBarry Smith . ts - timestep context 947bf11e45SBarry Smith 957bf11e45SBarry Smith Output Parameter: 96fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 97fb4a63b6SLois Curfman McInnes 98fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 997bf11e45SBarry Smith @*/ 1007bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 1017bf11e45SBarry Smith { 1027bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1037bf11e45SBarry Smith int ierr; 1047bf11e45SBarry Smith 1057bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 1067bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt, pseudo->dtctx); CHKERRQ(ierr); 1077bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 1087bf11e45SBarry Smith return 0; 1097bf11e45SBarry Smith } 1107bf11e45SBarry Smith 1117bf11e45SBarry Smith 1127bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 1137bf11e45SBarry Smith /*@C 1147bf11e45SBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify last timestep. 1157bf11e45SBarry Smith 1167bf11e45SBarry Smith Input Parameters: 1177bf11e45SBarry Smith . ts - the timestep context 1187bf11e45SBarry Smith . dtctx - unused timestep context 1197bf11e45SBarry Smith 1207bf11e45SBarry Smith Output Parameter: 1217bf11e45SBarry Smith . newdt - the timestep to use for the next step 1227bf11e45SBarry Smith 1237bf11e45SBarry Smith @*/ 1247bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void* dtctx,double* newdt,int *flag) 1257bf11e45SBarry Smith { 1267bf11e45SBarry Smith *flag = 1; 1277bf11e45SBarry Smith return 0; 1287bf11e45SBarry Smith } 1297bf11e45SBarry Smith 1307bf11e45SBarry Smith /*@ 1317bf11e45SBarry Smith TSPseudoSetVerifyTimeStep - Sets the user routine to verify quality of last timestep. 1327bf11e45SBarry Smith 1337bf11e45SBarry Smith Input Parameters: 1347bf11e45SBarry Smith . ts - timestep context 1357bf11e45SBarry Smith . dt - function to verify 1367bf11e45SBarry Smith . ctx - [optional] context required by function 1377bf11e45SBarry Smith 1387bf11e45SBarry Smith @*/ 1397bf11e45SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 1407bf11e45SBarry Smith { 1417bf11e45SBarry Smith TS_Pseudo *pseudo; 1427bf11e45SBarry Smith 1437bf11e45SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 1447bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 1457bf11e45SBarry Smith 1467bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 1477bf11e45SBarry Smith pseudo->verify = dt; 1487bf11e45SBarry Smith pseudo->verifyctx = ctx; 1497bf11e45SBarry Smith return 0; 1507bf11e45SBarry Smith } 1517bf11e45SBarry Smith 1527bf11e45SBarry Smith /*@ 153fb4a63b6SLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies that the last timestep was OK. 1547bf11e45SBarry Smith 155fb4a63b6SLois Curfman McInnes Input Parameters: 1567bf11e45SBarry Smith . ts - timestep context 1577bf11e45SBarry Smith . update - latest solution 1587bf11e45SBarry Smith 159fb4a63b6SLois Curfman McInnes Output Parameters: 160fb4a63b6SLois Curfman McInnes . dt - newly computed timestep (if it had to shrink) 1617bf11e45SBarry Smith . flag - indicates if current timestep was ok 1627bf11e45SBarry Smith 1637bf11e45SBarry Smith @*/ 1647bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 1657bf11e45SBarry Smith { 1667bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1677bf11e45SBarry Smith int ierr; 1687bf11e45SBarry Smith 1697bf11e45SBarry Smith if (!pseudo->verify) {*flag = 1; return 0;} 1707bf11e45SBarry Smith 1717bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 1727bf11e45SBarry Smith 1737bf11e45SBarry Smith return 0; 1747bf11e45SBarry Smith } 1757bf11e45SBarry Smith 1767bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1777bf11e45SBarry Smith 1787bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 1792d3f70b5SBarry Smith { 1802d3f70b5SBarry Smith Vec sol = ts->vec_sol; 1817bf11e45SBarry Smith int ierr,i,max_steps = ts->max_steps,its,ok; 1827bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1837bf11e45SBarry Smith double current_time_step; 1842d3f70b5SBarry Smith 1852d3f70b5SBarry Smith *steps = -ts->steps; 1862d3f70b5SBarry Smith 1877bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 1882d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 1897bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 1907bf11e45SBarry Smith current_time_step = ts->time_step; 1917bf11e45SBarry Smith while (1) { 1927bf11e45SBarry Smith ts->ptime += current_time_step; 1937bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 1947bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 1957bf11e45SBarry Smith if (ok) break; 1967bf11e45SBarry Smith ts->ptime -= current_time_step; 1977bf11e45SBarry Smith current_time_step = ts->time_step; 1987bf11e45SBarry Smith } 1997bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 2002d3f70b5SBarry Smith ts->steps++; 2012d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 2022d3f70b5SBarry Smith } 2032d3f70b5SBarry Smith 2042d3f70b5SBarry Smith *steps += ts->steps; 2052d3f70b5SBarry Smith *time = ts->ptime; 2062d3f70b5SBarry Smith return 0; 2072d3f70b5SBarry Smith } 2082d3f70b5SBarry Smith 2092d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2107bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj ) 2112d3f70b5SBarry Smith { 2122d3f70b5SBarry Smith TS ts = (TS) obj; 2137bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2142d3f70b5SBarry Smith int ierr; 2152d3f70b5SBarry Smith 2167bf11e45SBarry Smith ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 2177bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 2187bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 2192d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 2207bf11e45SBarry Smith PetscFree(pseudo); 2212d3f70b5SBarry Smith return 0; 2222d3f70b5SBarry Smith } 2232d3f70b5SBarry Smith 2242d3f70b5SBarry Smith 2252d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2262d3f70b5SBarry Smith /* 2272d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 2282d3f70b5SBarry Smith */ 2292d3f70b5SBarry Smith 2307bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 2312d3f70b5SBarry Smith { 2322d3f70b5SBarry Smith TS ts; 2332d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 2342d3f70b5SBarry Smith int ierr; 2352d3f70b5SBarry Smith 2362d3f70b5SBarry Smith MatShellGetContext(mat,(void **)&ts); 2372d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2382d3f70b5SBarry Smith 2392d3f70b5SBarry Smith /* apply user provided function */ 2402d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 2412d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2422d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 2432d3f70b5SBarry Smith return 0; 2442d3f70b5SBarry Smith } 2452d3f70b5SBarry Smith 2462d3f70b5SBarry Smith /* 2472d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2482d3f70b5SBarry Smith 2492d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2502d3f70b5SBarry Smith */ 2517bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2522d3f70b5SBarry Smith { 2532d3f70b5SBarry Smith TS ts = (TS) ctx; 2542d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2552d3f70b5SBarry Smith int ierr,i,n; 2562d3f70b5SBarry Smith 2572d3f70b5SBarry Smith /* apply user provided function */ 2582d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 2597bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2602d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 2612d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 2622d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 2632d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 2642d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 2652d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2662d3f70b5SBarry Smith } 2672d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 2682d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 2692d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 2702d3f70b5SBarry Smith 2712d3f70b5SBarry Smith return 0; 2722d3f70b5SBarry Smith } 2732d3f70b5SBarry Smith 2742d3f70b5SBarry Smith /* 2752d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2762d3f70b5SBarry Smith 2772d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2782d3f70b5SBarry Smith */ 2797bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2802d3f70b5SBarry Smith { 2812d3f70b5SBarry Smith TS ts = (TS) ctx; 2822d3f70b5SBarry Smith int ierr; 2832d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 2842d3f70b5SBarry Smith MatType mtype; 2852d3f70b5SBarry Smith 2862d3f70b5SBarry Smith /* construct users Jacobian */ 2872d3f70b5SBarry Smith if (ts->rhsjacobian) { 2882d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 2892d3f70b5SBarry Smith } 2902d3f70b5SBarry Smith 2912d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 2922d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 2932d3f70b5SBarry Smith if (mtype != MATSHELL) { 2942d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 2952d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 2962d3f70b5SBarry Smith } 2972d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 2982d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 2992d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 3002d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 3012d3f70b5SBarry Smith } 3022d3f70b5SBarry Smith 3032d3f70b5SBarry Smith return 0; 3042d3f70b5SBarry Smith } 3052d3f70b5SBarry Smith 3062d3f70b5SBarry Smith 3077bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 3082d3f70b5SBarry Smith { 3097bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3102d3f70b5SBarry Smith int ierr, M, m; 3112d3f70b5SBarry Smith 3127bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 3137bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 3147bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 3152d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 3162d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 3172d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 3182d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 3197bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 3202d3f70b5SBarry Smith } 3217bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 3222d3f70b5SBarry Smith return 0; 3232d3f70b5SBarry Smith } 3242d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3252d3f70b5SBarry Smith 3262d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 3272d3f70b5SBarry Smith { 3287bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3292d3f70b5SBarry Smith 3307bf11e45SBarry Smith PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 3312d3f70b5SBarry Smith return 0; 3322d3f70b5SBarry Smith } 3332d3f70b5SBarry Smith 3347bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3352d3f70b5SBarry Smith { 3362d3f70b5SBarry Smith int ierr,flg; 337*28aa8177SBarry Smith double inc; 3382d3f70b5SBarry Smith 3392d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 3402d3f70b5SBarry Smith 3412d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 3422d3f70b5SBarry Smith if (flg) { 343*28aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 344*28aa8177SBarry Smith } 345*28aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 346*28aa8177SBarry Smith if (flg) { 347*28aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 3482d3f70b5SBarry Smith } 3492d3f70b5SBarry Smith return 0; 3502d3f70b5SBarry Smith } 3512d3f70b5SBarry Smith 3527bf11e45SBarry Smith static int TSPrintHelp_Pseudo(TS ts) 3532d3f70b5SBarry Smith { 3542d3f70b5SBarry Smith 3552d3f70b5SBarry Smith return 0; 3562d3f70b5SBarry Smith } 3572d3f70b5SBarry Smith 3587bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer) 3592d3f70b5SBarry Smith { 3602d3f70b5SBarry Smith return 0; 3612d3f70b5SBarry Smith } 3622d3f70b5SBarry Smith 3632d3f70b5SBarry Smith /* ------------------------------------------------------------ */ 3647bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 3652d3f70b5SBarry Smith { 3667bf11e45SBarry Smith TS_Pseudo *pseudo; 3672d3f70b5SBarry Smith int ierr; 3682d3f70b5SBarry Smith MatType mtype; 3692d3f70b5SBarry Smith 3707bf11e45SBarry Smith ts->type = TS_PSEUDO; 3717bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 3727bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 3737bf11e45SBarry Smith ts->view = TSView_Pseudo; 3742d3f70b5SBarry Smith 3752d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 3767bf11e45SBarry Smith SETERRQ(1,"TSCreate_Pseudo:Only for nonlinear problems"); 3772d3f70b5SBarry Smith } 3782d3f70b5SBarry Smith if (!ts->A) { 3797bf11e45SBarry Smith SETERRQ(1,"TSCreate_Pseudo:Must set Jacobian"); 3802d3f70b5SBarry Smith } 3812d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 3822d3f70b5SBarry Smith if (mtype == MATSHELL) { 3832d3f70b5SBarry Smith ts->Ashell = ts->A; 3842d3f70b5SBarry Smith } 3857bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 3867bf11e45SBarry Smith ts->step = TSStep_Pseudo; 3877bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 3887bf11e45SBarry Smith 3897bf11e45SBarry Smith /* create the required nonlinear solver context */ 3902d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 3912d3f70b5SBarry Smith 3927bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 3937bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 3947bf11e45SBarry Smith ts->data = (void *) pseudo; 3952d3f70b5SBarry Smith 396*28aa8177SBarry Smith pseudo->dt_increment = 1.1; 397*28aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 3982d3f70b5SBarry Smith return 0; 3992d3f70b5SBarry Smith } 4002d3f70b5SBarry Smith 4012d3f70b5SBarry Smith 402*28aa8177SBarry Smith /*@ 403*28aa8177SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 404*28aa8177SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 405*28aa8177SBarry Smith 406*28aa8177SBarry Smith Input Parameters: 407*28aa8177SBarry Smith . ts - the timestep context 408*28aa8177SBarry Smith . inc - the scaling factor >= 1.0 409*28aa8177SBarry Smith 410*28aa8177SBarry Smith Options Database Key: 411*28aa8177SBarry Smith $ -ts_pseudo_increment <increment> 412*28aa8177SBarry Smith 413*28aa8177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 414*28aa8177SBarry Smith @*/ 415*28aa8177SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 416*28aa8177SBarry Smith { 417*28aa8177SBarry Smith TS_Pseudo *pseudo; 418*28aa8177SBarry Smith 419*28aa8177SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 420*28aa8177SBarry Smith if (ts->type != TS_PSEUDO) return 0; 421*28aa8177SBarry Smith 422*28aa8177SBarry Smith pseudo = (TS_Pseudo*) ts->data; 423*28aa8177SBarry Smith pseudo->dt_increment = inc; 424*28aa8177SBarry Smith return 0; 425*28aa8177SBarry Smith } 4262d3f70b5SBarry Smith 4272d3f70b5SBarry Smith 4282d3f70b5SBarry Smith 429