12d3f70b5SBarry Smith #ifndef lint 2*b5ed506cSLois Curfman McInnes static char vcid[] = "$Id: posindep.c,v 1.12 1997/01/21 18:42:10 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) */ 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 /* ------------------------------------------------------------------------------*/ 3258562591SBarry Smith #undef __FUNC__ 3358562591SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep" 342d3f70b5SBarry Smith /*@C 35fb4a63b6SLois Curfman McInnes TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 36fb4a63b6SLois Curfman McInnes Use with TSPseudoSetTimeStep(). 372d3f70b5SBarry Smith 382d3f70b5SBarry Smith Input Parameters: 392d3f70b5SBarry Smith . ts - the timestep context 407bf11e45SBarry Smith . dtctx - unused timestep context 412d3f70b5SBarry Smith 422d3f70b5SBarry Smith Output Parameter: 432d3f70b5SBarry Smith . newdt - the timestep to use for the next step 442d3f70b5SBarry Smith 45fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, default 46fb4a63b6SLois Curfman McInnes 47564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 482d3f70b5SBarry Smith @*/ 497bf11e45SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 502d3f70b5SBarry Smith { 517bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 52ca4b7087SBarry Smith double inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 532d3f70b5SBarry Smith int ierr; 542d3f70b5SBarry Smith 557bf11e45SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 567bf11e45SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 577bf11e45SBarry Smith if (pseudo->initial_fnorm == 0.0) { 582d3f70b5SBarry Smith /* first time through so compute initial function norm */ 597bf11e45SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 60ca4b7087SBarry Smith fnorm_previous = pseudo->fnorm; 612d3f70b5SBarry Smith } 62ca4b7087SBarry Smith if (pseudo->fnorm == 0.0) { 63ca4b7087SBarry Smith *newdt = 1.e12*inc*ts->time_step; 64ca4b7087SBarry Smith } 65ca4b7087SBarry Smith else if (pseudo->increment_dt_from_initial_dt) { 66ca4b7087SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 67ca4b7087SBarry Smith } else { 68ca4b7087SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 69ca4b7087SBarry Smith } 70ca4b7087SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 712d3f70b5SBarry Smith return 0; 722d3f70b5SBarry Smith } 732d3f70b5SBarry Smith 7458562591SBarry Smith #undef __FUNC__ 7558562591SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep" 762d3f70b5SBarry Smith /*@ 77fb4a63b6SLois Curfman McInnes TSPseudoSetTimeStep - Sets the user-defined routine to be 78fb4a63b6SLois Curfman McInnes called at each pseudo-timestep to update the timestep. 792d3f70b5SBarry Smith 802d3f70b5SBarry Smith Input Parameters: 812d3f70b5SBarry Smith . ts - timestep context 822d3f70b5SBarry Smith . dt - function to compute timestep 83564e8f4eSLois Curfman McInnes . ctx - [optional] user-defined context for private data 84564e8f4eSLois Curfman McInnes required by the function (may be PETSC_NULL) 85564e8f4eSLois Curfman McInnes 86564e8f4eSLois Curfman McInnes Calling sequence of func: 87564e8f4eSLois Curfman McInnes . func (TS ts,double *newdt,void *ctx); 88564e8f4eSLois Curfman McInnes 89564e8f4eSLois Curfman McInnes . newdt - the newly computed timestep 90564e8f4eSLois Curfman McInnes . ctx - [optional] timestep context 91564e8f4eSLois Curfman McInnes 92564e8f4eSLois Curfman McInnes Notes: 93564e8f4eSLois Curfman McInnes The routine set here will be called by TSPseudoComputeTimeStep() 94564e8f4eSLois Curfman McInnes during the timestepping process. 952d3f70b5SBarry Smith 96fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, set 97fb4a63b6SLois Curfman McInnes 98564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 992d3f70b5SBarry Smith @*/ 1007bf11e45SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 1012d3f70b5SBarry Smith { 1027bf11e45SBarry Smith TS_Pseudo *pseudo; 1032d3f70b5SBarry Smith 1042d3f70b5SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 1057bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 1067bf11e45SBarry Smith 1077bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 1087bf11e45SBarry Smith pseudo->dt = dt; 1097bf11e45SBarry Smith pseudo->dtctx = ctx; 1102d3f70b5SBarry Smith return 0; 1112d3f70b5SBarry Smith } 1122d3f70b5SBarry Smith 11358562591SBarry Smith #undef __FUNC__ 11458562591SBarry Smith #define __FUNC__ "TSPseudoComputeTimeStep" 1157bf11e45SBarry Smith /*@ 1167bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 117564e8f4eSLois Curfman McInnes pseudo-timestepping process. 1182d3f70b5SBarry Smith 1197bf11e45SBarry Smith Input Parameter: 1207bf11e45SBarry Smith . ts - timestep context 1217bf11e45SBarry Smith 1227bf11e45SBarry Smith Output Parameter: 123fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 124fb4a63b6SLois Curfman McInnes 125564e8f4eSLois Curfman McInnes 126564e8f4eSLois Curfman McInnes Notes: 127564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 128564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 129564e8f4eSLois Curfman McInnes 130fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 131564e8f4eSLois Curfman McInnes 132564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 1337bf11e45SBarry Smith @*/ 1347bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 1357bf11e45SBarry Smith { 1367bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1377bf11e45SBarry Smith int ierr; 1387bf11e45SBarry Smith 1397bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 1407bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr); 1417bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 1427bf11e45SBarry Smith return 0; 1437bf11e45SBarry Smith } 1447bf11e45SBarry Smith 1457bf11e45SBarry Smith 1467bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 14758562591SBarry Smith #undef __FUNC__ 14858562591SBarry Smith #define __FUNC__ "TSPseudoDefaultVerifyTimeStep" 1497bf11e45SBarry Smith /*@C 150639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 1517bf11e45SBarry Smith 1527bf11e45SBarry Smith Input Parameters: 1537bf11e45SBarry Smith . ts - the timestep context 1547bf11e45SBarry Smith . dtctx - unused timestep context 155564e8f4eSLois Curfman McInnes . update - latest solution vector 1567bf11e45SBarry Smith 157564e8f4eSLois Curfman McInnes Output Parameters: 1587bf11e45SBarry Smith . newdt - the timestep to use for the next step 159564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 1607bf11e45SBarry Smith 161564e8f4eSLois Curfman McInnes Note: 162564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 163564e8f4eSLois Curfman McInnes timestep. 164564e8f4eSLois Curfman McInnes 165564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 166564e8f4eSLois Curfman McInnes 167564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 1687bf11e45SBarry Smith @*/ 1697bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag) 1707bf11e45SBarry Smith { 1717bf11e45SBarry Smith *flag = 1; 1727bf11e45SBarry Smith return 0; 1737bf11e45SBarry Smith } 1747bf11e45SBarry Smith 1757bf11e45SBarry Smith /*@ 176564e8f4eSLois Curfman McInnes TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 177564e8f4eSLois Curfman McInnes last timestep. 1787bf11e45SBarry Smith 1797bf11e45SBarry Smith Input Parameters: 1807bf11e45SBarry Smith . ts - timestep context 181564e8f4eSLois Curfman McInnes . dt - user-defined function to verify timestep 182564e8f4eSLois Curfman McInnes . ctx - [optional] user-defined context for private data 183564e8f4eSLois Curfman McInnes for the timestep verification routine (may be PETSC_NULL) 1847bf11e45SBarry Smith 185564e8f4eSLois Curfman McInnes Calling sequence of func: 186564e8f4eSLois Curfman McInnes . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 187564e8f4eSLois Curfman McInnes 188564e8f4eSLois Curfman McInnes . update - latest solution vector 189564e8f4eSLois Curfman McInnes . ctx - [optional] timestep context 190564e8f4eSLois Curfman McInnes . newdt - the timestep to use for the next step 191564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 192564e8f4eSLois Curfman McInnes 193564e8f4eSLois Curfman McInnes Notes: 194564e8f4eSLois Curfman McInnes The routine set here will be called by TSPseudoVerifyTimeStep() 195564e8f4eSLois Curfman McInnes during the timestepping process. 196564e8f4eSLois Curfman McInnes 197564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, verify 198564e8f4eSLois Curfman McInnes 199564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 2007bf11e45SBarry Smith @*/ 2017bf11e45SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 2027bf11e45SBarry Smith { 2037bf11e45SBarry Smith TS_Pseudo *pseudo; 2047bf11e45SBarry Smith 2057bf11e45SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 2067bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 2077bf11e45SBarry Smith 2087bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 2097bf11e45SBarry Smith pseudo->verify = dt; 2107bf11e45SBarry Smith pseudo->verifyctx = ctx; 2117bf11e45SBarry Smith return 0; 2127bf11e45SBarry Smith } 2137bf11e45SBarry Smith 21458562591SBarry Smith #undef __FUNC__ 21558562591SBarry Smith #define __FUNC__ "TSPseudoVerifyTimeStep" 2167bf11e45SBarry Smith /*@ 217564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 2187bf11e45SBarry Smith 219fb4a63b6SLois Curfman McInnes Input Parameters: 2207bf11e45SBarry Smith . ts - timestep context 221564e8f4eSLois Curfman McInnes . update - latest solution vector 2227bf11e45SBarry Smith 223fb4a63b6SLois Curfman McInnes Output Parameters: 224fb4a63b6SLois Curfman McInnes . dt - newly computed timestep (if it had to shrink) 2257bf11e45SBarry Smith . flag - indicates if current timestep was ok 2267bf11e45SBarry Smith 227564e8f4eSLois Curfman McInnes Notes: 228564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 229564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 230564e8f4eSLois Curfman McInnes 231564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 232564e8f4eSLois Curfman McInnes 233564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 2347bf11e45SBarry Smith @*/ 2357bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 2367bf11e45SBarry Smith { 2377bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2387bf11e45SBarry Smith int ierr; 2397bf11e45SBarry Smith 2407bf11e45SBarry Smith if (!pseudo->verify) {*flag = 1; return 0;} 2417bf11e45SBarry Smith 2427bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 2437bf11e45SBarry Smith 2447bf11e45SBarry Smith return 0; 2457bf11e45SBarry Smith } 2467bf11e45SBarry Smith 2477bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 2487bf11e45SBarry Smith 24958562591SBarry Smith #undef __FUNC__ 25058562591SBarry Smith #define __FUNC__ "TSStep_Pseudo" 2517bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 2522d3f70b5SBarry Smith { 2532d3f70b5SBarry Smith Vec sol = ts->vec_sol; 2547bf11e45SBarry Smith int ierr,i,max_steps = ts->max_steps,its,ok; 2557bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2567bf11e45SBarry Smith double current_time_step; 2572d3f70b5SBarry Smith 2582d3f70b5SBarry Smith *steps = -ts->steps; 2592d3f70b5SBarry Smith 2607bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 2612d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 2627bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 2637bf11e45SBarry Smith current_time_step = ts->time_step; 2647bf11e45SBarry Smith while (1) { 2657bf11e45SBarry Smith ts->ptime += current_time_step; 2667bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 2677a00f4a9SLois Curfman McInnes ts->nonlinear_its += PetscAbsInt(its); 2687bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 2697bf11e45SBarry Smith if (ok) break; 2707bf11e45SBarry Smith ts->ptime -= current_time_step; 2717bf11e45SBarry Smith current_time_step = ts->time_step; 2727bf11e45SBarry Smith } 2737bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 2742d3f70b5SBarry Smith ts->steps++; 2752d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 2762d3f70b5SBarry Smith } 2772d3f70b5SBarry Smith 2782d3f70b5SBarry Smith *steps += ts->steps; 2792d3f70b5SBarry Smith *time = ts->ptime; 2802d3f70b5SBarry Smith return 0; 2812d3f70b5SBarry Smith } 2822d3f70b5SBarry Smith 2832d3f70b5SBarry Smith /*------------------------------------------------------------*/ 28458562591SBarry Smith #undef __FUNC__ 28558562591SBarry Smith #define __FUNC__ "TSDestroy_Pseudo" 2867bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj ) 2872d3f70b5SBarry Smith { 2882d3f70b5SBarry Smith TS ts = (TS) obj; 2897bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2902d3f70b5SBarry Smith int ierr; 2912d3f70b5SBarry Smith 2927bf11e45SBarry Smith ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 2937bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 2947bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 2952d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 2967bf11e45SBarry Smith PetscFree(pseudo); 2972d3f70b5SBarry Smith return 0; 2982d3f70b5SBarry Smith } 2992d3f70b5SBarry Smith 3002d3f70b5SBarry Smith 3012d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3022d3f70b5SBarry Smith /* 3032d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 3042d3f70b5SBarry Smith */ 3052d3f70b5SBarry Smith 30658562591SBarry Smith #undef __FUNC__ 30758562591SBarry Smith #define __FUNC__ "TSPseudoMatMult" 3087bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 3092d3f70b5SBarry Smith { 3102d3f70b5SBarry Smith TS ts; 3112d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 3122d3f70b5SBarry Smith int ierr; 3132d3f70b5SBarry Smith 3142d3f70b5SBarry Smith MatShellGetContext(mat,(void **)&ts); 3152d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 3162d3f70b5SBarry Smith 3172d3f70b5SBarry Smith /* apply user provided function */ 3182d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 3192d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 3202d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 3212d3f70b5SBarry Smith return 0; 3222d3f70b5SBarry Smith } 3232d3f70b5SBarry Smith 3242d3f70b5SBarry Smith /* 3252d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 3262d3f70b5SBarry Smith 3272d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 3282d3f70b5SBarry Smith */ 32958562591SBarry Smith #undef __FUNC__ 33058562591SBarry Smith #define __FUNC__ "TSPseudoFunction" 3317bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 3322d3f70b5SBarry Smith { 3332d3f70b5SBarry Smith TS ts = (TS) ctx; 3342d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 3352d3f70b5SBarry Smith int ierr,i,n; 3362d3f70b5SBarry Smith 3372d3f70b5SBarry Smith /* apply user provided function */ 3382d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 3397bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 3402d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 3412d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 3422d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 3432d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 3442d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 3452d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 3462d3f70b5SBarry Smith } 3472d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 3482d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 3492d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 3502d3f70b5SBarry Smith 3512d3f70b5SBarry Smith return 0; 3522d3f70b5SBarry Smith } 3532d3f70b5SBarry Smith 3542d3f70b5SBarry Smith /* 3552d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 3562d3f70b5SBarry Smith 3572d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 3582d3f70b5SBarry Smith */ 35958562591SBarry Smith #undef __FUNC__ 36058562591SBarry Smith #define __FUNC__ "TSPseudoJacobian" 3617bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 3622d3f70b5SBarry Smith { 3632d3f70b5SBarry Smith TS ts = (TS) ctx; 3642d3f70b5SBarry Smith int ierr; 3652d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 3662d3f70b5SBarry Smith MatType mtype; 3672d3f70b5SBarry Smith 3682d3f70b5SBarry Smith /* construct users Jacobian */ 3692d3f70b5SBarry Smith if (ts->rhsjacobian) { 3702d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 3712d3f70b5SBarry Smith } 3722d3f70b5SBarry Smith 3732d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 3742d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 3752d3f70b5SBarry Smith if (mtype != MATSHELL) { 3762d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 3772d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 3782d3f70b5SBarry Smith } 3792d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 3802d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 3812d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 3822d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 3832d3f70b5SBarry Smith } 3842d3f70b5SBarry Smith 3852d3f70b5SBarry Smith return 0; 3862d3f70b5SBarry Smith } 3872d3f70b5SBarry Smith 3882d3f70b5SBarry Smith 38958562591SBarry Smith #undef __FUNC__ 39058562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo" 3917bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 3922d3f70b5SBarry Smith { 3937bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3942d3f70b5SBarry Smith int ierr, M, m; 3952d3f70b5SBarry Smith 3967bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 3977bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 3987bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 3992d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 4002d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 4012d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 4022d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 4037bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 4042d3f70b5SBarry Smith } 4057bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 4062d3f70b5SBarry Smith return 0; 4072d3f70b5SBarry Smith } 4082d3f70b5SBarry Smith /*------------------------------------------------------------*/ 4092d3f70b5SBarry Smith 41058562591SBarry Smith #undef __FUNC__ 41158562591SBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor" 4122d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 4132d3f70b5SBarry Smith { 4147bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 4152d3f70b5SBarry Smith 4167bf11e45SBarry Smith PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 4172d3f70b5SBarry Smith return 0; 4182d3f70b5SBarry Smith } 4192d3f70b5SBarry Smith 42058562591SBarry Smith #undef __FUNC__ 42158562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo" 4227bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 4232d3f70b5SBarry Smith { 4242d3f70b5SBarry Smith int ierr,flg; 42528aa8177SBarry Smith double inc; 4262d3f70b5SBarry Smith 4272d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 4282d3f70b5SBarry Smith 4292d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 4302d3f70b5SBarry Smith if (flg) { 43128aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 43228aa8177SBarry Smith } 43328aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 43428aa8177SBarry Smith if (flg) { 43528aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 4362d3f70b5SBarry Smith } 437ca4b7087SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr); 438ca4b7087SBarry Smith if (flg) { 439ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr); 440ca4b7087SBarry Smith } 4412d3f70b5SBarry Smith return 0; 4422d3f70b5SBarry Smith } 4432d3f70b5SBarry Smith 44458562591SBarry Smith #undef __FUNC__ 44558562591SBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo" 446ca4b7087SBarry Smith static int TSPrintHelp_Pseudo(TS ts,char *p) 4472d3f70b5SBarry Smith { 448ca4b7087SBarry Smith PetscPrintf(ts->comm," Options for TS Pseudo timestepper:\n"); 449ca4b7087SBarry Smith PetscPrintf(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p); 450ca4b7087SBarry Smith PetscPrintf(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p); 451ca4b7087SBarry Smith PetscPrintf(ts->comm," initial fnorm/current fnorm to determine new timestep\n"); 452ca4b7087SBarry Smith PetscPrintf(ts->comm," default is current_dt * previous fnorm/current fnorm\n"); 4532d3f70b5SBarry Smith return 0; 4542d3f70b5SBarry Smith } 4552d3f70b5SBarry Smith 45658562591SBarry Smith #undef __FUNC__ 45758562591SBarry Smith #define __FUNC__ "TSView_Pseudo" 4587bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer) 4592d3f70b5SBarry Smith { 4602d3f70b5SBarry Smith return 0; 4612d3f70b5SBarry Smith } 4622d3f70b5SBarry Smith 4632d3f70b5SBarry Smith /* ------------------------------------------------------------ */ 46458562591SBarry Smith #undef __FUNC__ 46558562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo" 4667bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 4672d3f70b5SBarry Smith { 4687bf11e45SBarry Smith TS_Pseudo *pseudo; 4692d3f70b5SBarry Smith int ierr; 4702d3f70b5SBarry Smith MatType mtype; 4712d3f70b5SBarry Smith 4727bf11e45SBarry Smith ts->type = TS_PSEUDO; 4737bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 4747bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 4757bf11e45SBarry Smith ts->view = TSView_Pseudo; 4762d3f70b5SBarry Smith 4772d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 478e3372554SBarry Smith SETERRQ(1,0,"Only for nonlinear problems"); 4792d3f70b5SBarry Smith } 4802d3f70b5SBarry Smith if (!ts->A) { 481e3372554SBarry Smith SETERRQ(1,0,"Must set Jacobian"); 4822d3f70b5SBarry Smith } 4832d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 4842d3f70b5SBarry Smith if (mtype == MATSHELL) { 4852d3f70b5SBarry Smith ts->Ashell = ts->A; 4862d3f70b5SBarry Smith } 4877bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 4887bf11e45SBarry Smith ts->step = TSStep_Pseudo; 4897bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 4907bf11e45SBarry Smith 4917bf11e45SBarry Smith /* create the required nonlinear solver context */ 4922d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 4932d3f70b5SBarry Smith 4947bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 4957bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 4967bf11e45SBarry Smith ts->data = (void *) pseudo; 4972d3f70b5SBarry Smith 49828aa8177SBarry Smith pseudo->dt_increment = 1.1; 499ca4b7087SBarry Smith pseudo->increment_dt_from_initial_dt = 0; 50028aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 5012d3f70b5SBarry Smith return 0; 5022d3f70b5SBarry Smith } 5032d3f70b5SBarry Smith 5042d3f70b5SBarry Smith 50558562591SBarry Smith #undef __FUNC__ 50658562591SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement" 50728aa8177SBarry Smith /*@ 50828aa8177SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 50928aa8177SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 51028aa8177SBarry Smith 51128aa8177SBarry Smith Input Parameters: 51228aa8177SBarry Smith . ts - the timestep context 51328aa8177SBarry Smith . inc - the scaling factor >= 1.0 51428aa8177SBarry Smith 51528aa8177SBarry Smith Options Database Key: 51628aa8177SBarry Smith $ -ts_pseudo_increment <increment> 51728aa8177SBarry Smith 518564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, increment 519564e8f4eSLois Curfman McInnes 52028aa8177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 52128aa8177SBarry Smith @*/ 52228aa8177SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 52328aa8177SBarry Smith { 52428aa8177SBarry Smith TS_Pseudo *pseudo; 52528aa8177SBarry Smith 52628aa8177SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 52728aa8177SBarry Smith if (ts->type != TS_PSEUDO) return 0; 52828aa8177SBarry Smith 52928aa8177SBarry Smith pseudo = (TS_Pseudo*) ts->data; 53028aa8177SBarry Smith pseudo->dt_increment = inc; 53128aa8177SBarry Smith return 0; 53228aa8177SBarry Smith } 5332d3f70b5SBarry Smith 534ca4b7087SBarry Smith #undef __FUNC__ 535ca4b7087SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt" 536ca4b7087SBarry Smith /*@ 537*b5ed506cSLois Curfman McInnes TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 538*b5ed506cSLois Curfman McInnes is computed via the formula 539*b5ed506cSLois Curfman McInnes $ dt = initial_dt*initial_fnorm/current_fnorm 540*b5ed506cSLois Curfman McInnes rather than the default update, 541*b5ed506cSLois Curfman McInnes $ dt = current_dt*previous_fnorm/current_fnorm. 542ca4b7087SBarry Smith 543ca4b7087SBarry Smith Input Parameters: 544ca4b7087SBarry Smith . ts - the timestep context 545ca4b7087SBarry Smith 546ca4b7087SBarry Smith Options Database Key: 547ca4b7087SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 548ca4b7087SBarry Smith 549ca4b7087SBarry Smith .keywords: timestep, pseudo, set, increment 550ca4b7087SBarry Smith 551ca4b7087SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 552ca4b7087SBarry Smith @*/ 553ca4b7087SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts) 554ca4b7087SBarry Smith { 555ca4b7087SBarry Smith TS_Pseudo *pseudo; 556ca4b7087SBarry Smith 557ca4b7087SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 558ca4b7087SBarry Smith if (ts->type != TS_PSEUDO) return 0; 559ca4b7087SBarry Smith 560ca4b7087SBarry Smith pseudo = (TS_Pseudo*) ts->data; 561ca4b7087SBarry Smith pseudo->increment_dt_from_initial_dt = 1; 562ca4b7087SBarry Smith return 0; 563ca4b7087SBarry Smith } 5642d3f70b5SBarry Smith 5652d3f70b5SBarry Smith 566