12d3f70b5SBarry Smith #ifndef lint 2*639f9d9dSBarry Smith static char vcid[] = "$Id: posindep.c,v 1.5 1996/10/01 15:48:03 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) */ 2528aa8177SBarry Smith 2628aa8177SBarry 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 43564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 442d3f70b5SBarry Smith @*/ 457bf11e45SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 462d3f70b5SBarry Smith { 477bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 4828aa8177SBarry 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 } 5728aa8177SBarry Smith if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 5828aa8177SBarry 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 69564e8f4eSLois Curfman McInnes . ctx - [optional] user-defined context for private data 70564e8f4eSLois Curfman McInnes required by the function (may be PETSC_NULL) 71564e8f4eSLois Curfman McInnes 72564e8f4eSLois Curfman McInnes Calling sequence of func: 73564e8f4eSLois Curfman McInnes . func (TS ts,double *newdt,void *ctx); 74564e8f4eSLois Curfman McInnes 75564e8f4eSLois Curfman McInnes . newdt - the newly computed timestep 76564e8f4eSLois Curfman McInnes . ctx - [optional] timestep context 77564e8f4eSLois Curfman McInnes 78564e8f4eSLois Curfman McInnes Notes: 79564e8f4eSLois Curfman McInnes The routine set here will be called by TSPseudoComputeTimeStep() 80564e8f4eSLois Curfman McInnes during the timestepping process. 812d3f70b5SBarry Smith 82fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, set 83fb4a63b6SLois Curfman McInnes 84564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 852d3f70b5SBarry Smith @*/ 867bf11e45SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 872d3f70b5SBarry Smith { 887bf11e45SBarry Smith TS_Pseudo *pseudo; 892d3f70b5SBarry Smith 902d3f70b5SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 917bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 927bf11e45SBarry Smith 937bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 947bf11e45SBarry Smith pseudo->dt = dt; 957bf11e45SBarry Smith pseudo->dtctx = ctx; 962d3f70b5SBarry Smith return 0; 972d3f70b5SBarry Smith } 982d3f70b5SBarry Smith 997bf11e45SBarry Smith /*@ 1007bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 101564e8f4eSLois Curfman McInnes pseudo-timestepping process. 1022d3f70b5SBarry Smith 1037bf11e45SBarry Smith Input Parameter: 1047bf11e45SBarry Smith . ts - timestep context 1057bf11e45SBarry Smith 1067bf11e45SBarry Smith Output Parameter: 107fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 108fb4a63b6SLois Curfman McInnes 109564e8f4eSLois Curfman McInnes 110564e8f4eSLois Curfman McInnes Notes: 111564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 112564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 113564e8f4eSLois Curfman McInnes 114fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 115564e8f4eSLois Curfman McInnes 116564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 1177bf11e45SBarry Smith @*/ 1187bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 1197bf11e45SBarry Smith { 1207bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1217bf11e45SBarry Smith int ierr; 1227bf11e45SBarry Smith 1237bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 1247bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr); 1257bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 1267bf11e45SBarry Smith return 0; 1277bf11e45SBarry Smith } 1287bf11e45SBarry Smith 1297bf11e45SBarry Smith 1307bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 1317bf11e45SBarry Smith /*@C 132*639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 1337bf11e45SBarry Smith 1347bf11e45SBarry Smith Input Parameters: 1357bf11e45SBarry Smith . ts - the timestep context 1367bf11e45SBarry Smith . dtctx - unused timestep context 137564e8f4eSLois Curfman McInnes . update - latest solution vector 1387bf11e45SBarry Smith 139564e8f4eSLois Curfman McInnes Output Parameters: 1407bf11e45SBarry Smith . newdt - the timestep to use for the next step 141564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 1427bf11e45SBarry Smith 143564e8f4eSLois Curfman McInnes Note: 144564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 145564e8f4eSLois Curfman McInnes timestep. 146564e8f4eSLois Curfman McInnes 147564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 148564e8f4eSLois Curfman McInnes 149564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 1507bf11e45SBarry Smith @*/ 1517bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag) 1527bf11e45SBarry Smith { 1537bf11e45SBarry Smith *flag = 1; 1547bf11e45SBarry Smith return 0; 1557bf11e45SBarry Smith } 1567bf11e45SBarry Smith 1577bf11e45SBarry Smith /*@ 158564e8f4eSLois Curfman McInnes TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 159564e8f4eSLois Curfman McInnes last timestep. 1607bf11e45SBarry Smith 1617bf11e45SBarry Smith Input Parameters: 1627bf11e45SBarry Smith . ts - timestep context 163564e8f4eSLois Curfman McInnes . dt - user-defined function to verify timestep 164564e8f4eSLois Curfman McInnes . ctx - [optional] user-defined context for private data 165564e8f4eSLois Curfman McInnes for the timestep verification routine (may be PETSC_NULL) 1667bf11e45SBarry Smith 167564e8f4eSLois Curfman McInnes Calling sequence of func: 168564e8f4eSLois Curfman McInnes . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 169564e8f4eSLois Curfman McInnes 170564e8f4eSLois Curfman McInnes . update - latest solution vector 171564e8f4eSLois Curfman McInnes . ctx - [optional] timestep context 172564e8f4eSLois Curfman McInnes . newdt - the timestep to use for the next step 173564e8f4eSLois Curfman McInnes . flag - flag indicating whether the last time step was acceptable 174564e8f4eSLois Curfman McInnes 175564e8f4eSLois Curfman McInnes Notes: 176564e8f4eSLois Curfman McInnes The routine set here will be called by TSPseudoVerifyTimeStep() 177564e8f4eSLois Curfman McInnes during the timestepping process. 178564e8f4eSLois Curfman McInnes 179564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, verify 180564e8f4eSLois Curfman McInnes 181564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 1827bf11e45SBarry Smith @*/ 1837bf11e45SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 1847bf11e45SBarry Smith { 1857bf11e45SBarry Smith TS_Pseudo *pseudo; 1867bf11e45SBarry Smith 1877bf11e45SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 1887bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 1897bf11e45SBarry Smith 1907bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 1917bf11e45SBarry Smith pseudo->verify = dt; 1927bf11e45SBarry Smith pseudo->verifyctx = ctx; 1937bf11e45SBarry Smith return 0; 1947bf11e45SBarry Smith } 1957bf11e45SBarry Smith 1967bf11e45SBarry Smith /*@ 197564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1987bf11e45SBarry Smith 199fb4a63b6SLois Curfman McInnes Input Parameters: 2007bf11e45SBarry Smith . ts - timestep context 201564e8f4eSLois Curfman McInnes . update - latest solution vector 2027bf11e45SBarry Smith 203fb4a63b6SLois Curfman McInnes Output Parameters: 204fb4a63b6SLois Curfman McInnes . dt - newly computed timestep (if it had to shrink) 2057bf11e45SBarry Smith . flag - indicates if current timestep was ok 2067bf11e45SBarry Smith 207564e8f4eSLois Curfman McInnes Notes: 208564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 209564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 210564e8f4eSLois Curfman McInnes 211564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 212564e8f4eSLois Curfman McInnes 213564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 2147bf11e45SBarry Smith @*/ 2157bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 2167bf11e45SBarry Smith { 2177bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2187bf11e45SBarry Smith int ierr; 2197bf11e45SBarry Smith 2207bf11e45SBarry Smith if (!pseudo->verify) {*flag = 1; return 0;} 2217bf11e45SBarry Smith 2227bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 2237bf11e45SBarry Smith 2247bf11e45SBarry Smith return 0; 2257bf11e45SBarry Smith } 2267bf11e45SBarry Smith 2277bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 2287bf11e45SBarry Smith 2297bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 2302d3f70b5SBarry Smith { 2312d3f70b5SBarry Smith Vec sol = ts->vec_sol; 2327bf11e45SBarry Smith int ierr,i,max_steps = ts->max_steps,its,ok; 2337bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2347bf11e45SBarry Smith double current_time_step; 2352d3f70b5SBarry Smith 2362d3f70b5SBarry Smith *steps = -ts->steps; 2372d3f70b5SBarry Smith 2387bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 2392d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 2407bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 2417bf11e45SBarry Smith current_time_step = ts->time_step; 2427bf11e45SBarry Smith while (1) { 2437bf11e45SBarry Smith ts->ptime += current_time_step; 2447bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 2457bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 2467bf11e45SBarry Smith if (ok) break; 2477bf11e45SBarry Smith ts->ptime -= current_time_step; 2487bf11e45SBarry Smith current_time_step = ts->time_step; 2497bf11e45SBarry Smith } 2507bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 2512d3f70b5SBarry Smith ts->steps++; 2522d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 2532d3f70b5SBarry Smith } 2542d3f70b5SBarry Smith 2552d3f70b5SBarry Smith *steps += ts->steps; 2562d3f70b5SBarry Smith *time = ts->ptime; 2572d3f70b5SBarry Smith return 0; 2582d3f70b5SBarry Smith } 2592d3f70b5SBarry Smith 2602d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2617bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj ) 2622d3f70b5SBarry Smith { 2632d3f70b5SBarry Smith TS ts = (TS) obj; 2647bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2652d3f70b5SBarry Smith int ierr; 2662d3f70b5SBarry Smith 2677bf11e45SBarry Smith ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 2687bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 2697bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 2702d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 2717bf11e45SBarry Smith PetscFree(pseudo); 2722d3f70b5SBarry Smith return 0; 2732d3f70b5SBarry Smith } 2742d3f70b5SBarry Smith 2752d3f70b5SBarry Smith 2762d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2772d3f70b5SBarry Smith /* 2782d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 2792d3f70b5SBarry Smith */ 2802d3f70b5SBarry Smith 2817bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 2822d3f70b5SBarry Smith { 2832d3f70b5SBarry Smith TS ts; 2842d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 2852d3f70b5SBarry Smith int ierr; 2862d3f70b5SBarry Smith 2872d3f70b5SBarry Smith MatShellGetContext(mat,(void **)&ts); 2882d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2892d3f70b5SBarry Smith 2902d3f70b5SBarry Smith /* apply user provided function */ 2912d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 2922d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2932d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 2942d3f70b5SBarry Smith return 0; 2952d3f70b5SBarry Smith } 2962d3f70b5SBarry Smith 2972d3f70b5SBarry Smith /* 2982d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2992d3f70b5SBarry Smith 3002d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 3012d3f70b5SBarry Smith */ 3027bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 3032d3f70b5SBarry Smith { 3042d3f70b5SBarry Smith TS ts = (TS) ctx; 3052d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 3062d3f70b5SBarry Smith int ierr,i,n; 3072d3f70b5SBarry Smith 3082d3f70b5SBarry Smith /* apply user provided function */ 3092d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 3107bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 3112d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 3122d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 3132d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 3142d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 3152d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 3162d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 3172d3f70b5SBarry Smith } 3182d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 3192d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 3202d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 3212d3f70b5SBarry Smith 3222d3f70b5SBarry Smith return 0; 3232d3f70b5SBarry Smith } 3242d3f70b5SBarry Smith 3252d3f70b5SBarry Smith /* 3262d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 3272d3f70b5SBarry Smith 3282d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 3292d3f70b5SBarry Smith */ 3307bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 3312d3f70b5SBarry Smith { 3322d3f70b5SBarry Smith TS ts = (TS) ctx; 3332d3f70b5SBarry Smith int ierr; 3342d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 3352d3f70b5SBarry Smith MatType mtype; 3362d3f70b5SBarry Smith 3372d3f70b5SBarry Smith /* construct users Jacobian */ 3382d3f70b5SBarry Smith if (ts->rhsjacobian) { 3392d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 3402d3f70b5SBarry Smith } 3412d3f70b5SBarry Smith 3422d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 3432d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 3442d3f70b5SBarry Smith if (mtype != MATSHELL) { 3452d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 3462d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 3472d3f70b5SBarry Smith } 3482d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 3492d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 3502d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 3512d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 3522d3f70b5SBarry Smith } 3532d3f70b5SBarry Smith 3542d3f70b5SBarry Smith return 0; 3552d3f70b5SBarry Smith } 3562d3f70b5SBarry Smith 3572d3f70b5SBarry Smith 3587bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 3592d3f70b5SBarry Smith { 3607bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3612d3f70b5SBarry Smith int ierr, M, m; 3622d3f70b5SBarry Smith 3637bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 3647bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 3657bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 3662d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 3672d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 3682d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 3692d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 3707bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 3712d3f70b5SBarry Smith } 3727bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 3732d3f70b5SBarry Smith return 0; 3742d3f70b5SBarry Smith } 3752d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3762d3f70b5SBarry Smith 3772d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 3782d3f70b5SBarry Smith { 3797bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3802d3f70b5SBarry Smith 3817bf11e45SBarry Smith PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 3822d3f70b5SBarry Smith return 0; 3832d3f70b5SBarry Smith } 3842d3f70b5SBarry Smith 3857bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3862d3f70b5SBarry Smith { 3872d3f70b5SBarry Smith int ierr,flg; 38828aa8177SBarry Smith double inc; 3892d3f70b5SBarry Smith 3902d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 3912d3f70b5SBarry Smith 3922d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 3932d3f70b5SBarry Smith if (flg) { 39428aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 39528aa8177SBarry Smith } 39628aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 39728aa8177SBarry Smith if (flg) { 39828aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 3992d3f70b5SBarry Smith } 4002d3f70b5SBarry Smith return 0; 4012d3f70b5SBarry Smith } 4022d3f70b5SBarry Smith 4037bf11e45SBarry Smith static int TSPrintHelp_Pseudo(TS ts) 4042d3f70b5SBarry Smith { 4052d3f70b5SBarry Smith return 0; 4062d3f70b5SBarry Smith } 4072d3f70b5SBarry Smith 4087bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer) 4092d3f70b5SBarry Smith { 4102d3f70b5SBarry Smith return 0; 4112d3f70b5SBarry Smith } 4122d3f70b5SBarry Smith 4132d3f70b5SBarry Smith /* ------------------------------------------------------------ */ 4147bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 4152d3f70b5SBarry Smith { 4167bf11e45SBarry Smith TS_Pseudo *pseudo; 4172d3f70b5SBarry Smith int ierr; 4182d3f70b5SBarry Smith MatType mtype; 4192d3f70b5SBarry Smith 4207bf11e45SBarry Smith ts->type = TS_PSEUDO; 4217bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 4227bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 4237bf11e45SBarry Smith ts->view = TSView_Pseudo; 4242d3f70b5SBarry Smith 4252d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 4267bf11e45SBarry Smith SETERRQ(1,"TSCreate_Pseudo:Only for nonlinear problems"); 4272d3f70b5SBarry Smith } 4282d3f70b5SBarry Smith if (!ts->A) { 4297bf11e45SBarry Smith SETERRQ(1,"TSCreate_Pseudo:Must set Jacobian"); 4302d3f70b5SBarry Smith } 4312d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 4322d3f70b5SBarry Smith if (mtype == MATSHELL) { 4332d3f70b5SBarry Smith ts->Ashell = ts->A; 4342d3f70b5SBarry Smith } 4357bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 4367bf11e45SBarry Smith ts->step = TSStep_Pseudo; 4377bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 4387bf11e45SBarry Smith 4397bf11e45SBarry Smith /* create the required nonlinear solver context */ 4402d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 4412d3f70b5SBarry Smith 4427bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 4437bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 4447bf11e45SBarry Smith ts->data = (void *) pseudo; 4452d3f70b5SBarry Smith 44628aa8177SBarry Smith pseudo->dt_increment = 1.1; 44728aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 4482d3f70b5SBarry Smith return 0; 4492d3f70b5SBarry Smith } 4502d3f70b5SBarry Smith 4512d3f70b5SBarry Smith 45228aa8177SBarry Smith /*@ 45328aa8177SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 45428aa8177SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 45528aa8177SBarry Smith 45628aa8177SBarry Smith Input Parameters: 45728aa8177SBarry Smith . ts - the timestep context 45828aa8177SBarry Smith . inc - the scaling factor >= 1.0 45928aa8177SBarry Smith 46028aa8177SBarry Smith Options Database Key: 46128aa8177SBarry Smith $ -ts_pseudo_increment <increment> 46228aa8177SBarry Smith 463564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, set, increment 464564e8f4eSLois Curfman McInnes 46528aa8177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 46628aa8177SBarry Smith @*/ 46728aa8177SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 46828aa8177SBarry Smith { 46928aa8177SBarry Smith TS_Pseudo *pseudo; 47028aa8177SBarry Smith 47128aa8177SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 47228aa8177SBarry Smith if (ts->type != TS_PSEUDO) return 0; 47328aa8177SBarry Smith 47428aa8177SBarry Smith pseudo = (TS_Pseudo*) ts->data; 47528aa8177SBarry Smith pseudo->dt_increment = inc; 47628aa8177SBarry Smith return 0; 47728aa8177SBarry Smith } 4782d3f70b5SBarry Smith 4792d3f70b5SBarry Smith 4802d3f70b5SBarry Smith 481