12d3f70b5SBarry Smith #ifndef lint 22d3f70b5SBarry Smith static char vcid[] = "$Id: posindep.c,v 1.5 1996/09/14 12:37:20 bsmith Exp bsmith $"; 32d3f70b5SBarry Smith #endif 42d3f70b5SBarry Smith /* 52d3f70b5SBarry Smith Code for Time Stepping 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 19*7bf11e45SBarry Smith int (*dt)(TS,double*,void*); /* compute next timestep, and related context */ 202d3f70b5SBarry Smith void *dtctx; 21*7bf11e45SBarry Smith int (*verify)(TS,Vec,void*,double*,int*); /* verify previous timestep and related context */ 22*7bf11e45SBarry Smith void *verifyctx; 232d3f70b5SBarry Smith 24*7bf11e45SBarry Smith double initial_fnorm,fnorm; /* original and current norm of F(u) */ 25*7bf11e45SBarry Smith } TS_Pseudo; 262d3f70b5SBarry Smith 272d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 282d3f70b5SBarry Smith /*@C 29*7bf11e45SBarry Smith TSPseudoDefaultTimeStep - Default code to compute 30*7bf11e45SBarry Smith pseudo timestepping. Use with TSPseudoSetTimeStep(). 312d3f70b5SBarry Smith 322d3f70b5SBarry Smith Input Parameters: 332d3f70b5SBarry Smith . ts - the time step context 34*7bf11e45SBarry Smith . dtctx - unused timestep context 352d3f70b5SBarry Smith 362d3f70b5SBarry Smith Output Parameter: 372d3f70b5SBarry Smith . newdt - the time step to use for the next step 382d3f70b5SBarry Smith 392d3f70b5SBarry Smith @*/ 40*7bf11e45SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 412d3f70b5SBarry Smith { 42*7bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 432d3f70b5SBarry Smith int ierr; 442d3f70b5SBarry Smith 45*7bf11e45SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 46*7bf11e45SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 47*7bf11e45SBarry Smith if (pseudo->initial_fnorm == 0.0) { 482d3f70b5SBarry Smith /* first time through so compute initial function norm */ 49*7bf11e45SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 502d3f70b5SBarry Smith } 51*7bf11e45SBarry Smith if (pseudo->fnorm == 0.0) *newdt = 1.e12*1.1*ts->time_step; 52*7bf11e45SBarry Smith else *newdt = 1.1*ts->time_step*pseudo->initial_fnorm/pseudo->fnorm; 532d3f70b5SBarry Smith return 0; 542d3f70b5SBarry Smith } 552d3f70b5SBarry Smith 562d3f70b5SBarry Smith /*@ 57*7bf11e45SBarry Smith TSPseudoSetTimeStep - Sets the user routine to be 582d3f70b5SBarry Smith called at each pseudo-time-step to update the time-step. 592d3f70b5SBarry Smith 602d3f70b5SBarry Smith Input Parameters: 612d3f70b5SBarry Smith . ts - time step context 622d3f70b5SBarry Smith . dt - function to compute timestep 632d3f70b5SBarry Smith . ctx - [optional] context required by function 642d3f70b5SBarry Smith 652d3f70b5SBarry Smith @*/ 66*7bf11e45SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 672d3f70b5SBarry Smith { 68*7bf11e45SBarry Smith TS_Pseudo *pseudo; 692d3f70b5SBarry Smith 702d3f70b5SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 71*7bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 72*7bf11e45SBarry Smith 73*7bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 74*7bf11e45SBarry Smith pseudo->dt = dt; 75*7bf11e45SBarry Smith pseudo->dtctx = ctx; 762d3f70b5SBarry Smith return 0; 772d3f70b5SBarry Smith } 782d3f70b5SBarry Smith 79*7bf11e45SBarry Smith /*@ 80*7bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 81*7bf11e45SBarry Smith pseudo-timestepping. 822d3f70b5SBarry Smith 83*7bf11e45SBarry Smith Input Parameter: 84*7bf11e45SBarry Smith . ts - time step context 85*7bf11e45SBarry Smith 86*7bf11e45SBarry Smith Output Parameter: 87*7bf11e45SBarry Smith . dt - newly computed time-step 88*7bf11e45SBarry Smith @*/ 89*7bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt) 90*7bf11e45SBarry Smith { 91*7bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 92*7bf11e45SBarry Smith int ierr; 93*7bf11e45SBarry Smith 94*7bf11e45SBarry Smith PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 95*7bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt, pseudo->dtctx); CHKERRQ(ierr); 96*7bf11e45SBarry Smith PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 97*7bf11e45SBarry Smith return 0; 98*7bf11e45SBarry Smith } 99*7bf11e45SBarry Smith 100*7bf11e45SBarry Smith 101*7bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 102*7bf11e45SBarry Smith /*@C 103*7bf11e45SBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify last timestep. 104*7bf11e45SBarry Smith 105*7bf11e45SBarry Smith Input Parameters: 106*7bf11e45SBarry Smith . ts - the time step context 107*7bf11e45SBarry Smith . dtctx - unused timestep context 108*7bf11e45SBarry Smith 109*7bf11e45SBarry Smith Output Parameter: 110*7bf11e45SBarry Smith . newdt - the time step to use for the next step 111*7bf11e45SBarry Smith 112*7bf11e45SBarry Smith @*/ 113*7bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void* dtctx,double* newdt,int *flag) 114*7bf11e45SBarry Smith { 115*7bf11e45SBarry Smith *flag = 1; 116*7bf11e45SBarry Smith return 0; 117*7bf11e45SBarry Smith } 118*7bf11e45SBarry Smith 119*7bf11e45SBarry Smith /*@ 120*7bf11e45SBarry Smith TSPseudoSetVerifyTimeStep - Sets the user routine to verify quality of last time step. 121*7bf11e45SBarry Smith 122*7bf11e45SBarry Smith Input Parameters: 123*7bf11e45SBarry Smith . ts - time step context 124*7bf11e45SBarry Smith . dt - function to verify 125*7bf11e45SBarry Smith . ctx - [optional] context required by function 126*7bf11e45SBarry Smith 127*7bf11e45SBarry Smith @*/ 128*7bf11e45SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 129*7bf11e45SBarry Smith { 130*7bf11e45SBarry Smith TS_Pseudo *pseudo; 131*7bf11e45SBarry Smith 132*7bf11e45SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 133*7bf11e45SBarry Smith if (ts->type != TS_PSEUDO) return 0; 134*7bf11e45SBarry Smith 135*7bf11e45SBarry Smith pseudo = (TS_Pseudo*) ts->data; 136*7bf11e45SBarry Smith pseudo->verify = dt; 137*7bf11e45SBarry Smith pseudo->verifyctx = ctx; 138*7bf11e45SBarry Smith return 0; 139*7bf11e45SBarry Smith } 140*7bf11e45SBarry Smith 141*7bf11e45SBarry Smith /*@ 142*7bf11e45SBarry Smith TSPseudoVerifyTimeStep - Verifies that the last time step was ok. 143*7bf11e45SBarry Smith 144*7bf11e45SBarry Smith Input Parameter: 145*7bf11e45SBarry Smith . ts - time step context 146*7bf11e45SBarry Smith . update - latest solution 147*7bf11e45SBarry Smith 148*7bf11e45SBarry Smith Output Parameter: 149*7bf11e45SBarry Smith . dt - newly computed time-step (if it had to shrink) 150*7bf11e45SBarry Smith . flag - indicates if current timestep was ok 151*7bf11e45SBarry Smith 152*7bf11e45SBarry Smith @*/ 153*7bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 154*7bf11e45SBarry Smith { 155*7bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 156*7bf11e45SBarry Smith int ierr; 157*7bf11e45SBarry Smith 158*7bf11e45SBarry Smith if (!pseudo->verify) {*flag = 1; return 0;} 159*7bf11e45SBarry Smith 160*7bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 161*7bf11e45SBarry Smith 162*7bf11e45SBarry Smith return 0; 163*7bf11e45SBarry Smith } 164*7bf11e45SBarry Smith 165*7bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 166*7bf11e45SBarry Smith 167*7bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time) 1682d3f70b5SBarry Smith { 1692d3f70b5SBarry Smith Vec sol = ts->vec_sol; 170*7bf11e45SBarry Smith int ierr,i,max_steps = ts->max_steps,its,ok; 171*7bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 172*7bf11e45SBarry Smith double current_time_step; 1732d3f70b5SBarry Smith 1742d3f70b5SBarry Smith *steps = -ts->steps; 1752d3f70b5SBarry Smith 176*7bf11e45SBarry Smith ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 1772d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 178*7bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 179*7bf11e45SBarry Smith current_time_step = ts->time_step; 180*7bf11e45SBarry Smith while (1) { 181*7bf11e45SBarry Smith ts->ptime += current_time_step; 182*7bf11e45SBarry Smith ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 183*7bf11e45SBarry Smith ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 184*7bf11e45SBarry Smith if (ok) break; 185*7bf11e45SBarry Smith ts->ptime -= current_time_step; 186*7bf11e45SBarry Smith current_time_step = ts->time_step; 187*7bf11e45SBarry Smith } 188*7bf11e45SBarry Smith ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 1892d3f70b5SBarry Smith ts->steps++; 1902d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1912d3f70b5SBarry Smith } 1922d3f70b5SBarry Smith 1932d3f70b5SBarry Smith *steps += ts->steps; 1942d3f70b5SBarry Smith *time = ts->ptime; 1952d3f70b5SBarry Smith return 0; 1962d3f70b5SBarry Smith } 1972d3f70b5SBarry Smith 1982d3f70b5SBarry Smith /*------------------------------------------------------------*/ 199*7bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj ) 2002d3f70b5SBarry Smith { 2012d3f70b5SBarry Smith TS ts = (TS) obj; 202*7bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2032d3f70b5SBarry Smith int ierr; 2042d3f70b5SBarry Smith 205*7bf11e45SBarry Smith ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 206*7bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 207*7bf11e45SBarry Smith if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 2082d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 209*7bf11e45SBarry Smith PetscFree(pseudo); 2102d3f70b5SBarry Smith return 0; 2112d3f70b5SBarry Smith } 2122d3f70b5SBarry Smith 2132d3f70b5SBarry Smith 2142d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2152d3f70b5SBarry Smith /* 2162d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 2172d3f70b5SBarry Smith */ 2182d3f70b5SBarry Smith 219*7bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y) 2202d3f70b5SBarry Smith { 2212d3f70b5SBarry Smith TS ts; 2222d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 2232d3f70b5SBarry Smith int ierr; 2242d3f70b5SBarry Smith 2252d3f70b5SBarry Smith MatShellGetContext(mat,(void **)&ts); 2262d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 2272d3f70b5SBarry Smith 2282d3f70b5SBarry Smith /* apply user provided function */ 2292d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 2302d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 2312d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 2322d3f70b5SBarry Smith return 0; 2332d3f70b5SBarry Smith } 2342d3f70b5SBarry Smith 2352d3f70b5SBarry Smith /* 2362d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 2372d3f70b5SBarry Smith 2382d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 2392d3f70b5SBarry Smith */ 240*7bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 2412d3f70b5SBarry Smith { 2422d3f70b5SBarry Smith TS ts = (TS) ctx; 2432d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 2442d3f70b5SBarry Smith int ierr,i,n; 2452d3f70b5SBarry Smith 2462d3f70b5SBarry Smith /* apply user provided function */ 2472d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 248*7bf11e45SBarry Smith /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 2492d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 2502d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 2512d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 2522d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 2532d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 2542d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 2552d3f70b5SBarry Smith } 2562d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 2572d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 2582d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 2592d3f70b5SBarry Smith 2602d3f70b5SBarry Smith return 0; 2612d3f70b5SBarry Smith } 2622d3f70b5SBarry Smith 2632d3f70b5SBarry Smith /* 2642d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 2652d3f70b5SBarry Smith 2662d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 2672d3f70b5SBarry Smith */ 268*7bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 2692d3f70b5SBarry Smith { 2702d3f70b5SBarry Smith TS ts = (TS) ctx; 2712d3f70b5SBarry Smith int ierr; 2722d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 2732d3f70b5SBarry Smith MatType mtype; 2742d3f70b5SBarry Smith 2752d3f70b5SBarry Smith /* construct users Jacobian */ 2762d3f70b5SBarry Smith if (ts->rhsjacobian) { 2772d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 2782d3f70b5SBarry Smith } 2792d3f70b5SBarry Smith 2802d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 2812d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 2822d3f70b5SBarry Smith if (mtype != MATSHELL) { 2832d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 2842d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 2852d3f70b5SBarry Smith } 2862d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 2872d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 2882d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 2892d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 2902d3f70b5SBarry Smith } 2912d3f70b5SBarry Smith 2922d3f70b5SBarry Smith return 0; 2932d3f70b5SBarry Smith } 2942d3f70b5SBarry Smith 2952d3f70b5SBarry Smith 296*7bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts) 2972d3f70b5SBarry Smith { 298*7bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 2992d3f70b5SBarry Smith int ierr, M, m; 3002d3f70b5SBarry Smith 301*7bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 302*7bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 303*7bf11e45SBarry Smith ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 3042d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 3052d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 3062d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 3072d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 308*7bf11e45SBarry Smith ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 3092d3f70b5SBarry Smith } 310*7bf11e45SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 3112d3f70b5SBarry Smith return 0; 3122d3f70b5SBarry Smith } 3132d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3142d3f70b5SBarry Smith 3152d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 3162d3f70b5SBarry Smith { 317*7bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3182d3f70b5SBarry Smith 319*7bf11e45SBarry Smith PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 3202d3f70b5SBarry Smith return 0; 3212d3f70b5SBarry Smith } 3222d3f70b5SBarry Smith 323*7bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts) 3242d3f70b5SBarry Smith { 3252d3f70b5SBarry Smith int ierr,flg; 3262d3f70b5SBarry Smith 3272d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 3282d3f70b5SBarry Smith 3292d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 3302d3f70b5SBarry Smith if (flg) { 3312d3f70b5SBarry Smith TSSetMonitor(ts,TSPseudoDefaultMonitor,0); 3322d3f70b5SBarry Smith } 3332d3f70b5SBarry Smith return 0; 3342d3f70b5SBarry Smith } 3352d3f70b5SBarry Smith 336*7bf11e45SBarry Smith static int TSPrintHelp_Pseudo(TS ts) 3372d3f70b5SBarry Smith { 3382d3f70b5SBarry Smith 3392d3f70b5SBarry Smith return 0; 3402d3f70b5SBarry Smith } 3412d3f70b5SBarry Smith 342*7bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer) 3432d3f70b5SBarry Smith { 3442d3f70b5SBarry Smith return 0; 3452d3f70b5SBarry Smith } 3462d3f70b5SBarry Smith 3472d3f70b5SBarry Smith /* ------------------------------------------------------------ */ 348*7bf11e45SBarry Smith int TSCreate_Pseudo(TS ts ) 3492d3f70b5SBarry Smith { 350*7bf11e45SBarry Smith TS_Pseudo *pseudo; 3512d3f70b5SBarry Smith int ierr; 3522d3f70b5SBarry Smith MatType mtype; 3532d3f70b5SBarry Smith 354*7bf11e45SBarry Smith ts->type = TS_PSEUDO; 355*7bf11e45SBarry Smith ts->destroy = TSDestroy_Pseudo; 356*7bf11e45SBarry Smith ts->printhelp = TSPrintHelp_Pseudo; 357*7bf11e45SBarry Smith ts->view = TSView_Pseudo; 3582d3f70b5SBarry Smith 3592d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 360*7bf11e45SBarry Smith SETERRQ(1,"TSCreate_Pseudo:Only for nonlinear problems"); 3612d3f70b5SBarry Smith } 3622d3f70b5SBarry Smith if (!ts->A) { 363*7bf11e45SBarry Smith SETERRQ(1,"TSCreate_Pseudo:Must set Jacobian"); 3642d3f70b5SBarry Smith } 3652d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 3662d3f70b5SBarry Smith if (mtype == MATSHELL) { 3672d3f70b5SBarry Smith ts->Ashell = ts->A; 3682d3f70b5SBarry Smith } 369*7bf11e45SBarry Smith ts->setup = TSSetUp_Pseudo; 370*7bf11e45SBarry Smith ts->step = TSStep_Pseudo; 371*7bf11e45SBarry Smith ts->setfromoptions = TSSetFromOptions_Pseudo; 372*7bf11e45SBarry Smith 373*7bf11e45SBarry Smith /* create the required nonlinear solver context */ 3742d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 3752d3f70b5SBarry Smith 376*7bf11e45SBarry Smith pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 377*7bf11e45SBarry Smith PetscMemzero(pseudo,sizeof(TS_Pseudo)); 378*7bf11e45SBarry Smith ts->data = (void *) pseudo; 3792d3f70b5SBarry Smith 3802d3f70b5SBarry Smith return 0; 3812d3f70b5SBarry Smith } 3822d3f70b5SBarry Smith 3832d3f70b5SBarry Smith 3842d3f70b5SBarry Smith 3852d3f70b5SBarry Smith 3862d3f70b5SBarry Smith 387