1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*606d414cSSatish Balay static char vcid[] = "$Id: posindep.c,v 1.35 1999/06/08 22:58:09 balay Exp balay $"; 32d3f70b5SBarry Smith #endif 42d3f70b5SBarry Smith /* 5fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 62d3f70b5SBarry Smith */ 72d3f70b5SBarry Smith #include "src/ts/tsimpl.h" /*I "ts.h" I*/ 82d3f70b5SBarry Smith 92d3f70b5SBarry Smith typedef struct { 102d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 112d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 122d3f70b5SBarry Smith Vec rhs; /* work vector for RHS; vec_sol/dt */ 132d3f70b5SBarry Smith 142d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 152d3f70b5SBarry Smith 167bf11e45SBarry Smith int (*dt)(TS,double*,void*); /* compute next timestep, and related context */ 172d3f70b5SBarry Smith void *dtctx; 187bf11e45SBarry Smith int (*verify)(TS,Vec,void*,double*,int*); /* verify previous timestep and related context */ 197bf11e45SBarry Smith void *verifyctx; 202d3f70b5SBarry Smith 217bf11e45SBarry Smith double initial_fnorm,fnorm; /* original and current norm of F(u) */ 22ca4b7087SBarry Smith double fnorm_previous; 2328aa8177SBarry Smith 2428aa8177SBarry Smith double dt_increment; /* scaling that dt is incremented each time-step */ 25ca4b7087SBarry Smith int increment_dt_from_initial_dt; 267bf11e45SBarry Smith } TS_Pseudo; 272d3f70b5SBarry Smith 282d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 292d3f70b5SBarry Smith 3058562591SBarry Smith #undef __FUNC__ 3158562591SBarry Smith #define __FUNC__ "TSPseudoComputeTimeStep" 327bf11e45SBarry Smith /*@ 337bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 34564e8f4eSLois Curfman McInnes pseudo-timestepping process. 352d3f70b5SBarry Smith 3615091d37SBarry Smith Collective on TS 3715091d37SBarry Smith 387bf11e45SBarry Smith Input Parameter: 397bf11e45SBarry Smith . ts - timestep context 407bf11e45SBarry Smith 417bf11e45SBarry Smith Output Parameter: 42fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 43fb4a63b6SLois Curfman McInnes 4415091d37SBarry Smith Level: advanced 45564e8f4eSLois Curfman McInnes 46564e8f4eSLois Curfman McInnes Notes: 47564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 48564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 49564e8f4eSLois Curfman McInnes 50fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 51564e8f4eSLois Curfman McInnes 52564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 537bf11e45SBarry Smith @*/ 547bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 557bf11e45SBarry Smith { 567bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 577bf11e45SBarry Smith int ierr; 587bf11e45SBarry Smith 593a40ed3dSBarry Smith PetscFunctionBegin; 607bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 617bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 627bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 633a40ed3dSBarry Smith PetscFunctionReturn(0); 647bf11e45SBarry Smith } 657bf11e45SBarry Smith 667bf11e45SBarry Smith 677bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 6858562591SBarry Smith #undef __FUNC__ 6958562591SBarry Smith #define __FUNC__ "TSPseudoDefaultVerifyTimeStep" 707bf11e45SBarry Smith /*@C 71639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 727bf11e45SBarry Smith 7315091d37SBarry Smith Collective on TS 7415091d37SBarry Smith 757bf11e45SBarry Smith Input Parameters: 7615091d37SBarry Smith + ts - the timestep context 777bf11e45SBarry Smith . dtctx - unused timestep context 7815091d37SBarry Smith - update - latest solution vector 797bf11e45SBarry Smith 80564e8f4eSLois Curfman McInnes Output Parameters: 8115091d37SBarry Smith + newdt - the timestep to use for the next step 8215091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 837bf11e45SBarry Smith 8415091d37SBarry Smith Level: advanced 85fee21e36SBarry Smith 86564e8f4eSLois Curfman McInnes Note: 87564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 88564e8f4eSLois Curfman McInnes timestep. 89564e8f4eSLois Curfman McInnes 90564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 91564e8f4eSLois Curfman McInnes 92564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 937bf11e45SBarry Smith @*/ 947bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag) 957bf11e45SBarry Smith { 963a40ed3dSBarry Smith PetscFunctionBegin; 977bf11e45SBarry Smith *flag = 1; 983a40ed3dSBarry Smith PetscFunctionReturn(0); 997bf11e45SBarry Smith } 1007bf11e45SBarry Smith 1017bf11e45SBarry Smith 10258562591SBarry Smith #undef __FUNC__ 10358562591SBarry Smith #define __FUNC__ "TSPseudoVerifyTimeStep" 1047bf11e45SBarry Smith /*@ 105564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1067bf11e45SBarry Smith 10715091d37SBarry Smith Collective on TS 10815091d37SBarry Smith 109fb4a63b6SLois Curfman McInnes Input Parameters: 11015091d37SBarry Smith + ts - timestep context 11115091d37SBarry Smith - update - latest solution vector 1127bf11e45SBarry Smith 113fb4a63b6SLois Curfman McInnes Output Parameters: 11415091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11515091d37SBarry Smith - flag - indicates if current timestep was ok 1167bf11e45SBarry Smith 11715091d37SBarry Smith Level: advanced 118fee21e36SBarry Smith 119564e8f4eSLois Curfman McInnes Notes: 120564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 121564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 122564e8f4eSLois Curfman McInnes 123564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 124564e8f4eSLois Curfman McInnes 125564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1267bf11e45SBarry Smith @*/ 1277bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 1287bf11e45SBarry Smith { 1297bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1307bf11e45SBarry Smith int ierr; 1317bf11e45SBarry Smith 1323a40ed3dSBarry Smith PetscFunctionBegin; 1333a40ed3dSBarry Smith if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);} 1347bf11e45SBarry Smith 1357bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag );CHKERRQ(ierr); 1367bf11e45SBarry Smith 1373a40ed3dSBarry Smith PetscFunctionReturn(0); 1387bf11e45SBarry Smith } 1397bf11e45SBarry Smith 1407bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1417bf11e45SBarry Smith 14258562591SBarry Smith #undef __FUNC__ 14358562591SBarry Smith #define __FUNC__ "TSStep_Pseudo" 1447bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 1452d3f70b5SBarry Smith { 1462d3f70b5SBarry Smith Vec sol = ts->vec_sol; 147f2267985SLois Curfman McInnes int ierr,i,max_steps = ts->max_steps,its,ok,lits; 1487bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1497bf11e45SBarry Smith double current_time_step; 1502d3f70b5SBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 1522d3f70b5SBarry Smith *steps = -ts->steps; 1532d3f70b5SBarry Smith 1547bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr); 1552d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 1567bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 1577bf11e45SBarry Smith current_time_step = ts->time_step; 1587bf11e45SBarry Smith while (1) { 1597bf11e45SBarry Smith ts->ptime += current_time_step; 1607bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its);CHKERRQ(ierr); 161f2267985SLois Curfman McInnes ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr); 162f2267985SLois Curfman McInnes ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits; 1637bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr); 1647bf11e45SBarry Smith if (ok) break; 1657bf11e45SBarry Smith ts->ptime -= current_time_step; 1667bf11e45SBarry Smith current_time_step = ts->time_step; 1677bf11e45SBarry Smith } 1687bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr); 1692d3f70b5SBarry Smith ts->steps++; 1702d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1712d3f70b5SBarry Smith } 1722d3f70b5SBarry Smith 1732d3f70b5SBarry Smith *steps += ts->steps; 1742d3f70b5SBarry Smith *time = ts->ptime; 1753a40ed3dSBarry Smith PetscFunctionReturn(0); 1762d3f70b5SBarry Smith } 1772d3f70b5SBarry Smith 1782d3f70b5SBarry Smith /*------------------------------------------------------------*/ 17958562591SBarry Smith #undef __FUNC__ 180d4bb536fSBarry Smith #define __FUNC__ "TSDestroy_Pseudo" 181e1311b90SBarry Smith static int TSDestroy_Pseudo(TS ts ) 1822d3f70b5SBarry Smith { 1837bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 1842d3f70b5SBarry Smith int ierr; 1852d3f70b5SBarry Smith 1863a40ed3dSBarry Smith PetscFunctionBegin; 187e4ef1970SSatish Balay if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 1887bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1897bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 1902d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 191*606d414cSSatish Balay ierr = PetscFree(pseudo);CHKERRQ(ierr); 1923a40ed3dSBarry Smith PetscFunctionReturn(0); 1932d3f70b5SBarry Smith } 1942d3f70b5SBarry Smith 1952d3f70b5SBarry Smith 1962d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1972d3f70b5SBarry Smith /* 1982d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 1992d3f70b5SBarry Smith */ 2002d3f70b5SBarry Smith 20158562591SBarry Smith #undef __FUNC__ 20258562591SBarry Smith #define __FUNC__ "TSPseudoMatMult" 2037bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 2042d3f70b5SBarry Smith { 2052d3f70b5SBarry Smith TS ts; 2062d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 2072d3f70b5SBarry Smith int ierr; 2082d3f70b5SBarry Smith 2093a40ed3dSBarry Smith PetscFunctionBegin; 2103a40ed3dSBarry Smith ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr); 2112d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2122d3f70b5SBarry Smith 2132d3f70b5SBarry Smith /* apply user provided function */ 2142d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y);CHKERRQ(ierr); 2152d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2162d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y);CHKERRQ(ierr); 2173a40ed3dSBarry Smith PetscFunctionReturn(0); 2182d3f70b5SBarry Smith } 2192d3f70b5SBarry Smith 2202d3f70b5SBarry Smith /* 2212d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2222d3f70b5SBarry Smith 2232d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2242d3f70b5SBarry Smith */ 22558562591SBarry Smith #undef __FUNC__ 22658562591SBarry Smith #define __FUNC__ "TSPseudoFunction" 2277bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2282d3f70b5SBarry Smith { 2292d3f70b5SBarry Smith TS ts = (TS) ctx; 2302d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2312d3f70b5SBarry Smith int ierr,i,n; 2322d3f70b5SBarry Smith 2333a40ed3dSBarry Smith PetscFunctionBegin; 2342d3f70b5SBarry Smith /* apply user provided function */ 2352d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr); 2367bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2372d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr); 2382d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1);CHKERRQ(ierr); 2392d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr); 2402d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 2412d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 2422d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2432d3f70b5SBarry Smith } 244888f2ed8SSatish Balay ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr); 245888f2ed8SSatish Balay ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr); 246888f2ed8SSatish Balay ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr); 2472d3f70b5SBarry Smith 2483a40ed3dSBarry Smith PetscFunctionReturn(0); 2492d3f70b5SBarry Smith } 2502d3f70b5SBarry Smith 2512d3f70b5SBarry Smith /* 2522d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2532d3f70b5SBarry Smith 2542d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2552d3f70b5SBarry Smith */ 25658562591SBarry Smith #undef __FUNC__ 25758562591SBarry Smith #define __FUNC__ "TSPseudoJacobian" 2587bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2592d3f70b5SBarry Smith { 2602d3f70b5SBarry Smith TS ts = (TS) ctx; 2612d3f70b5SBarry Smith int ierr; 2622d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 2632d3f70b5SBarry Smith MatType mtype; 2642d3f70b5SBarry Smith 2653a40ed3dSBarry Smith PetscFunctionBegin; 2662d3f70b5SBarry Smith /* construct users Jacobian */ 2672d3f70b5SBarry Smith if (ts->rhsjacobian) { 2682d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 2692d3f70b5SBarry Smith } 2702d3f70b5SBarry Smith 2712d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 272888f2ed8SSatish Balay ierr = MatGetType(*AA,&mtype,PETSC_NULL);CHKERRQ(ierr); 2732d3f70b5SBarry Smith if (mtype != MATSHELL) { 2742d3f70b5SBarry Smith ierr = MatScale(&mone,*AA);CHKERRQ(ierr); 2752d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA);CHKERRQ(ierr); 2762d3f70b5SBarry Smith } 277888f2ed8SSatish Balay ierr = MatGetType(*BB,&mtype,PETSC_NULL);CHKERRQ(ierr); 2782d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 2792d3f70b5SBarry Smith ierr = MatScale(&mone,*BB);CHKERRQ(ierr); 2802d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB);CHKERRQ(ierr); 2812d3f70b5SBarry Smith } 2822d3f70b5SBarry Smith 2833a40ed3dSBarry Smith PetscFunctionReturn(0); 2842d3f70b5SBarry Smith } 2852d3f70b5SBarry Smith 2862d3f70b5SBarry Smith 28758562591SBarry Smith #undef __FUNC__ 28858562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo" 2897bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 2902d3f70b5SBarry Smith { 2917bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2922d3f70b5SBarry Smith int ierr, M, m; 2932d3f70b5SBarry Smith 2943a40ed3dSBarry Smith PetscFunctionBegin; 2957bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 2967bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 2977bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 2982d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 2992d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M);CHKERRQ(ierr); 3002d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m);CHKERRQ(ierr); 3012d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A);CHKERRQ(ierr); 3027bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 3032d3f70b5SBarry Smith } 3047bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 3053a40ed3dSBarry Smith PetscFunctionReturn(0); 3062d3f70b5SBarry Smith } 3072d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3082d3f70b5SBarry Smith 30958562591SBarry Smith #undef __FUNC__ 310d4bb536fSBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor" 3112d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 3122d3f70b5SBarry Smith { 3137bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 314d132466eSBarry Smith int ierr; 3152d3f70b5SBarry Smith 3163a40ed3dSBarry Smith PetscFunctionBegin; 317d132466eSBarry Smith ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm);CHKERRQ(ierr); 3183a40ed3dSBarry Smith PetscFunctionReturn(0); 3192d3f70b5SBarry Smith } 3202d3f70b5SBarry Smith 32158562591SBarry Smith #undef __FUNC__ 32258562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo" 3237bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3242d3f70b5SBarry Smith { 3252d3f70b5SBarry Smith int ierr,flg; 32628aa8177SBarry Smith double inc; 3272d3f70b5SBarry Smith 3283a40ed3dSBarry Smith PetscFunctionBegin; 3292d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 3302d3f70b5SBarry Smith 3312d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg);CHKERRQ(ierr); 3322d3f70b5SBarry Smith if (flg) { 33328aa8177SBarry Smith ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0);CHKERRQ(ierr); 33428aa8177SBarry Smith } 33528aa8177SBarry Smith ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg);CHKERRQ(ierr); 33628aa8177SBarry Smith if (flg) { 33728aa8177SBarry Smith ierr = TSPseudoSetTimeStepIncrement(ts,inc);CHKERRQ(ierr); 3382d3f70b5SBarry Smith } 339ca4b7087SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr); 340ca4b7087SBarry Smith if (flg) { 341ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 342ca4b7087SBarry Smith } 3433a40ed3dSBarry Smith PetscFunctionReturn(0); 3442d3f70b5SBarry Smith } 3452d3f70b5SBarry Smith 34658562591SBarry Smith #undef __FUNC__ 347d4bb536fSBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo" 348ca4b7087SBarry Smith static int TSPrintHelp_Pseudo(TS ts,char *p) 3492d3f70b5SBarry Smith { 350d132466eSBarry Smith int ierr; 351d132466eSBarry Smith 3523a40ed3dSBarry Smith PetscFunctionBegin; 353d132466eSBarry Smith ierr = (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n");CHKERRQ(ierr); 354d132466eSBarry Smith ierr = (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p);CHKERRQ(ierr); 355d132466eSBarry Smith ierr = (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p);CHKERRQ(ierr); 356d132466eSBarry Smith ierr = (*PetscHelpPrintf)(ts->comm," initial fnorm/current fnorm to determine new timestep\n");CHKERRQ(ierr); 357d132466eSBarry Smith ierr = (*PetscHelpPrintf)(ts->comm," default is current_dt * previous fnorm/current fnorm\n");CHKERRQ(ierr); 3583a40ed3dSBarry Smith PetscFunctionReturn(0); 3592d3f70b5SBarry Smith } 3602d3f70b5SBarry Smith 36158562591SBarry Smith #undef __FUNC__ 362d4bb536fSBarry Smith #define __FUNC__ "TSView_Pseudo" 363e1311b90SBarry Smith static int TSView_Pseudo(TS ts,Viewer viewer) 3642d3f70b5SBarry Smith { 3653a40ed3dSBarry Smith PetscFunctionBegin; 3663a40ed3dSBarry Smith PetscFunctionReturn(0); 3672d3f70b5SBarry Smith } 3682d3f70b5SBarry Smith 36982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 37082bf6240SBarry Smith #undef __FUNC__ 37182bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep" 37282bf6240SBarry Smith /*@ 37382bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 37482bf6240SBarry Smith last timestep. 37582bf6240SBarry Smith 37615091d37SBarry Smith Collective on TS 37715091d37SBarry Smith 37882bf6240SBarry Smith Input Parameters: 37915091d37SBarry Smith + ts - timestep context 38082bf6240SBarry Smith . dt - user-defined function to verify timestep 38115091d37SBarry Smith - ctx - [optional] user-defined context for private data 38282bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 38382bf6240SBarry Smith 38415091d37SBarry Smith Level: advanced 385fee21e36SBarry Smith 38682bf6240SBarry Smith Calling sequence of func: 38782bf6240SBarry Smith . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 38882bf6240SBarry Smith 38982bf6240SBarry Smith . update - latest solution vector 39082bf6240SBarry Smith . ctx - [optional] timestep context 39182bf6240SBarry Smith . newdt - the timestep to use for the next step 39282bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 39382bf6240SBarry Smith 39482bf6240SBarry Smith Notes: 39582bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 39682bf6240SBarry Smith during the timestepping process. 39782bf6240SBarry Smith 39882bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 39982bf6240SBarry Smith 40082bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 40182bf6240SBarry Smith @*/ 40282bf6240SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 40382bf6240SBarry Smith { 40482bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,Vec,void*,double *,int *),void *); 40582bf6240SBarry Smith 40682bf6240SBarry Smith PetscFunctionBegin; 40782bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 40882bf6240SBarry Smith 40994d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void **)&f);CHKERRQ(ierr); 41082bf6240SBarry Smith if (f) { 41182bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 41282bf6240SBarry Smith } 41382bf6240SBarry Smith PetscFunctionReturn(0); 41482bf6240SBarry Smith } 41582bf6240SBarry Smith 41682bf6240SBarry Smith #undef __FUNC__ 41782bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement" 41882bf6240SBarry Smith /*@ 41982bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 42082bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 42182bf6240SBarry Smith 422fee21e36SBarry Smith Collective on TS 423fee21e36SBarry Smith 42415091d37SBarry Smith Input Parameters: 42515091d37SBarry Smith + ts - the timestep context 42615091d37SBarry Smith - inc - the scaling factor >= 1.0 42715091d37SBarry Smith 42882bf6240SBarry Smith Options Database Key: 42982bf6240SBarry Smith $ -ts_pseudo_increment <increment> 43082bf6240SBarry Smith 43115091d37SBarry Smith Level: advanced 43215091d37SBarry Smith 43382bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 43482bf6240SBarry Smith 43582bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 43682bf6240SBarry Smith @*/ 43782bf6240SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc) 43882bf6240SBarry Smith { 43982bf6240SBarry Smith int ierr, (*f)(TS,double); 44082bf6240SBarry Smith 44182bf6240SBarry Smith PetscFunctionBegin; 44282bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 44382bf6240SBarry Smith 44494d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void **)&f);CHKERRQ(ierr); 44582bf6240SBarry Smith if (f) { 44682bf6240SBarry Smith ierr = (*f)(ts,inc);CHKERRQ(ierr); 44782bf6240SBarry Smith } 44882bf6240SBarry Smith PetscFunctionReturn(0); 44982bf6240SBarry Smith } 45082bf6240SBarry Smith 45182bf6240SBarry Smith #undef __FUNC__ 45282bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt" 45382bf6240SBarry Smith /*@ 45482bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 45582bf6240SBarry Smith is computed via the formula 45682bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 45782bf6240SBarry Smith rather than the default update, 45882bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 45982bf6240SBarry Smith 46015091d37SBarry Smith Collective on TS 46115091d37SBarry Smith 46282bf6240SBarry Smith Input Parameter: 46382bf6240SBarry Smith . ts - the timestep context 46482bf6240SBarry Smith 46582bf6240SBarry Smith Options Database Key: 46682bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 46782bf6240SBarry Smith 46815091d37SBarry Smith Level: advanced 46915091d37SBarry Smith 47082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 47182bf6240SBarry Smith 47282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 47382bf6240SBarry Smith @*/ 47482bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts) 47582bf6240SBarry Smith { 47682bf6240SBarry Smith int ierr, (*f)(TS); 47782bf6240SBarry Smith 47882bf6240SBarry Smith PetscFunctionBegin; 47982bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 48082bf6240SBarry Smith 48194d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void **)&f);CHKERRQ(ierr); 48282bf6240SBarry Smith if (f) { 48382bf6240SBarry Smith ierr = (*f)(ts);CHKERRQ(ierr); 48482bf6240SBarry Smith } 48582bf6240SBarry Smith PetscFunctionReturn(0); 48682bf6240SBarry Smith } 48782bf6240SBarry Smith 48882bf6240SBarry Smith 48982bf6240SBarry Smith #undef __FUNC__ 49082bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep" 49182bf6240SBarry Smith /*@ 49282bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 49382bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 49482bf6240SBarry Smith 49515091d37SBarry Smith Collective on TS 49615091d37SBarry Smith 49782bf6240SBarry Smith Input Parameters: 49815091d37SBarry Smith + ts - timestep context 49982bf6240SBarry Smith . dt - function to compute timestep 50015091d37SBarry Smith - ctx - [optional] user-defined context for private data 50182bf6240SBarry Smith required by the function (may be PETSC_NULL) 50282bf6240SBarry Smith 50315091d37SBarry Smith Level: intermediate 504fee21e36SBarry Smith 50582bf6240SBarry Smith Calling sequence of func: 50682bf6240SBarry Smith . func (TS ts,double *newdt,void *ctx); 50782bf6240SBarry Smith 50882bf6240SBarry Smith . newdt - the newly computed timestep 50982bf6240SBarry Smith . ctx - [optional] timestep context 51082bf6240SBarry Smith 51182bf6240SBarry Smith Notes: 51282bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 51382bf6240SBarry Smith during the timestepping process. 51482bf6240SBarry Smith 51582bf6240SBarry Smith .keywords: timestep, pseudo, set 51682bf6240SBarry Smith 51782bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 51882bf6240SBarry Smith @*/ 51982bf6240SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 52082bf6240SBarry Smith { 52182bf6240SBarry Smith int ierr, (*f)(TS,int (*)(TS,double *,void *),void *); 52282bf6240SBarry Smith 52382bf6240SBarry Smith PetscFunctionBegin; 52482bf6240SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 52582bf6240SBarry Smith 52694d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void **)&f);CHKERRQ(ierr); 52782bf6240SBarry Smith if (f) { 52882bf6240SBarry Smith ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 52982bf6240SBarry Smith } 53082bf6240SBarry Smith PetscFunctionReturn(0); 53182bf6240SBarry Smith } 53282bf6240SBarry Smith 53382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 53482bf6240SBarry Smith 535fb2e594dSBarry Smith EXTERN_C_BEGIN 53682bf6240SBarry Smith #undef __FUNC__ 53782bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep_Pseudo" 53882bf6240SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 53982bf6240SBarry Smith { 54082bf6240SBarry Smith TS_Pseudo *pseudo; 54182bf6240SBarry Smith 54282bf6240SBarry Smith PetscFunctionBegin; 54382bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 54482bf6240SBarry Smith pseudo->verify = dt; 54582bf6240SBarry Smith pseudo->verifyctx = ctx; 54682bf6240SBarry Smith PetscFunctionReturn(0); 54782bf6240SBarry Smith } 548fb2e594dSBarry Smith EXTERN_C_END 54982bf6240SBarry Smith 550fb2e594dSBarry Smith EXTERN_C_BEGIN 55182bf6240SBarry Smith #undef __FUNC__ 55282bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement_Pseudo" 55382bf6240SBarry Smith int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,double inc) 55482bf6240SBarry Smith { 55582bf6240SBarry Smith TS_Pseudo *pseudo; 55682bf6240SBarry Smith 55782bf6240SBarry Smith PetscFunctionBegin; 55882bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 55982bf6240SBarry Smith pseudo->dt_increment = inc; 56082bf6240SBarry Smith PetscFunctionReturn(0); 56182bf6240SBarry Smith } 562fb2e594dSBarry Smith EXTERN_C_END 56382bf6240SBarry Smith 564fb2e594dSBarry Smith EXTERN_C_BEGIN 56582bf6240SBarry Smith #undef __FUNC__ 56682bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 56782bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 56882bf6240SBarry Smith { 56982bf6240SBarry Smith TS_Pseudo *pseudo; 57082bf6240SBarry Smith 57182bf6240SBarry Smith PetscFunctionBegin; 57282bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 57382bf6240SBarry Smith pseudo->increment_dt_from_initial_dt = 1; 57482bf6240SBarry Smith PetscFunctionReturn(0); 57582bf6240SBarry Smith } 576fb2e594dSBarry Smith EXTERN_C_END 57782bf6240SBarry Smith 578fb2e594dSBarry Smith EXTERN_C_BEGIN 57982bf6240SBarry Smith #undef __FUNC__ 58082bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep_Pseudo" 58182bf6240SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,double*,void*),void* ctx) 58282bf6240SBarry Smith { 58382bf6240SBarry Smith TS_Pseudo *pseudo; 58482bf6240SBarry Smith 58582bf6240SBarry Smith PetscFunctionBegin; 58682bf6240SBarry Smith pseudo = (TS_Pseudo*) ts->data; 58782bf6240SBarry Smith pseudo->dt = dt; 58882bf6240SBarry Smith pseudo->dtctx = ctx; 58982bf6240SBarry Smith PetscFunctionReturn(0); 59082bf6240SBarry Smith } 591fb2e594dSBarry Smith EXTERN_C_END 59282bf6240SBarry Smith 59382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 59482bf6240SBarry Smith 595fb2e594dSBarry Smith EXTERN_C_BEGIN 59658562591SBarry Smith #undef __FUNC__ 59758562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo" 5987bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 5992d3f70b5SBarry Smith { 6007bf11e45SBarry Smith TS_Pseudo *pseudo; 6012d3f70b5SBarry Smith int ierr; 6022d3f70b5SBarry Smith MatType mtype; 6032d3f70b5SBarry Smith 6043a40ed3dSBarry Smith PetscFunctionBegin; 6057bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 6067bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 6077bf11e45SBarry Smith ts->view = TSView_Pseudo; 6082d3f70b5SBarry Smith 6092d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 610a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems"); 6112d3f70b5SBarry Smith } 6122d3f70b5SBarry Smith if (!ts->A) { 613a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian"); 6142d3f70b5SBarry Smith } 615888f2ed8SSatish Balay ierr = MatGetType(ts->A,&mtype,PETSC_NULL);CHKERRQ(ierr); 6162d3f70b5SBarry Smith if (mtype == MATSHELL) { 6172d3f70b5SBarry Smith ts->Ashell = ts->A; 6182d3f70b5SBarry Smith } 6197bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 6207bf11e45SBarry Smith ts->step = TSStep_Pseudo; 6217bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 6227bf11e45SBarry Smith 6237bf11e45SBarry Smith /* create the required nonlinear solver context */ 6242d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 6252d3f70b5SBarry Smith 6267bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo);CHKPTRQ(pseudo); 627eed86810SBarry Smith PLogObjectMemory(ts,sizeof(TS_Pseudo)); 628eed86810SBarry Smith 629549d3d68SSatish Balay ierr = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr); 6307bf11e45SBarry Smith ts->data = (void *) pseudo; 6312d3f70b5SBarry Smith 63228aa8177SBarry Smith pseudo->dt_increment = 1.1; 633ca4b7087SBarry Smith pseudo->increment_dt_from_initial_dt = 0; 63428aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 63582bf6240SBarry Smith 63694d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 637e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 638e1311b90SBarry Smith (void*)TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 63994d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 640e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 641e1311b90SBarry Smith (void*)TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 64294d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 643e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 644e1311b90SBarry Smith (void*)TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 64594d884c6SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C","TSPseudoSetTimeStep_Pseudo", 646e1311b90SBarry Smith (void*)TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6473a40ed3dSBarry Smith PetscFunctionReturn(0); 6482d3f70b5SBarry Smith } 649fb2e594dSBarry Smith EXTERN_C_END 6502d3f70b5SBarry Smith 65158562591SBarry Smith #undef __FUNC__ 65282bf6240SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep" 65382bf6240SBarry Smith /*@C 65482bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 65582bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 65628aa8177SBarry Smith 65715091d37SBarry Smith Collective on TS 65815091d37SBarry Smith 65928aa8177SBarry Smith Input Parameters: 66028aa8177SBarry Smith . ts - the timestep context 66182bf6240SBarry Smith . dtctx - unused timestep context 66228aa8177SBarry Smith 66382bf6240SBarry Smith Output Parameter: 66482bf6240SBarry Smith . newdt - the timestep to use for the next step 66528aa8177SBarry Smith 66615091d37SBarry Smith Level: advanced 66715091d37SBarry Smith 66882bf6240SBarry Smith .keywords: timestep, pseudo, default 669564e8f4eSLois Curfman McInnes 67082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 67128aa8177SBarry Smith @*/ 67282bf6240SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 67328aa8177SBarry Smith { 67482bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 67582bf6240SBarry Smith double inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 67682bf6240SBarry Smith int ierr; 67728aa8177SBarry Smith 6783a40ed3dSBarry Smith PetscFunctionBegin; 67982bf6240SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 68082bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 68182bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 68282bf6240SBarry Smith /* first time through so compute initial function norm */ 68382bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 68482bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 68582bf6240SBarry Smith } 68682bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 68782bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 68882bf6240SBarry Smith } 68982bf6240SBarry Smith else if (pseudo->increment_dt_from_initial_dt) { 69082bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 69182bf6240SBarry Smith } else { 69282bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 69382bf6240SBarry Smith } 69482bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6953a40ed3dSBarry Smith PetscFunctionReturn(0); 69628aa8177SBarry Smith } 6972d3f70b5SBarry Smith 698