xref: /petsc/src/ts/impls/pseudo/posindep.c (revision fb4a63b695def8baf7b53fd7b7ab95f1c9876fc4)
12d3f70b5SBarry Smith #ifndef lint
2*fb4a63b6SLois Curfman McInnes static char vcid[] = "$Id: posindep.c,v 1.2 1996/09/29 14:39:45 bsmith Exp curfman $";
32d3f70b5SBarry Smith #endif
42d3f70b5SBarry Smith /*
5*fb4a63b6SLois 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) */
257bf11e45SBarry Smith } TS_Pseudo;
262d3f70b5SBarry Smith 
272d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
282d3f70b5SBarry Smith /*@C
29*fb4a63b6SLois Curfman McInnes    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
30*fb4a63b6SLois Curfman McInnes    Use with TSPseudoSetTimeStep().
312d3f70b5SBarry Smith 
322d3f70b5SBarry Smith    Input Parameters:
332d3f70b5SBarry Smith .  ts - the timestep context
347bf11e45SBarry Smith .  dtctx - unused timestep context
352d3f70b5SBarry Smith 
362d3f70b5SBarry Smith    Output Parameter:
372d3f70b5SBarry Smith .  newdt - the timestep to use for the next step
382d3f70b5SBarry Smith 
39*fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, default
40*fb4a63b6SLois Curfman McInnes 
41*fb4a63b6SLois Curfman McInnes .seealso: TSPseudoSetTimeStep()
422d3f70b5SBarry Smith @*/
437bf11e45SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx)
442d3f70b5SBarry Smith {
457bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
462d3f70b5SBarry Smith   int       ierr;
472d3f70b5SBarry Smith 
487bf11e45SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
497bf11e45SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr);
507bf11e45SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
512d3f70b5SBarry Smith     /* first time through so compute initial function norm */
527bf11e45SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
532d3f70b5SBarry Smith   }
547bf11e45SBarry Smith   if (pseudo->fnorm == 0.0) *newdt = 1.e12*1.1*ts->time_step;
557bf11e45SBarry Smith   else                      *newdt = 1.1*ts->time_step*pseudo->initial_fnorm/pseudo->fnorm;
562d3f70b5SBarry Smith   return 0;
572d3f70b5SBarry Smith }
582d3f70b5SBarry Smith 
592d3f70b5SBarry Smith /*@
60*fb4a63b6SLois Curfman McInnes    TSPseudoSetTimeStep - Sets the user-defined routine to be
61*fb4a63b6SLois Curfman McInnes    called at each pseudo-timestep to update the timestep.
622d3f70b5SBarry Smith 
632d3f70b5SBarry Smith    Input Parameters:
642d3f70b5SBarry Smith .  ts - timestep context
652d3f70b5SBarry Smith .  dt - function to compute timestep
662d3f70b5SBarry Smith .  ctx - [optional] context required by function
672d3f70b5SBarry Smith 
68*fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, set
69*fb4a63b6SLois Curfman McInnes 
70*fb4a63b6SLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep()
712d3f70b5SBarry Smith @*/
727bf11e45SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx)
732d3f70b5SBarry Smith {
747bf11e45SBarry Smith   TS_Pseudo *pseudo;
752d3f70b5SBarry Smith 
762d3f70b5SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
777bf11e45SBarry Smith   if (ts->type != TS_PSEUDO) return 0;
787bf11e45SBarry Smith 
797bf11e45SBarry Smith   pseudo          = (TS_Pseudo*) ts->data;
807bf11e45SBarry Smith   pseudo->dt      = dt;
817bf11e45SBarry Smith   pseudo->dtctx   = ctx;
822d3f70b5SBarry Smith   return 0;
832d3f70b5SBarry Smith }
842d3f70b5SBarry Smith 
857bf11e45SBarry Smith /*@
867bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
877bf11e45SBarry Smith     pseudo-timestepping.
882d3f70b5SBarry Smith 
897bf11e45SBarry Smith     Input Parameter:
907bf11e45SBarry Smith .   ts - timestep context
917bf11e45SBarry Smith 
927bf11e45SBarry Smith     Output Parameter:
93*fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
94*fb4a63b6SLois Curfman McInnes 
95*fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
967bf11e45SBarry Smith @*/
977bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt)
987bf11e45SBarry Smith {
997bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
1007bf11e45SBarry Smith   int       ierr;
1017bf11e45SBarry Smith 
1027bf11e45SBarry Smith   PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
1037bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt, pseudo->dtctx); CHKERRQ(ierr);
1047bf11e45SBarry Smith   PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
1057bf11e45SBarry Smith   return 0;
1067bf11e45SBarry Smith }
1077bf11e45SBarry Smith 
1087bf11e45SBarry Smith 
1097bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
1107bf11e45SBarry Smith /*@C
1117bf11e45SBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify last timestep.
1127bf11e45SBarry Smith 
1137bf11e45SBarry Smith    Input Parameters:
1147bf11e45SBarry Smith .  ts - the timestep context
1157bf11e45SBarry Smith .  dtctx - unused timestep context
1167bf11e45SBarry Smith 
1177bf11e45SBarry Smith    Output Parameter:
1187bf11e45SBarry Smith .  newdt - the timestep to use for the next step
1197bf11e45SBarry Smith 
1207bf11e45SBarry Smith @*/
1217bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void* dtctx,double* newdt,int *flag)
1227bf11e45SBarry Smith {
1237bf11e45SBarry Smith   *flag = 1;
1247bf11e45SBarry Smith   return 0;
1257bf11e45SBarry Smith }
1267bf11e45SBarry Smith 
1277bf11e45SBarry Smith /*@
1287bf11e45SBarry Smith    TSPseudoSetVerifyTimeStep - Sets the user routine to verify quality of last timestep.
1297bf11e45SBarry Smith 
1307bf11e45SBarry Smith    Input Parameters:
1317bf11e45SBarry Smith .  ts - timestep context
1327bf11e45SBarry Smith .  dt - function to verify
1337bf11e45SBarry Smith .  ctx - [optional] context required by function
1347bf11e45SBarry Smith 
1357bf11e45SBarry Smith @*/
1367bf11e45SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx)
1377bf11e45SBarry Smith {
1387bf11e45SBarry Smith   TS_Pseudo *pseudo;
1397bf11e45SBarry Smith 
1407bf11e45SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1417bf11e45SBarry Smith   if (ts->type != TS_PSEUDO) return 0;
1427bf11e45SBarry Smith 
1437bf11e45SBarry Smith   pseudo              = (TS_Pseudo*) ts->data;
1447bf11e45SBarry Smith   pseudo->verify      = dt;
1457bf11e45SBarry Smith   pseudo->verifyctx   = ctx;
1467bf11e45SBarry Smith   return 0;
1477bf11e45SBarry Smith }
1487bf11e45SBarry Smith 
1497bf11e45SBarry Smith /*@
150*fb4a63b6SLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies that the last timestep was OK.
1517bf11e45SBarry Smith 
152*fb4a63b6SLois Curfman McInnes     Input Parameters:
1537bf11e45SBarry Smith .   ts - timestep context
1547bf11e45SBarry Smith .   update - latest solution
1557bf11e45SBarry Smith 
156*fb4a63b6SLois Curfman McInnes     Output Parameters:
157*fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep (if it had to shrink)
1587bf11e45SBarry Smith .   flag - indicates if current timestep was ok
1597bf11e45SBarry Smith 
1607bf11e45SBarry Smith @*/
1617bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag)
1627bf11e45SBarry Smith {
1637bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
1647bf11e45SBarry Smith   int       ierr;
1657bf11e45SBarry Smith 
1667bf11e45SBarry Smith   if (!pseudo->verify) {*flag = 1; return 0;}
1677bf11e45SBarry Smith 
1687bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr);
1697bf11e45SBarry Smith 
1707bf11e45SBarry Smith   return 0;
1717bf11e45SBarry Smith }
1727bf11e45SBarry Smith 
1737bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1747bf11e45SBarry Smith 
1757bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time)
1762d3f70b5SBarry Smith {
1772d3f70b5SBarry Smith   Vec       sol = ts->vec_sol;
1787bf11e45SBarry Smith   int       ierr,i,max_steps = ts->max_steps,its,ok;
1797bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
1807bf11e45SBarry Smith   double    current_time_step;
1812d3f70b5SBarry Smith 
1822d3f70b5SBarry Smith   *steps = -ts->steps;
1832d3f70b5SBarry Smith 
1847bf11e45SBarry Smith   ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr);
1852d3f70b5SBarry Smith   for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) {
1867bf11e45SBarry Smith     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr);
1877bf11e45SBarry Smith     current_time_step = ts->time_step;
1887bf11e45SBarry Smith     while (1) {
1897bf11e45SBarry Smith       ts->ptime  += current_time_step;
1907bf11e45SBarry Smith       ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr);
1917bf11e45SBarry Smith       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr);
1927bf11e45SBarry Smith       if (ok) break;
1937bf11e45SBarry Smith       ts->ptime        -= current_time_step;
1947bf11e45SBarry Smith       current_time_step = ts->time_step;
1957bf11e45SBarry Smith     }
1967bf11e45SBarry Smith     ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr);
1972d3f70b5SBarry Smith     ts->steps++;
1982d3f70b5SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1992d3f70b5SBarry Smith   }
2002d3f70b5SBarry Smith 
2012d3f70b5SBarry Smith   *steps += ts->steps;
2022d3f70b5SBarry Smith   *time  = ts->ptime;
2032d3f70b5SBarry Smith   return 0;
2042d3f70b5SBarry Smith }
2052d3f70b5SBarry Smith 
2062d3f70b5SBarry Smith /*------------------------------------------------------------*/
2077bf11e45SBarry Smith static int TSDestroy_Pseudo(PetscObject obj )
2082d3f70b5SBarry Smith {
2092d3f70b5SBarry Smith   TS        ts = (TS) obj;
2107bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
2112d3f70b5SBarry Smith   int       ierr;
2122d3f70b5SBarry Smith 
2137bf11e45SBarry Smith   ierr = VecDestroy(pseudo->update); CHKERRQ(ierr);
2147bf11e45SBarry Smith   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
2157bf11e45SBarry Smith   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
2162d3f70b5SBarry Smith   if (ts->Ashell)   {ierr = MatDestroy(ts->A); CHKERRQ(ierr);}
2177bf11e45SBarry Smith   PetscFree(pseudo);
2182d3f70b5SBarry Smith   return 0;
2192d3f70b5SBarry Smith }
2202d3f70b5SBarry Smith 
2212d3f70b5SBarry Smith 
2222d3f70b5SBarry Smith /*------------------------------------------------------------*/
2232d3f70b5SBarry Smith /*
2242d3f70b5SBarry Smith     This matrix shell multiply where user provided Shell matrix
2252d3f70b5SBarry Smith */
2262d3f70b5SBarry Smith 
2277bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y)
2282d3f70b5SBarry Smith {
2292d3f70b5SBarry Smith   TS     ts;
2302d3f70b5SBarry Smith   Scalar mdt,mone = -1.0;
2312d3f70b5SBarry Smith   int    ierr;
2322d3f70b5SBarry Smith 
2332d3f70b5SBarry Smith   MatShellGetContext(mat,(void **)&ts);
2342d3f70b5SBarry Smith   mdt = 1.0/ts->time_step;
2352d3f70b5SBarry Smith 
2362d3f70b5SBarry Smith   /* apply user provided function */
2372d3f70b5SBarry Smith   ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr);
2382d3f70b5SBarry Smith   /* shift and scale by 1/dt - F */
2392d3f70b5SBarry Smith   ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr);
2402d3f70b5SBarry Smith   return 0;
2412d3f70b5SBarry Smith }
2422d3f70b5SBarry Smith 
2432d3f70b5SBarry Smith /*
2442d3f70b5SBarry Smith     This defines the nonlinear equation that is to be solved with SNES
2452d3f70b5SBarry Smith 
2462d3f70b5SBarry Smith               (U^{n+1} - U^{n})/dt - F(U^{n+1})
2472d3f70b5SBarry Smith */
2487bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
2492d3f70b5SBarry Smith {
2502d3f70b5SBarry Smith   TS     ts = (TS) ctx;
2512d3f70b5SBarry Smith   Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
2522d3f70b5SBarry Smith   int    ierr,i,n;
2532d3f70b5SBarry Smith 
2542d3f70b5SBarry Smith   /* apply user provided function */
2552d3f70b5SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr);
2567bf11e45SBarry Smith   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
2572d3f70b5SBarry Smith   ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr);
2582d3f70b5SBarry Smith   ierr = VecGetArray(x,&unp1); CHKERRQ(ierr);
2592d3f70b5SBarry Smith   ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr);
2602d3f70b5SBarry Smith   ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
2612d3f70b5SBarry Smith   for ( i=0; i<n; i++ ) {
2622d3f70b5SBarry Smith     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
2632d3f70b5SBarry Smith   }
2642d3f70b5SBarry Smith   ierr = VecRestoreArray(ts->vec_sol,&un);
2652d3f70b5SBarry Smith   ierr = VecRestoreArray(x,&unp1);
2662d3f70b5SBarry Smith   ierr = VecRestoreArray(y,&Funp1);
2672d3f70b5SBarry Smith 
2682d3f70b5SBarry Smith   return 0;
2692d3f70b5SBarry Smith }
2702d3f70b5SBarry Smith 
2712d3f70b5SBarry Smith /*
2722d3f70b5SBarry Smith    This constructs the Jacobian needed for SNES
2732d3f70b5SBarry Smith 
2742d3f70b5SBarry Smith              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
2752d3f70b5SBarry Smith */
2767bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
2772d3f70b5SBarry Smith {
2782d3f70b5SBarry Smith   TS      ts = (TS) ctx;
2792d3f70b5SBarry Smith   int     ierr;
2802d3f70b5SBarry Smith   Scalar  mone = -1.0, mdt = 1.0/ts->time_step;
2812d3f70b5SBarry Smith   MatType mtype;
2822d3f70b5SBarry Smith 
2832d3f70b5SBarry Smith   /* construct users Jacobian */
2842d3f70b5SBarry Smith   if (ts->rhsjacobian) {
2852d3f70b5SBarry Smith     ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr);
2862d3f70b5SBarry Smith   }
2872d3f70b5SBarry Smith 
2882d3f70b5SBarry Smith   /* shift and scale Jacobian, if not a shell matrix */
2892d3f70b5SBarry Smith   ierr = MatGetType(*AA,&mtype,PETSC_NULL);
2902d3f70b5SBarry Smith   if (mtype != MATSHELL) {
2912d3f70b5SBarry Smith     ierr = MatScale(&mone,*AA); CHKERRQ(ierr);
2922d3f70b5SBarry Smith     ierr = MatShift(&mdt,*AA); CHKERRQ(ierr);
2932d3f70b5SBarry Smith   }
2942d3f70b5SBarry Smith   ierr = MatGetType(*BB,&mtype,PETSC_NULL);
2952d3f70b5SBarry Smith   if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) {
2962d3f70b5SBarry Smith     ierr = MatScale(&mone,*BB); CHKERRQ(ierr);
2972d3f70b5SBarry Smith     ierr = MatShift(&mdt,*BB); CHKERRQ(ierr);
2982d3f70b5SBarry Smith   }
2992d3f70b5SBarry Smith 
3002d3f70b5SBarry Smith   return 0;
3012d3f70b5SBarry Smith }
3022d3f70b5SBarry Smith 
3032d3f70b5SBarry Smith 
3047bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts)
3052d3f70b5SBarry Smith {
3067bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
3072d3f70b5SBarry Smith   int         ierr, M, m;
3082d3f70b5SBarry Smith 
3097bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr);
3107bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr);
3117bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
3122d3f70b5SBarry Smith   if (ts->Ashell) { /* construct new shell matrix */
3132d3f70b5SBarry Smith     ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr);
3142d3f70b5SBarry Smith     ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr);
3152d3f70b5SBarry Smith     ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr);
3167bf11e45SBarry Smith     ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr);
3172d3f70b5SBarry Smith   }
3187bf11e45SBarry Smith   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
3192d3f70b5SBarry Smith   return 0;
3202d3f70b5SBarry Smith }
3212d3f70b5SBarry Smith /*------------------------------------------------------------*/
3222d3f70b5SBarry Smith 
3232d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
3242d3f70b5SBarry Smith {
3257bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
3262d3f70b5SBarry Smith 
3277bf11e45SBarry Smith   PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm);
3282d3f70b5SBarry Smith   return 0;
3292d3f70b5SBarry Smith }
3302d3f70b5SBarry Smith 
3317bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts)
3322d3f70b5SBarry Smith {
3332d3f70b5SBarry Smith   int ierr,flg;
3342d3f70b5SBarry Smith 
3352d3f70b5SBarry Smith   ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr);
3362d3f70b5SBarry Smith 
3372d3f70b5SBarry Smith   ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr);
3382d3f70b5SBarry Smith   if (flg) {
3392d3f70b5SBarry Smith     TSSetMonitor(ts,TSPseudoDefaultMonitor,0);
3402d3f70b5SBarry Smith   }
3412d3f70b5SBarry Smith   return 0;
3422d3f70b5SBarry Smith }
3432d3f70b5SBarry Smith 
3447bf11e45SBarry Smith static int TSPrintHelp_Pseudo(TS ts)
3452d3f70b5SBarry Smith {
3462d3f70b5SBarry Smith 
3472d3f70b5SBarry Smith   return 0;
3482d3f70b5SBarry Smith }
3492d3f70b5SBarry Smith 
3507bf11e45SBarry Smith static int TSView_Pseudo(PetscObject obj,Viewer viewer)
3512d3f70b5SBarry Smith {
3522d3f70b5SBarry Smith   return 0;
3532d3f70b5SBarry Smith }
3542d3f70b5SBarry Smith 
3552d3f70b5SBarry Smith /* ------------------------------------------------------------ */
3567bf11e45SBarry Smith int TSCreate_Pseudo(TS ts )
3572d3f70b5SBarry Smith {
3587bf11e45SBarry Smith   TS_Pseudo *pseudo;
3592d3f70b5SBarry Smith   int       ierr;
3602d3f70b5SBarry Smith   MatType   mtype;
3612d3f70b5SBarry Smith 
3627bf11e45SBarry Smith   ts->type 	      = TS_PSEUDO;
3637bf11e45SBarry Smith   ts->destroy         = TSDestroy_Pseudo;
3647bf11e45SBarry Smith   ts->printhelp       = TSPrintHelp_Pseudo;
3657bf11e45SBarry Smith   ts->view            = TSView_Pseudo;
3662d3f70b5SBarry Smith 
3672d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
3687bf11e45SBarry Smith     SETERRQ(1,"TSCreate_Pseudo:Only for nonlinear problems");
3692d3f70b5SBarry Smith   }
3702d3f70b5SBarry Smith   if (!ts->A) {
3717bf11e45SBarry Smith     SETERRQ(1,"TSCreate_Pseudo:Must set Jacobian");
3722d3f70b5SBarry Smith   }
3732d3f70b5SBarry Smith   ierr = MatGetType(ts->A,&mtype,PETSC_NULL);
3742d3f70b5SBarry Smith   if (mtype == MATSHELL) {
3752d3f70b5SBarry Smith     ts->Ashell = ts->A;
3762d3f70b5SBarry Smith   }
3777bf11e45SBarry Smith   ts->setup           = TSSetUp_Pseudo;
3787bf11e45SBarry Smith   ts->step            = TSStep_Pseudo;
3797bf11e45SBarry Smith   ts->setfromoptions  = TSSetFromOptions_Pseudo;
3807bf11e45SBarry Smith 
3817bf11e45SBarry Smith   /* create the required nonlinear solver context */
3822d3f70b5SBarry Smith   ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr);
3832d3f70b5SBarry Smith 
3847bf11e45SBarry Smith   pseudo   = PetscNew(TS_Pseudo); CHKPTRQ(pseudo);
3857bf11e45SBarry Smith   PetscMemzero(pseudo,sizeof(TS_Pseudo));
3867bf11e45SBarry Smith   ts->data = (void *) pseudo;
3872d3f70b5SBarry Smith 
3882d3f70b5SBarry Smith   return 0;
3892d3f70b5SBarry Smith }
3902d3f70b5SBarry Smith 
3912d3f70b5SBarry Smith 
3922d3f70b5SBarry Smith 
3932d3f70b5SBarry Smith 
3942d3f70b5SBarry Smith 
395