12d3f70b5SBarry Smith #ifndef lint 2*7a00f4a9SLois Curfman McInnes static char vcid[] = "$Id: posindep.c,v 1.10 1997/01/06 20:42:03 bsmith Exp curfman $"; 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) */ 2528aa8177SBarry Smith 2628aa8177SBarry Smith double dt_increment; /* scaling that dt is incremented each time-step */ 277bf11e45SBarry Smith } TS_Pseudo; 282d3f70b5SBarry Smith 292d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 3058562591SBarry Smith #undef __FUNC__ 3158562591SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep" 322d3f70b5SBarry Smith /*@C 33fb4a63b6SLois Curfman McInnes TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 34fb4a63b6SLois Curfman McInnes Use with TSPseudoSetTimeStep(). 352d3f70b5SBarry Smith 362d3f70b5SBarry Smith Input Parameters: 372d3f70b5SBarry Smith . ts - the timestep context 387bf11e45SBarry Smith . dtctx - unused timestep context 392d3f70b5SBarry Smith 402d3f70b5SBarry Smith Output Parameter: 412d3f70b5SBarry Smith . newdt - the timestep to use for the next step 422d3f70b5SBarry Smith 43fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, default 44fb4a63b6SLois Curfman McInnes 45564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 462d3f70b5SBarry Smith @*/ 477bf11e45SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 482d3f70b5SBarry Smith { 497bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 5028aa8177SBarry Smith double inc = pseudo->dt_increment; 512d3f70b5SBarry Smith int ierr; 522d3f70b5SBarry Smith 537bf11e45SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 547bf11e45SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 557bf11e45SBarry Smith if (pseudo->initial_fnorm == 0.0) { 562d3f70b5SBarry Smith /* first time through so compute initial function norm */ 577bf11e45SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 582d3f70b5SBarry Smith } 5928aa8177SBarry Smith if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 6028aa8177SBarry Smith else *newdt = inc*ts->time_step*pseudo->initial_fnorm/pseudo->fnorm; 612d3f70b5SBarry Smith return 0; 622d3f70b5SBarry Smith } 632d3f70b5SBarry Smith 6458562591SBarry Smith #undef __FUNC__ 6558562591SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep" 662d3f70b5SBarry Smith /*@ 67fb4a63b6SLois Curfman McInnes TSPseudoSetTimeStep - Sets the user-defined routine to be 68fb4a63b6SLois Curfman McInnes called at each pseudo-timestep to update the timestep. 692d3f70b5SBarry Smith 702d3f70b5SBarry Smith Input Parameters: 712d3f70b5SBarry Smith . ts - timestep context 722d3f70b5SBarry Smith . dt - function to compute timestep 73564e8f4eSLois Curfman McInnes . ctx - [optional] user-defined context for private data 74564e8f4eSLois Curfman McInnes required by the function (may be PETSC_NULL) 75564e8f4eSLois Curfman McInnes 76564e8f4eSLois Curfman McInnes Calling sequence of func: 77564e8f4eSLois Curfman McInnes . func (TS ts,double *newdt,void *ctx); 78564e8f4eSLois Curfman McInnes 79564e8f4eSLois Curfman McInnes . newdt - the newly computed timestep 80564e8f4eSLois Curfman McInnes . ctx - [optional] timestep context 81564e8f4eSLois Curfman McInnes 82564e8f4eSLois Curfman McInnes Notes: 83564e8f4eSLois Curfman McInnes The routine set here will be called by TSPseudoComputeTimeStep() 84564e8f4eSLois Curfman McInnes during the timestepping process. 852d3f70b5SBarry Smith 86fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, set 87fb4a63b6SLois Curfman McInnes 88564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 892d3f70b5SBarry Smith @*/ 907bf11e45SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 912d3f70b5SBarry Smith { 927bf11e45SBarry Smith TS_Pseudo *pseudo; 932d3f70b5SBarry Smith 942d3f70b5SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 957bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 967bf11e45SBarry Smith 977bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 987bf11e45SBarry Smith pseudo->dt = dt; 997bf11e45SBarry Smith pseudo->dtctx = ctx; 1002d3f70b5SBarry Smith return 0; 1012d3f70b5SBarry Smith } 1022d3f70b5SBarry Smith 10358562591SBarry Smith #undef __FUNC__ 10458562591SBarry Smith #define __FUNC__ "TSPseudoComputeTimeStep" 1057bf11e45SBarry Smith /*@ 1067bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 107564e8f4eSLois Curfman McInnes pseudo-timestepping process. 1082d3f70b5SBarry Smith 1097bf11e45SBarry Smith Input Parameter: 1107bf11e45SBarry Smith . ts - timestep context 1117bf11e45SBarry Smith 1127bf11e45SBarry Smith Output Parameter: 113fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 114fb4a63b6SLois Curfman McInnes 115564e8f4eSLois Curfman McInnes 116564e8f4eSLois Curfman McInnes Notes: 117564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 118564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 119564e8f4eSLois Curfman McInnes 120fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 121564e8f4eSLois Curfman McInnes 122564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 1237bf11e45SBarry Smith @*/ 1247bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 1257bf11e45SBarry Smith { 1267bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1277bf11e45SBarry Smith int ierr; 1287bf11e45SBarry Smith 1297bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 1307bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr); 1317bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 1327bf11e45SBarry Smith return 0; 1337bf11e45SBarry Smith } 1347bf11e45SBarry Smith 1357bf11e45SBarry Smith 1367bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 13758562591SBarry Smith #undef __FUNC__ 13858562591SBarry Smith #define __FUNC__ "TSPseudoDefaultVerifyTimeStep" 1397bf11e45SBarry Smith /*@C 140639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 1417bf11e45SBarry Smith 1427bf11e45SBarry Smith Input Parameters: 1437bf11e45SBarry Smith . ts - the timestep context 1447bf11e45SBarry Smith . dtctx - unused timestep context 145564e8f4eSLois Curfman McInnes . update - latest solution vector 1467bf11e45SBarry Smith 147564e8f4eSLois Curfman McInnes Output Parameters: 1487bf11e45SBarry Smith . newdt - the timestep to use for the next step 149564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 1507bf11e45SBarry Smith 151564e8f4eSLois Curfman McInnes Note: 152564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 153564e8f4eSLois Curfman McInnes timestep. 154564e8f4eSLois Curfman McInnes 155564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 156564e8f4eSLois Curfman McInnes 157564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 1587bf11e45SBarry Smith @*/ 1597bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag) 1607bf11e45SBarry Smith { 1617bf11e45SBarry Smith *flag = 1; 1627bf11e45SBarry Smith return 0; 1637bf11e45SBarry Smith } 1647bf11e45SBarry Smith 1657bf11e45SBarry Smith /*@ 166564e8f4eSLois Curfman McInnes TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 167564e8f4eSLois Curfman McInnes last timestep. 1687bf11e45SBarry Smith 1697bf11e45SBarry Smith Input Parameters: 1707bf11e45SBarry Smith . ts - timestep context 171564e8f4eSLois Curfman McInnes . dt - user-defined function to verify timestep 172564e8f4eSLois Curfman McInnes . ctx - [optional] user-defined context for private data 173564e8f4eSLois Curfman McInnes for the timestep verification routine (may be PETSC_NULL) 1747bf11e45SBarry Smith 175564e8f4eSLois Curfman McInnes Calling sequence of func: 176564e8f4eSLois Curfman McInnes . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 177564e8f4eSLois Curfman McInnes 178564e8f4eSLois Curfman McInnes . update - latest solution vector 179564e8f4eSLois Curfman McInnes . ctx - [optional] timestep context 180564e8f4eSLois Curfman McInnes . newdt - the timestep to use for the next step 181564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 182564e8f4eSLois Curfman McInnes 183564e8f4eSLois Curfman McInnes Notes: 184564e8f4eSLois Curfman McInnes The routine set here will be called by TSPseudoVerifyTimeStep() 185564e8f4eSLois Curfman McInnes during the timestepping process. 186564e8f4eSLois Curfman McInnes 187564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, verify 188564e8f4eSLois Curfman McInnes 189564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 1907bf11e45SBarry Smith @*/ 1917bf11e45SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 1927bf11e45SBarry Smith { 1937bf11e45SBarry Smith TS_Pseudo *pseudo; 1947bf11e45SBarry Smith 1957bf11e45SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 1967bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 1977bf11e45SBarry Smith 1987bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 1997bf11e45SBarry Smith pseudo->verify = dt; 2007bf11e45SBarry Smith pseudo->verifyctx = ctx; 2017bf11e45SBarry Smith return 0; 2027bf11e45SBarry Smith } 2037bf11e45SBarry Smith 20458562591SBarry Smith #undef __FUNC__ 20558562591SBarry Smith #define __FUNC__ "TSPseudoVerifyTimeStep" 2067bf11e45SBarry Smith /*@ 207564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 2087bf11e45SBarry Smith 209fb4a63b6SLois Curfman McInnes Input Parameters: 2107bf11e45SBarry Smith . ts - timestep context 211564e8f4eSLois Curfman McInnes . update - latest solution vector 2127bf11e45SBarry Smith 213fb4a63b6SLois Curfman McInnes Output Parameters: 214fb4a63b6SLois Curfman McInnes . dt - newly computed timestep (if it had to shrink) 2157bf11e45SBarry Smith . flag - indicates if current timestep was ok 2167bf11e45SBarry Smith 217564e8f4eSLois Curfman McInnes Notes: 218564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 219564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 220564e8f4eSLois Curfman McInnes 221564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 222564e8f4eSLois Curfman McInnes 223564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 2247bf11e45SBarry Smith @*/ 2257bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 2267bf11e45SBarry Smith { 2277bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2287bf11e45SBarry Smith int ierr; 2297bf11e45SBarry Smith 2307bf11e45SBarry Smith if (!pseudo->verify) {*flag = 1; return 0;} 2317bf11e45SBarry Smith 2327bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 2337bf11e45SBarry Smith 2347bf11e45SBarry Smith return 0; 2357bf11e45SBarry Smith } 2367bf11e45SBarry Smith 2377bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 2387bf11e45SBarry Smith 23958562591SBarry Smith #undef __FUNC__ 24058562591SBarry Smith #define __FUNC__ "TSStep_Pseudo" 2417bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 2422d3f70b5SBarry Smith { 2432d3f70b5SBarry Smith Vec sol = ts->vec_sol; 2447bf11e45SBarry Smith int ierr,i,max_steps = ts->max_steps,its,ok; 2457bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2467bf11e45SBarry Smith double current_time_step; 2472d3f70b5SBarry Smith 2482d3f70b5SBarry Smith *steps = -ts->steps; 2492d3f70b5SBarry Smith 2507bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 2512d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 2527bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 2537bf11e45SBarry Smith current_time_step = ts->time_step; 2547bf11e45SBarry Smith while (1) { 2557bf11e45SBarry Smith ts->ptime += current_time_step; 2567bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 257*7a00f4a9SLois Curfman McInnes ts->nonlinear_its += PetscAbsInt(its); 2587bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 2597bf11e45SBarry Smith if (ok) break; 2607bf11e45SBarry Smith ts->ptime -= current_time_step; 2617bf11e45SBarry Smith current_time_step = ts->time_step; 2627bf11e45SBarry Smith } 2637bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 2642d3f70b5SBarry Smith ts->steps++; 2652d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 2662d3f70b5SBarry Smith } 2672d3f70b5SBarry Smith 2682d3f70b5SBarry Smith *steps += ts->steps; 2692d3f70b5SBarry Smith *time = ts->ptime; 2702d3f70b5SBarry Smith return 0; 2712d3f70b5SBarry Smith } 2722d3f70b5SBarry Smith 2732d3f70b5SBarry Smith /*------------------------------------------------------------*/ 27458562591SBarry Smith #undef __FUNC__ 27558562591SBarry Smith #define __FUNC__ "TSDestroy_Pseudo" 2767bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj ) 2772d3f70b5SBarry Smith { 2782d3f70b5SBarry Smith TS ts = (TS) obj; 2797bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2802d3f70b5SBarry Smith int ierr; 2812d3f70b5SBarry Smith 2827bf11e45SBarry Smith ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 2837bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 2847bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 2852d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 2867bf11e45SBarry Smith PetscFree(pseudo); 2872d3f70b5SBarry Smith return 0; 2882d3f70b5SBarry Smith } 2892d3f70b5SBarry Smith 2902d3f70b5SBarry Smith 2912d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2922d3f70b5SBarry Smith /* 2932d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 2942d3f70b5SBarry Smith */ 2952d3f70b5SBarry Smith 29658562591SBarry Smith #undef __FUNC__ 29758562591SBarry Smith #define __FUNC__ "TSPseudoMatMult" 2987bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 2992d3f70b5SBarry Smith { 3002d3f70b5SBarry Smith TS ts; 3012d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 3022d3f70b5SBarry Smith int ierr; 3032d3f70b5SBarry Smith 3042d3f70b5SBarry Smith MatShellGetContext(mat,(void **)&ts); 3052d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 3062d3f70b5SBarry Smith 3072d3f70b5SBarry Smith /* apply user provided function */ 3082d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 3092d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 3102d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 3112d3f70b5SBarry Smith return 0; 3122d3f70b5SBarry Smith } 3132d3f70b5SBarry Smith 3142d3f70b5SBarry Smith /* 3152d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 3162d3f70b5SBarry Smith 3172d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 3182d3f70b5SBarry Smith */ 31958562591SBarry Smith #undef __FUNC__ 32058562591SBarry Smith #define __FUNC__ "TSPseudoFunction" 3217bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 3222d3f70b5SBarry Smith { 3232d3f70b5SBarry Smith TS ts = (TS) ctx; 3242d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 3252d3f70b5SBarry Smith int ierr,i,n; 3262d3f70b5SBarry Smith 3272d3f70b5SBarry Smith /* apply user provided function */ 3282d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 3297bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 3302d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 3312d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 3322d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 3332d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 3342d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 3352d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 3362d3f70b5SBarry Smith } 3372d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 3382d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 3392d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 3402d3f70b5SBarry Smith 3412d3f70b5SBarry Smith return 0; 3422d3f70b5SBarry Smith } 3432d3f70b5SBarry Smith 3442d3f70b5SBarry Smith /* 3452d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 3462d3f70b5SBarry Smith 3472d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 3482d3f70b5SBarry Smith */ 34958562591SBarry Smith #undef __FUNC__ 35058562591SBarry Smith #define __FUNC__ "TSPseudoJacobian" 3517bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 3522d3f70b5SBarry Smith { 3532d3f70b5SBarry Smith TS ts = (TS) ctx; 3542d3f70b5SBarry Smith int ierr; 3552d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 3562d3f70b5SBarry Smith MatType mtype; 3572d3f70b5SBarry Smith 3582d3f70b5SBarry Smith /* construct users Jacobian */ 3592d3f70b5SBarry Smith if (ts->rhsjacobian) { 3602d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 3612d3f70b5SBarry Smith } 3622d3f70b5SBarry Smith 3632d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 3642d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 3652d3f70b5SBarry Smith if (mtype != MATSHELL) { 3662d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 3672d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 3682d3f70b5SBarry Smith } 3692d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 3702d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 3712d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 3722d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 3732d3f70b5SBarry Smith } 3742d3f70b5SBarry Smith 3752d3f70b5SBarry Smith return 0; 3762d3f70b5SBarry Smith } 3772d3f70b5SBarry Smith 3782d3f70b5SBarry Smith 37958562591SBarry Smith #undef __FUNC__ 38058562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo" 3817bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 3822d3f70b5SBarry Smith { 3837bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3842d3f70b5SBarry Smith int ierr, M, m; 3852d3f70b5SBarry Smith 3867bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 3877bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 3887bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 3892d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 3902d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 3912d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 3922d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 3937bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 3942d3f70b5SBarry Smith } 3957bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 3962d3f70b5SBarry Smith return 0; 3972d3f70b5SBarry Smith } 3982d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3992d3f70b5SBarry Smith 40058562591SBarry Smith #undef __FUNC__ 40158562591SBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor" 4022d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 4032d3f70b5SBarry Smith { 4047bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 4052d3f70b5SBarry Smith 4067bf11e45SBarry Smith PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 4072d3f70b5SBarry Smith return 0; 4082d3f70b5SBarry Smith } 4092d3f70b5SBarry Smith 41058562591SBarry Smith #undef __FUNC__ 41158562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo" 4127bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 4132d3f70b5SBarry Smith { 4142d3f70b5SBarry Smith int ierr,flg; 41528aa8177SBarry Smith double inc; 4162d3f70b5SBarry Smith 4172d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 4182d3f70b5SBarry Smith 4192d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 4202d3f70b5SBarry Smith if (flg) { 42128aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 42228aa8177SBarry Smith } 42328aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 42428aa8177SBarry Smith if (flg) { 42528aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 4262d3f70b5SBarry Smith } 4272d3f70b5SBarry Smith return 0; 4282d3f70b5SBarry Smith } 4292d3f70b5SBarry Smith 43058562591SBarry Smith #undef __FUNC__ 43158562591SBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo" 4327bf11e45SBarry Smith static int TSPrintHelp_Pseudo(TS ts) 4332d3f70b5SBarry Smith { 4342d3f70b5SBarry Smith return 0; 4352d3f70b5SBarry Smith } 4362d3f70b5SBarry Smith 43758562591SBarry Smith #undef __FUNC__ 43858562591SBarry Smith #define __FUNC__ "TSView_Pseudo" 4397bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer) 4402d3f70b5SBarry Smith { 4412d3f70b5SBarry Smith return 0; 4422d3f70b5SBarry Smith } 4432d3f70b5SBarry Smith 4442d3f70b5SBarry Smith /* ------------------------------------------------------------ */ 44558562591SBarry Smith #undef __FUNC__ 44658562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo" 4477bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 4482d3f70b5SBarry Smith { 4497bf11e45SBarry Smith TS_Pseudo *pseudo; 4502d3f70b5SBarry Smith int ierr; 4512d3f70b5SBarry Smith MatType mtype; 4522d3f70b5SBarry Smith 4537bf11e45SBarry Smith ts->type = TS_PSEUDO; 4547bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 4557bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 4567bf11e45SBarry Smith ts->view = TSView_Pseudo; 4572d3f70b5SBarry Smith 4582d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 459e3372554SBarry Smith SETERRQ(1,0,"Only for nonlinear problems"); 4602d3f70b5SBarry Smith } 4612d3f70b5SBarry Smith if (!ts->A) { 462e3372554SBarry Smith SETERRQ(1,0,"Must set Jacobian"); 4632d3f70b5SBarry Smith } 4642d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 4652d3f70b5SBarry Smith if (mtype == MATSHELL) { 4662d3f70b5SBarry Smith ts->Ashell = ts->A; 4672d3f70b5SBarry Smith } 4687bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 4697bf11e45SBarry Smith ts->step = TSStep_Pseudo; 4707bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 4717bf11e45SBarry Smith 4727bf11e45SBarry Smith /* create the required nonlinear solver context */ 4732d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 4742d3f70b5SBarry Smith 4757bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 4767bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 4777bf11e45SBarry Smith ts->data = (void *) pseudo; 4782d3f70b5SBarry Smith 47928aa8177SBarry Smith pseudo->dt_increment = 1.1; 48028aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 4812d3f70b5SBarry Smith return 0; 4822d3f70b5SBarry Smith } 4832d3f70b5SBarry Smith 4842d3f70b5SBarry Smith 48558562591SBarry Smith #undef __FUNC__ 48658562591SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement" 48728aa8177SBarry Smith /*@ 48828aa8177SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 48928aa8177SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 49028aa8177SBarry Smith 49128aa8177SBarry Smith Input Parameters: 49228aa8177SBarry Smith . ts - the timestep context 49328aa8177SBarry Smith . inc - the scaling factor >= 1.0 49428aa8177SBarry Smith 49528aa8177SBarry Smith Options Database Key: 49628aa8177SBarry Smith $ -ts_pseudo_increment <increment> 49728aa8177SBarry Smith 498564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, increment 499564e8f4eSLois Curfman McInnes 50028aa8177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 50128aa8177SBarry Smith @*/ 50228aa8177SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 50328aa8177SBarry Smith { 50428aa8177SBarry Smith TS_Pseudo *pseudo; 50528aa8177SBarry Smith 50628aa8177SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 50728aa8177SBarry Smith if (ts->type != TS_PSEUDO) return 0; 50828aa8177SBarry Smith 50928aa8177SBarry Smith pseudo = (TS_Pseudo*) ts->data; 51028aa8177SBarry Smith pseudo->dt_increment = inc; 51128aa8177SBarry Smith return 0; 51228aa8177SBarry Smith } 5132d3f70b5SBarry Smith 5142d3f70b5SBarry Smith 5152d3f70b5SBarry Smith 516