Actual source code: rosw.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: /*
  2:   Code for timestepping with Rosenbrock W methods

  4:   Notes:
  5:   The general system is written as

  7:   F(t,U,Udot) = G(t,U)

  9:   where F represents the stiff part of the physics and G represents the non-stiff part.
 10:   This method is designed to be linearly implicit on F and can use an approximate and lagged Jacobian.

 12: */
 13:  #include <petsc/private/tsimpl.h>
 14:  #include <petscdm.h>

 16:  #include <petsc/private/kernels/blockinvert.h>

 18: static TSRosWType TSRosWDefault = TSROSWRA34PW2;
 19: static PetscBool  TSRosWRegisterAllCalled;
 20: static PetscBool  TSRosWPackageInitialized;

 22: typedef struct _RosWTableau *RosWTableau;
 23: struct _RosWTableau {
 24:   char      *name;
 25:   PetscInt  order;              /* Classical approximation order of the method */
 26:   PetscInt  s;                  /* Number of stages */
 27:   PetscInt  pinterp;            /* Interpolation order */
 28:   PetscReal *A;                 /* Propagation table, strictly lower triangular */
 29:   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
 30:   PetscBool *GammaZeroDiag;     /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */
 31:   PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/
 32:   PetscReal *b;                 /* Step completion table */
 33:   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
 34:   PetscReal *ASum;              /* Row sum of A */
 35:   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
 36:   PetscReal *At;                /* Propagation table in transformed variables */
 37:   PetscReal *bt;                /* Step completion table in transformed variables */
 38:   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
 39:   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
 40:   PetscReal ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
 41:   PetscReal *binterpt;          /* Dense output formula */
 42: };
 43: typedef struct _RosWTableauLink *RosWTableauLink;
 44: struct _RosWTableauLink {
 45:   struct _RosWTableau tab;
 46:   RosWTableauLink     next;
 47: };
 48: static RosWTableauLink RosWTableauList;

 50: typedef struct {
 51:   RosWTableau  tableau;
 52:   Vec          *Y;               /* States computed during the step, used to complete the step */
 53:   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
 54:   Vec          Ystage;           /* Work vector for the state value at each stage */
 55:   Vec          Zdot;             /* Ydot = Zdot + shift*Y */
 56:   Vec          Zstage;           /* Y = Zstage + Y */
 57:   Vec          vec_sol_prev;     /* Solution from the previous step (used for interpolation and rollback)*/
 58:   PetscScalar  *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
 59:   PetscReal    scoeff;           /* shift = scoeff/dt */
 60:   PetscReal    stage_time;
 61:   PetscReal    stage_explicit;     /* Flag indicates that the current stage is explicit */
 62:   PetscBool    recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
 63:   TSStepStatus status;
 64: } TS_RosW;

 66: /*MC
 67:      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).

 69:      Only an approximate Jacobian is needed.

 71:      Level: intermediate

 73: .seealso: TSROSW
 74: M*/

 76: /*MC
 77:      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).

 79:      Only an approximate Jacobian is needed.

 81:      Level: intermediate

 83: .seealso: TSROSW
 84: M*/

 86: /*MC
 87:      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.

 89:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.

 91:      Level: intermediate

 93: .seealso: TSROSW
 94: M*/

 96: /*MC
 97:      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.

 99:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.

101:      Level: intermediate

103: .seealso: TSROSW
104: M*/

106: /*MC
107:      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.

109:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.

111:      This is strongly A-stable with R(infty) = 0.73. The embedded method of order 2 is strongly A-stable with R(infty) = 0.73.

113:      References:
114: .  1. -   Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.

116:      Level: intermediate

118: .seealso: TSROSW
119: M*/

121: /*MC
122:      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.

124:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.

126:      This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48.

128:      References:
129: .  1. -   Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.

131:      Level: intermediate

133: .seealso: TSROSW
134: M*/

136: /*MC
137:      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme

139:      By default, the Jacobian is only recomputed once per step.

141:      Both the third order and embedded second order methods are stiffly accurate and L-stable.

143:      References:
144: .  1. -   Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.

146:      Level: intermediate

148: .seealso: TSROSW, TSROSWSANDU3
149: M*/

151: /*MC
152:      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme

154:      By default, the Jacobian is only recomputed once per step.

156:      The third order method is L-stable, but not stiffly accurate.
157:      The second order embedded method is strongly A-stable with R(infty) = 0.5.
158:      The internal stages are L-stable.
159:      This method is called ROS3 in the paper.

161:      References:
162: .  1. -   Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.

164:      Level: intermediate

166: .seealso: TSROSW, TSROSWRODAS3
167: M*/

169: /*MC
170:      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages

172:      By default, the Jacobian is only recomputed once per step.

174:      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)

176:      References:
177: .     Emil Constantinescu

179:      Level: intermediate

181: .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP
182: M*/

184: /*MC
185:      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages

187:      By default, the Jacobian is only recomputed once per step.

189:      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)

191:      References:
192: .     Emil Constantinescu

194:      Level: intermediate

196: .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP
197: M*/

199: /*MC
200:      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages

202:      By default, the Jacobian is only recomputed once per step.

204:      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)

206:      References:
207: .     Emil Constantinescu

209:      Level: intermediate

211: .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP
212: M*/

214: /*MC
215:      TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop

217:      By default, the Jacobian is only recomputed once per step.

219:      A(89.3 degrees)-stable, |R(infty)| = 0.454.

221:      This method does not provide a dense output formula.

223:      References:
224: +   1. -  Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979.
225: -   2. -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

227:      Hairer's code ros4.f

229:      Level: intermediate

231: .seealso: TSROSW, TSROSWSHAMP4, TSROSWVELDD4, TSROSW4L
232: M*/

234: /*MC
235:      TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine

237:      By default, the Jacobian is only recomputed once per step.

239:      A-stable, |R(infty)| = 1/3.

241:      This method does not provide a dense output formula.

243:      References:
244: +   1. -  Shampine, Implementation of Rosenbrock methods, 1982.
245: -   2. -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

247:      Hairer's code ros4.f

249:      Level: intermediate

251: .seealso: TSROSW, TSROSWGRK4T, TSROSWVELDD4, TSROSW4L
252: M*/

254: /*MC
255:      TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen

257:      By default, the Jacobian is only recomputed once per step.

259:      A(89.5 degrees)-stable, |R(infty)| = 0.24.

261:      This method does not provide a dense output formula.

263:      References:
264: +   1. -  van Veldhuizen, D stability and Kaps Rentrop methods, 1984.
265: -   2. -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

267:      Hairer's code ros4.f

269:      Level: intermediate

271: .seealso: TSROSW, TSROSWGRK4T, TSROSWSHAMP4, TSROSW4L
272: M*/

274: /*MC
275:      TSROSW4L - four stage, fourth order Rosenbrock (not W) method

277:      By default, the Jacobian is only recomputed once per step.

279:      A-stable and L-stable

281:      This method does not provide a dense output formula.

283:      References:
284: .  1. -   Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

286:      Hairer's code ros4.f

288:      Level: intermediate

290: .seealso: TSROSW, TSROSWGRK4T, TSROSWSHAMP4, TSROSW4L
291: M*/

293: /*@C
294:   TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in TSRosW

296:   Not Collective, but should be called by all processes which will need the schemes to be registered

298:   Level: advanced

300: .keywords: TS, TSRosW, register, all

302: .seealso:  TSRosWRegisterDestroy()
303: @*/
304: PetscErrorCode TSRosWRegisterAll(void)
305: {

309:   if (TSRosWRegisterAllCalled) return(0);
310:   TSRosWRegisterAllCalled = PETSC_TRUE;

312:   {
313:     const PetscReal A = 0;
314:     const PetscReal Gamma = 1;
315:     const PetscReal b = 1;
316:     const PetscReal binterpt=1;

318:     TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,NULL,1,&binterpt);
319:   }

321:   {
322:     const PetscReal A = 0;
323:     const PetscReal Gamma = 0.5;
324:     const PetscReal b = 1;
325:     const PetscReal binterpt=1;

327:     TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,NULL,1,&binterpt);
328:   }

330:   {
331:     /*const PetscReal g = 1. + 1./PetscSqrtReal(2.0);   Direct evaluation: 1.707106781186547524401. Used for setting up arrays of values known at compile time below. */
332:     const PetscReal
333:       A[2][2]     = {{0,0}, {1.,0}},
334:       Gamma[2][2] = {{1.707106781186547524401,0}, {-2.*1.707106781186547524401,1.707106781186547524401}},
335:       b[2]        = {0.5,0.5},
336:       b1[2]       = {1.0,0.0};
337:     PetscReal binterpt[2][2];
338:     binterpt[0][0] = 1.707106781186547524401 - 1.0;
339:     binterpt[1][0] = 2.0 - 1.707106781186547524401;
340:     binterpt[0][1] = 1.707106781186547524401 - 1.5;
341:     binterpt[1][1] = 1.5 - 1.707106781186547524401;

343:     TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]);
344:   }
345:   {
346:     /*const PetscReal g = 1. - 1./PetscSqrtReal(2.0);   Direct evaluation: 0.2928932188134524755992. Used for setting up arrays of values known at compile time below. */
347:     const PetscReal
348:       A[2][2]     = {{0,0}, {1.,0}},
349:       Gamma[2][2] = {{0.2928932188134524755992,0}, {-2.*0.2928932188134524755992,0.2928932188134524755992}},
350:       b[2]        = {0.5,0.5},
351:       b1[2]       = {1.0,0.0};
352:     PetscReal binterpt[2][2];
353:     binterpt[0][0] = 0.2928932188134524755992 - 1.0;
354:     binterpt[1][0] = 2.0 - 0.2928932188134524755992;
355:     binterpt[0][1] = 0.2928932188134524755992 - 1.5;
356:     binterpt[1][1] = 1.5 - 0.2928932188134524755992;

358:     TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]);
359:   }
360:   {
361:     /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */
362:     PetscReal binterpt[3][2];
363:     const PetscReal
364:       A[3][3] = {{0,0,0},
365:                  {1.5773502691896257e+00,0,0},
366:                  {0.5,0,0}},
367:       Gamma[3][3] = {{7.8867513459481287e-01,0,0},
368:                      {-1.5773502691896257e+00,7.8867513459481287e-01,0},
369:                      {-6.7075317547305480e-01,-1.7075317547305482e-01,7.8867513459481287e-01}},
370:       b[3]  = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
371:       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};

373:       binterpt[0][0] = -0.8094010767585034;
374:       binterpt[1][0] = -0.5;
375:       binterpt[2][0] = 2.3094010767585034;
376:       binterpt[0][1] = 0.9641016151377548;
377:       binterpt[1][1] = 0.5;
378:       binterpt[2][1] = -1.4641016151377548;

380:       TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);
381:   }
382:   {
383:     PetscReal  binterpt[4][3];
384:     /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */
385:     const PetscReal
386:       A[4][4] = {{0,0,0,0},
387:                  {8.7173304301691801e-01,0,0,0},
388:                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
389:                  {0,0,1.,0}},
390:       Gamma[4][4] = {{4.3586652150845900e-01,0,0,0},
391:                      {-8.7173304301691801e-01,4.3586652150845900e-01,0,0},
392:                      {-9.0338057013044082e-01,5.4180672388095326e-02,4.3586652150845900e-01,0},
393:                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,4.3586652150845900e-01}},
394:       b[4]  = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
395:       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};

397:     binterpt[0][0]=1.0564298455794094;
398:     binterpt[1][0]=2.296429974281067;
399:     binterpt[2][0]=-1.307599564525376;
400:     binterpt[3][0]=-1.045260255335102;
401:     binterpt[0][1]=-1.3864882699759573;
402:     binterpt[1][1]=-8.262611700275677;
403:     binterpt[2][1]=7.250979895056055;
404:     binterpt[3][1]=2.398120075195581;
405:     binterpt[0][2]=0.5721822314575016;
406:     binterpt[1][2]=4.742931142090097;
407:     binterpt[2][2]=-4.398120075195578;
408:     binterpt[3][2]=-0.9169932983520199;

410:     TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);
411:   }
412:   {
413:     /* const PetscReal g = 0.5;       Directly written in-place below */
414:     const PetscReal
415:       A[4][4] = {{0,0,0,0},
416:                  {0,0,0,0},
417:                  {1.,0,0,0},
418:                  {0.75,-0.25,0.5,0}},
419:       Gamma[4][4] = {{0.5,0,0,0},
420:                      {1.,0.5,0,0},
421:                      {-0.25,-0.25,0.5,0},
422:                      {1./12,1./12,-2./3,0.5}},
423:       b[4]  = {5./6,-1./6,-1./6,0.5},
424:       b2[4] = {0.75,-0.25,0.5,0};

426:     TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2,0,NULL);
427:   }
428:   {
429:     /*const PetscReal g = 0.43586652150845899941601945119356;       Directly written in-place below */
430:     const PetscReal
431:       A[3][3] = {{0,0,0},
432:                  {0.43586652150845899941601945119356,0,0},
433:                  {0.43586652150845899941601945119356,0,0}},
434:       Gamma[3][3] = {{0.43586652150845899941601945119356,0,0},
435:                      {-0.19294655696029095575009695436041,0.43586652150845899941601945119356,0},
436:                      {0,1.74927148125794685173529749738960,0.43586652150845899941601945119356}},
437:       b[3]  = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
438:       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};

440:     PetscReal binterpt[3][2];
441:     binterpt[0][0] = 3.793692883777660870425141387941;
442:     binterpt[1][0] = -2.918692883777660870425141387941;
443:     binterpt[2][0] = 0.125;
444:     binterpt[0][1] = -0.725741064379812106687651020584;
445:     binterpt[1][1] = 0.559074397713145440020984353917;
446:     binterpt[2][1] = 0.16666666666666666666666666666667;

448:     TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);
449:   }
450:   {
451:     /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
452:      * Direct evaluation: s3 = 1.732050807568877293527;
453:      *                     g = 0.7886751345948128822546;
454:      * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */
455:     const PetscReal
456:       A[3][3] = {{0,0,0},
457:                  {1,0,0},
458:                  {0.25,0.25,0}},
459:       Gamma[3][3] = {{0,0,0},
460:                      {(-3.0-1.732050807568877293527)/6.0,0.7886751345948128822546,0},
461:                      {(-3.0-1.732050807568877293527)/24.0,(-3.0-1.732050807568877293527)/8.0,0.7886751345948128822546}},
462:       b[3]  = {1./6.,1./6.,2./3.},
463:       b2[3] = {1./4.,1./4.,1./2.};
464:     PetscReal binterpt[3][2];

466:     binterpt[0][0]=0.089316397477040902157517886164709;
467:     binterpt[1][0]=-0.91068360252295909784248211383529;
468:     binterpt[2][0]=1.8213672050459181956849642276706;
469:     binterpt[0][1]=0.077350269189625764509148780501957;
470:     binterpt[1][1]=1.077350269189625764509148780502;
471:     binterpt[2][1]=-1.1547005383792515290182975610039;

473:     TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);
474:   }

476:   {
477:     const PetscReal
478:       A[4][4] = {{0,0,0,0},
479:                  {1./2.,0,0,0},
480:                  {1./2.,1./2.,0,0},
481:                  {1./6.,1./6.,1./6.,0}},
482:       Gamma[4][4] = {{1./2.,0,0,0},
483:                      {0.0,1./4.,0,0},
484:                      {-2.,-2./3.,2./3.,0},
485:                      {1./2.,5./36.,-2./9,0}},
486:       b[4]  = {1./6.,1./6.,1./6.,1./2.},
487:       b2[4] = {1./8.,3./4.,1./8.,0};
488:     PetscReal binterpt[4][3];

490:     binterpt[0][0]=6.25;
491:     binterpt[1][0]=-30.25;
492:     binterpt[2][0]=1.75;
493:     binterpt[3][0]=23.25;
494:     binterpt[0][1]=-9.75;
495:     binterpt[1][1]=58.75;
496:     binterpt[2][1]=-3.25;
497:     binterpt[3][1]=-45.75;
498:     binterpt[0][2]=3.6666666666666666666666666666667;
499:     binterpt[1][2]=-28.333333333333333333333333333333;
500:     binterpt[2][2]=1.6666666666666666666666666666667;
501:     binterpt[3][2]=23.;

503:     TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);
504:   }

506:   {
507:     const PetscReal
508:       A[4][4] = {{0,0,0,0},
509:                  {1./2.,0,0,0},
510:                  {1./2.,1./2.,0,0},
511:                  {1./6.,1./6.,1./6.,0}},
512:       Gamma[4][4] = {{1./2.,0,0,0},
513:                      {0.0,3./4.,0,0},
514:                      {-2./3.,-23./9.,2./9.,0},
515:                      {1./18.,65./108.,-2./27,0}},
516:       b[4]  = {1./6.,1./6.,1./6.,1./2.},
517:       b2[4] = {3./16.,10./16.,3./16.,0};
518:     PetscReal binterpt[4][3];

520:     binterpt[0][0]=1.6911764705882352941176470588235;
521:     binterpt[1][0]=3.6813725490196078431372549019608;
522:     binterpt[2][0]=0.23039215686274509803921568627451;
523:     binterpt[3][0]=-4.6029411764705882352941176470588;
524:     binterpt[0][1]=-0.95588235294117647058823529411765;
525:     binterpt[1][1]=-6.2401960784313725490196078431373;
526:     binterpt[2][1]=-0.31862745098039215686274509803922;
527:     binterpt[3][1]=7.5147058823529411764705882352941;
528:     binterpt[0][2]=-0.56862745098039215686274509803922;
529:     binterpt[1][2]=2.7254901960784313725490196078431;
530:     binterpt[2][2]=0.25490196078431372549019607843137;
531:     binterpt[3][2]=-2.4117647058823529411764705882353;

533:     TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);
534:   }

536:   {
537:     PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
538:     PetscReal binterpt[4][3];

540:     Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
541:     Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
542:     Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
543:     Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
544:     Gamma[1][2]=0; Gamma[1][3]=0;
545:     Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
546:     Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
547:     Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
548:     Gamma[2][3]=0;
549:     Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
550:     Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
551:     Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
552:     Gamma[3][3]=0;

554:     A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
555:     A[1][0]=0.8717330430169179988320388950590125027645343373957631;
556:     A[1][1]=0; A[1][2]=0; A[1][3]=0;
557:     A[2][0]=0.5275890119763004115618079766722914408876108660811028;
558:     A[2][1]=0.07241098802369958843819203208518599088698057726988732;
559:     A[2][2]=0; A[2][3]=0;
560:     A[3][0]=0.3990960076760701320627260685975778145384666450351314;
561:     A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
562:     A[3][2]=1.038461646937449311660120300601880176655352737312713;
563:     A[3][3]=0;

565:     b[0]=0.1876410243467238251612921333138006734899663569186926;
566:     b[1]=-0.5952974735769549480478230473706443582188442040780541;
567:     b[2]=0.9717899277217721234705114616271378792182450260943198;
568:     b[3]=0.4358665215084589994160194475295062513822671686978816;

570:     b2[0]=0.2147402862233891404862383521089097657790734483804460;
571:     b2[1]=-0.4851622638849390928209050538171743017757490232519684;
572:     b2[2]=0.8687250025203875511662123688667549217531982787600080;
573:     b2[3]=0.4016969751411624011684543450940068201770721128357014;

575:     binterpt[0][0]=2.2565812720167954547104627844105;
576:     binterpt[1][0]=1.349166413351089573796243820819;
577:     binterpt[2][0]=-2.4695174540533503758652847586647;
578:     binterpt[3][0]=-0.13623023131453465264142184656474;
579:     binterpt[0][1]=-3.0826699111559187902922463354557;
580:     binterpt[1][1]=-2.4689115685996042534544925650515;
581:     binterpt[2][1]=5.7428279814696677152129332773553;
582:     binterpt[3][1]=-0.19124650171414467146619437684812;
583:     binterpt[0][2]=1.0137296634858471607430756831148;
584:     binterpt[1][2]=0.52444768167155973161042570784064;
585:     binterpt[2][2]=-2.3015205996945452158771370439586;
586:     binterpt[3][2]=0.76334325453713832352363565300308;

588:     TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);
589:   }
590:   TSRosWRegisterRos4(TSROSWGRK4T,0.231,PETSC_DEFAULT,PETSC_DEFAULT,0,-0.1282612945269037e+01);
591:   TSRosWRegisterRos4(TSROSWSHAMP4,0.5,PETSC_DEFAULT,PETSC_DEFAULT,0,125./108.);
592:   TSRosWRegisterRos4(TSROSWVELDD4,0.22570811482256823492,PETSC_DEFAULT,PETSC_DEFAULT,0,-1.355958941201148);
593:   TSRosWRegisterRos4(TSROSW4L,0.57282,PETSC_DEFAULT,PETSC_DEFAULT,0,-1.093502252409163);
594:   return(0);
595: }



599: /*@C
600:    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().

602:    Not Collective

604:    Level: advanced

606: .keywords: TSRosW, register, destroy
607: .seealso: TSRosWRegister(), TSRosWRegisterAll()
608: @*/
609: PetscErrorCode TSRosWRegisterDestroy(void)
610: {
611:   PetscErrorCode  ierr;
612:   RosWTableauLink link;

615:   while ((link = RosWTableauList)) {
616:     RosWTableau t = &link->tab;
617:     RosWTableauList = link->next;
618:     PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);
619:     PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);
620:     PetscFree2(t->bembed,t->bembedt);
621:     PetscFree(t->binterpt);
622:     PetscFree(t->name);
623:     PetscFree(link);
624:   }
625:   TSRosWRegisterAllCalled = PETSC_FALSE;
626:   return(0);
627: }

629: /*@C
630:   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
631:   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
632:   when using static libraries.

634:   Level: developer

636: .keywords: TS, TSRosW, initialize, package
637: .seealso: PetscInitialize()
638: @*/
639: PetscErrorCode TSRosWInitializePackage(void)
640: {

644:   if (TSRosWPackageInitialized) return(0);
645:   TSRosWPackageInitialized = PETSC_TRUE;
646:   TSRosWRegisterAll();
647:   PetscRegisterFinalize(TSRosWFinalizePackage);
648:   return(0);
649: }

651: /*@C
652:   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
653:   called from PetscFinalize().

655:   Level: developer

657: .keywords: Petsc, destroy, package
658: .seealso: PetscFinalize()
659: @*/
660: PetscErrorCode TSRosWFinalizePackage(void)
661: {

665:   TSRosWPackageInitialized = PETSC_FALSE;
666:   TSRosWRegisterDestroy();
667:   return(0);
668: }

670: /*@C
671:    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation

673:    Not Collective, but the same schemes should be registered on all processes on which they will be used

675:    Input Parameters:
676: +  name - identifier for method
677: .  order - approximation order of method
678: .  s - number of stages, this is the dimension of the matrices below
679: .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
680: .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
681: .  b - Step completion table (dimension s)
682: .  bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
683: .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
684: -  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)

686:    Notes:
687:    Several Rosenbrock W methods are provided, this function is only needed to create new methods.

689:    Level: advanced

691: .keywords: TS, register

693: .seealso: TSRosW
694: @*/
695: PetscErrorCode TSRosWRegister(TSRosWType name,PetscInt order,PetscInt s,const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[],
696:                               PetscInt pinterp,const PetscReal binterpt[])
697: {
698:   PetscErrorCode  ierr;
699:   RosWTableauLink link;
700:   RosWTableau     t;
701:   PetscInt        i,j,k;
702:   PetscScalar     *GammaInv;


711:   PetscCalloc1(1,&link);
712:   t        = &link->tab;
713:   PetscStrallocpy(name,&t->name);
714:   t->order = order;
715:   t->s     = s;
716:   PetscMalloc5(s*s,&t->A,s*s,&t->Gamma,s,&t->b,s,&t->ASum,s,&t->GammaSum);
717:   PetscMalloc5(s*s,&t->At,s,&t->bt,s*s,&t->GammaInv,s,&t->GammaZeroDiag,s*s,&t->GammaExplicitCorr);
718:   PetscMemcpy(t->A,A,s*s*sizeof(A[0]));
719:   PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));
720:   PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));
721:   PetscMemcpy(t->b,b,s*sizeof(b[0]));
722:   if (bembed) {
723:     PetscMalloc2(s,&t->bembed,s,&t->bembedt);
724:     PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));
725:   }
726:   for (i=0; i<s; i++) {
727:     t->ASum[i]     = 0;
728:     t->GammaSum[i] = 0;
729:     for (j=0; j<s; j++) {
730:       t->ASum[i]     += A[i*s+j];
731:       t->GammaSum[i] += Gamma[i*s+j];
732:     }
733:   }
734:   PetscMalloc1(s*s,&GammaInv); /* Need to use Scalar for inverse, then convert back to Real */
735:   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
736:   for (i=0; i<s; i++) {
737:     if (Gamma[i*s+i] == 0.0) {
738:       GammaInv[i*s+i] = 1.0;
739:       t->GammaZeroDiag[i] = PETSC_TRUE;
740:     } else {
741:       t->GammaZeroDiag[i] = PETSC_FALSE;
742:     }
743:   }

745:   switch (s) {
746:   case 1: GammaInv[0] = 1./GammaInv[0]; break;
747:   case 2: PetscKernel_A_gets_inverse_A_2(GammaInv,0,PETSC_FALSE,NULL); break;
748:   case 3: PetscKernel_A_gets_inverse_A_3(GammaInv,0,PETSC_FALSE,NULL); break;
749:   case 4: PetscKernel_A_gets_inverse_A_4(GammaInv,0,PETSC_FALSE,NULL); break;
750:   case 5: {
751:     PetscInt  ipvt5[5];
752:     MatScalar work5[5*5];
753:     PetscKernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0,PETSC_FALSE,NULL); break;
754:   }
755:   case 6: PetscKernel_A_gets_inverse_A_6(GammaInv,0,PETSC_FALSE,NULL); break;
756:   case 7: PetscKernel_A_gets_inverse_A_7(GammaInv,0,PETSC_FALSE,NULL); break;
757:   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
758:   }
759:   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
760:   PetscFree(GammaInv);

762:   for (i=0; i<s; i++) {
763:     for (k=0; k<i+1; k++) {
764:       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
765:       for (j=k+1; j<i+1; j++) {
766:         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
767:       }
768:     }
769:   }

771:   for (i=0; i<s; i++) {
772:     for (j=0; j<s; j++) {
773:       t->At[i*s+j] = 0;
774:       for (k=0; k<s; k++) {
775:         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
776:       }
777:     }
778:     t->bt[i] = 0;
779:     for (j=0; j<s; j++) {
780:       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
781:     }
782:     if (bembed) {
783:       t->bembedt[i] = 0;
784:       for (j=0; j<s; j++) {
785:         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
786:       }
787:     }
788:   }
789:   t->ccfl = 1.0;                /* Fix this */

791:   t->pinterp = pinterp;
792:   PetscMalloc1(s*pinterp,&t->binterpt);
793:   PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));
794:   link->next = RosWTableauList;
795:   RosWTableauList = link;
796:   return(0);
797: }

799: /*@C
800:    TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing paramter choices

802:    Not Collective, but the same schemes should be registered on all processes on which they will be used

804:    Input Parameters:
805: +  name - identifier for method
806: .  gamma - leading coefficient (diagonal entry)
807: .  a2 - design parameter, see Table 7.2 of Hairer&Wanner
808: .  a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
809: .  b3 - design parameter, see Table 7.2 of Hairer&Wanner
810: .  beta43 - design parameter or PETSC_DEFAULT to use Equation 7.21 of Hairer&Wanner
811: .  e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer

813:    Notes:
814:    This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2.
815:    It is used here to implement several methods from the book and can be used to experiment with new methods.
816:    It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions.

818:    Level: developer

820: .keywords: TS, register

822: .seealso: TSRosW, TSRosWRegister()
823: @*/
824: PetscErrorCode TSRosWRegisterRos4(TSRosWType name,PetscReal gamma,PetscReal a2,PetscReal a3,PetscReal b3,PetscReal e4)
825: {
827:   /* Declare numeric constants so they can be quad precision without being truncated at double */
828:   const PetscReal one = 1,two = 2,three = 3,four = 4,five = 5,six = 6,eight = 8,twelve = 12,twenty = 20,twentyfour = 24,
829:     p32 = one/six - gamma + gamma*gamma,
830:     p42 = one/eight - gamma/three,
831:     p43 = one/twelve - gamma/three,
832:     p44 = one/twentyfour - gamma/two + three/two*gamma*gamma - gamma*gamma*gamma,
833:     p56 = one/twenty - gamma/four;
834:   PetscReal   a4,a32,a42,a43,b1,b2,b4,beta2p,beta3p,beta4p,beta32,beta42,beta43,beta32beta2p,beta4jbetajp;
835:   PetscReal   A[4][4],Gamma[4][4],b[4],bm[4];
836:   PetscScalar M[3][3],rhs[3];

839:   /* Step 1: choose Gamma (input) */
840:   /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */
841:   if (a3 == PETSC_DEFAULT) a3 = (one/five - a2/four)/(one/four - a2/three); /* Eq 7.22 */
842:   a4 = a3;                                                  /* consequence of 7.20 */

844:   /* Solve order conditions 7.15a, 7.15c, 7.15e */
845:   M[0][0] = one; M[0][1] = one;      M[0][2] = one;      /* 7.15a */
846:   M[1][0] = 0.0; M[1][1] = a2*a2;    M[1][2] = a4*a4;    /* 7.15c */
847:   M[2][0] = 0.0; M[2][1] = a2*a2*a2; M[2][2] = a4*a4*a4; /* 7.15e */
848:   rhs[0]  = one - b3;
849:   rhs[1]  = one/three - a3*a3*b3;
850:   rhs[2]  = one/four - a3*a3*a3*b3;
851:   PetscKernel_A_gets_inverse_A_3(&M[0][0],0,PETSC_FALSE,NULL);
852:   b1      = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
853:   b2      = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
854:   b4      = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);

856:   /* Step 3 */
857:   beta43       = (p56 - a2*p43) / (b4*a3*a3*(a3 - a2)); /* 7.21 */
858:   beta32beta2p =  p44 / (b4*beta43);                    /* 7.15h */
859:   beta4jbetajp = (p32 - b3*beta32beta2p) / b4;
860:   M[0][0]      = b2;                                    M[0][1] = b3;                 M[0][2] = b4;
861:   M[1][0]      = a4*a4*beta32beta2p-a3*a3*beta4jbetajp; M[1][1] = a2*a2*beta4jbetajp; M[1][2] = -a2*a2*beta32beta2p;
862:   M[2][0]      = b4*beta43*a3*a3-p43;                   M[2][1] = -b4*beta43*a2*a2;   M[2][2] = 0;
863:   rhs[0]       = one/two - gamma; rhs[1] = 0; rhs[2] = -a2*a2*p32;
864:   PetscKernel_A_gets_inverse_A_3(&M[0][0],0,PETSC_FALSE,NULL);
865:   beta2p       = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
866:   beta3p       = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
867:   beta4p       = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);

869:   /* Step 4: back-substitute */
870:   beta32 = beta32beta2p / beta2p;
871:   beta42 = (beta4jbetajp - beta43*beta3p) / beta2p;

873:   /* Step 5: 7.15f and 7.20, then 7.16 */
874:   a43 = 0;
875:   a32 = p42 / (b3*a3*beta2p + b4*a4*beta2p);
876:   a42 = a32;

878:   A[0][0]     = 0;          A[0][1] = 0;   A[0][2] = 0;   A[0][3] = 0;
879:   A[1][0]     = a2;         A[1][1] = 0;   A[1][2] = 0;   A[1][3] = 0;
880:   A[2][0]     = a3-a32;     A[2][1] = a32; A[2][2] = 0;   A[2][3] = 0;
881:   A[3][0]     = a4-a43-a42; A[3][1] = a42; A[3][2] = a43; A[3][3] = 0;
882:   Gamma[0][0] = gamma;                        Gamma[0][1] = 0;              Gamma[0][2] = 0;              Gamma[0][3] = 0;
883:   Gamma[1][0] = beta2p-A[1][0];               Gamma[1][1] = gamma;          Gamma[1][2] = 0;              Gamma[1][3] = 0;
884:   Gamma[2][0] = beta3p-beta32-A[2][0];        Gamma[2][1] = beta32-A[2][1]; Gamma[2][2] = gamma;          Gamma[2][3] = 0;
885:   Gamma[3][0] = beta4p-beta42-beta43-A[3][0]; Gamma[3][1] = beta42-A[3][1]; Gamma[3][2] = beta43-A[3][2]; Gamma[3][3] = gamma;
886:   b[0] = b1; b[1] = b2; b[2] = b3; b[3] = b4;

888:   /* Construct embedded formula using given e4. We are solving Equation 7.18. */
889:   bm[3] = b[3] - e4*gamma;                                          /* using definition of E4 */
890:   bm[2] = (p32 - beta4jbetajp*bm[3]) / (beta32*beta2p);             /* fourth row of 7.18 */
891:   bm[1] = (one/two - gamma - beta3p*bm[2] - beta4p*bm[3]) / beta2p; /* second row */
892:   bm[0] = one - bm[1] - bm[2] - bm[3];                              /* first row */

894:   {
895:     const PetscReal misfit = a2*a2*bm[1] + a3*a3*bm[2] + a4*a4*bm[3] - one/three;
896:     if (PetscAbs(misfit) > PETSC_SMALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Assumptions violated, could not construct a third order embedded method");
897:   }
898:   TSRosWRegister(name,4,4,&A[0][0],&Gamma[0][0],b,bm,0,NULL);
899:   return(0);
900: }

902: /*
903:  The step completion formula is

905:  x1 = x0 + b^T Y

907:  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
908:  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write

910:  x1e = x0 + be^T Y
911:      = x1 - b^T Y + be^T Y
912:      = x1 + (be - b)^T Y

914:  so we can evaluate the method of different order even after the step has been optimistically completed.
915: */
916: static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec U,PetscBool *done)
917: {
918:   TS_RosW        *ros = (TS_RosW*)ts->data;
919:   RosWTableau    tab  = ros->tableau;
920:   PetscScalar    *w   = ros->work;
921:   PetscInt       i;

925:   if (order == tab->order) {
926:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
927:       VecCopy(ts->vec_sol,U);
928:       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
929:       VecMAXPY(U,tab->s,w,ros->Y);
930:     } else {VecCopy(ts->vec_sol,U);}
931:     if (done) *done = PETSC_TRUE;
932:     return(0);
933:   } else if (order == tab->order-1) {
934:     if (!tab->bembedt) goto unavailable;
935:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
936:       VecCopy(ts->vec_sol,U);
937:       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
938:       VecMAXPY(U,tab->s,w,ros->Y);
939:     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
940:       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
941:       VecCopy(ts->vec_sol,U);
942:       VecMAXPY(U,tab->s,w,ros->Y);
943:     }
944:     if (done) *done = PETSC_TRUE;
945:     return(0);
946:   }
947:   unavailable:
948:   if (done) *done = PETSC_FALSE;
949:   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Rosenbrock-W '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
950:   return(0);
951: }

953: static PetscErrorCode TSRollBack_RosW(TS ts)
954: {
955:   TS_RosW        *ros = (TS_RosW*)ts->data;

959:   VecCopy(ros->vec_sol_prev,ts->vec_sol);
960:   return(0);
961: }

963: static PetscErrorCode TSStep_RosW(TS ts)
964: {
965:   TS_RosW         *ros = (TS_RosW*)ts->data;
966:   RosWTableau     tab  = ros->tableau;
967:   const PetscInt  s    = tab->s;
968:   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
969:   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
970:   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
971:   PetscScalar     *w   = ros->work;
972:   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
973:   SNES            snes;
974:   TSAdapt         adapt;
975:   PetscInt        i,j,its,lits;
976:   PetscInt        rejections = 0;
977:   PetscBool       stageok,accept = PETSC_TRUE;
978:   PetscReal       next_time_step = ts->time_step;
979:   PetscErrorCode  ierr;

982:   if (!ts->steprollback) {
983:     VecCopy(ts->vec_sol,ros->vec_sol_prev);
984:   }

986:   ros->status = TS_STEP_INCOMPLETE;
987:   while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
988:     const PetscReal h = ts->time_step;
989:     for (i=0; i<s; i++) {
990:       ros->stage_time = ts->ptime + h*ASum[i];
991:       TSPreStage(ts,ros->stage_time);
992:       if (GammaZeroDiag[i]) {
993:         ros->stage_explicit = PETSC_TRUE;
994:         ros->scoeff         = 1.;
995:       } else {
996:         ros->stage_explicit = PETSC_FALSE;
997:         ros->scoeff         = 1./Gamma[i*s+i];
998:       }

1000:       VecCopy(ts->vec_sol,Zstage);
1001:       for (j=0; j<i; j++) w[j] = At[i*s+j];
1002:       VecMAXPY(Zstage,i,w,Y);

1004:       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
1005:       VecZeroEntries(Zdot);
1006:       VecMAXPY(Zdot,i,w,Y);

1008:       /* Initial guess taken from last stage */
1009:       VecZeroEntries(Y[i]);

1011:       if (!ros->stage_explicit) {
1012:         TSGetSNES(ts,&snes);
1013:         if (!ros->recompute_jacobian && !i) {
1014:           SNESSetLagJacobian(snes,-2); /* Recompute the Jacobian on this solve, but not again */
1015:         }
1016:         SNESSolve(snes,NULL,Y[i]);
1017:         SNESGetIterationNumber(snes,&its);
1018:         SNESGetLinearSolveIterations(snes,&lits);
1019:         ts->snes_its += its; ts->ksp_its += lits;
1020:       } else {
1021:         Mat J,Jp;
1022:         VecZeroEntries(Ydot); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
1023:         TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);
1024:         VecScale(Y[i],-1.0);
1025:         VecAXPY(Y[i],-1.0,Zdot); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/

1027:         VecZeroEntries(Zstage); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
1028:         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
1029:         VecMAXPY(Zstage,i,w,Y);

1031:         /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
1032:         TSGetIJacobian(ts,&J,&Jp,NULL,NULL);
1033:         TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,J,Jp,PETSC_FALSE);
1034:         MatMult(J,Zstage,Zdot);
1035:         VecAXPY(Y[i],-1.0,Zdot);
1036:         ts->ksp_its += 1;

1038:         VecScale(Y[i],h);
1039:       }
1040:       TSPostStage(ts,ros->stage_time,i,Y);
1041:       TSGetAdapt(ts,&adapt);
1042:       TSAdaptCheckStage(adapt,ts,ros->stage_time,Y[i],&stageok);
1043:       if (!stageok) goto reject_step;
1044:     }

1046:     ros->status = TS_STEP_INCOMPLETE;
1047:     TSEvaluateStep_RosW(ts,tab->order,ts->vec_sol,NULL);
1048:     ros->status = TS_STEP_PENDING;
1049:     TSGetAdapt(ts,&adapt);
1050:     TSAdaptCandidatesClear(adapt);
1051:     TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);
1052:     TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);
1053:     ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1054:     if (!accept) { /* Roll back the current step */
1055:       TSRollBack_RosW(ts);
1056:       ts->time_step = next_time_step;
1057:       goto reject_step;
1058:     }

1060:     ts->ptime += ts->time_step;
1061:     ts->time_step = next_time_step;
1062:     break;

1064:   reject_step:
1065:     ts->reject++; accept = PETSC_FALSE;
1066:     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1067:       ts->reason = TS_DIVERGED_STEP_REJECTED;
1068:       PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);
1069:     }
1070:   }
1071:   return(0);
1072: }

1074: static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec U)
1075: {
1076:   TS_RosW         *ros = (TS_RosW*)ts->data;
1077:   PetscInt        s    = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j;
1078:   PetscReal       h;
1079:   PetscReal       tt,t;
1080:   PetscScalar     *bt;
1081:   const PetscReal *Bt = ros->tableau->binterpt;
1082:   PetscErrorCode  ierr;
1083:   const PetscReal *GammaInv = ros->tableau->GammaInv;
1084:   PetscScalar     *w        = ros->work;
1085:   Vec             *Y        = ros->Y;

1088:   if (!Bt) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);

1090:   switch (ros->status) {
1091:   case TS_STEP_INCOMPLETE:
1092:   case TS_STEP_PENDING:
1093:     h = ts->time_step;
1094:     t = (itime - ts->ptime)/h;
1095:     break;
1096:   case TS_STEP_COMPLETE:
1097:     h = ts->ptime - ts->ptime_prev;
1098:     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
1099:     break;
1100:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
1101:   }
1102:   PetscMalloc1(s,&bt);
1103:   for (i=0; i<s; i++) bt[i] = 0;
1104:   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
1105:     for (i=0; i<s; i++) {
1106:       bt[i] += Bt[i*pinterp+j] * tt;
1107:     }
1108:   }

1110:   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1111:   /* U <- 0*/
1112:   VecZeroEntries(U);
1113:   /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
1114:   for (j=0; j<s; j++) w[j] = 0;
1115:   for (j=0; j<s; j++) {
1116:     for (i=j; i<s; i++) {
1117:       w[j] +=  bt[i]*GammaInv[i*s+j];
1118:     }
1119:   }
1120:   VecMAXPY(U,i,w,Y);
1121:   /* U <- y(t) + U */
1122:   VecAXPY(U,1,ros->vec_sol_prev);

1124:   PetscFree(bt);
1125:   return(0);
1126: }

1128: /*------------------------------------------------------------*/

1130: static PetscErrorCode TSRosWTableauReset(TS ts)
1131: {
1132:   TS_RosW        *ros = (TS_RosW*)ts->data;
1133:   RosWTableau    tab  = ros->tableau;

1137:   if (!tab) return(0);
1138:   VecDestroyVecs(tab->s,&ros->Y);
1139:   PetscFree(ros->work);
1140:   return(0);
1141: }

1143: static PetscErrorCode TSReset_RosW(TS ts)
1144: {
1145:   TS_RosW        *ros = (TS_RosW*)ts->data;

1149:   TSRosWTableauReset(ts);
1150:   VecDestroy(&ros->Ydot);
1151:   VecDestroy(&ros->Ystage);
1152:   VecDestroy(&ros->Zdot);
1153:   VecDestroy(&ros->Zstage);
1154:   VecDestroy(&ros->vec_sol_prev);
1155:   return(0);
1156: }

1158: static PetscErrorCode TSRosWGetVecs(TS ts,DM dm,Vec *Ydot,Vec *Zdot,Vec *Ystage,Vec *Zstage)
1159: {
1160:   TS_RosW        *rw = (TS_RosW*)ts->data;

1164:   if (Ydot) {
1165:     if (dm && dm != ts->dm) {
1166:       DMGetNamedGlobalVector(dm,"TSRosW_Ydot",Ydot);
1167:     } else *Ydot = rw->Ydot;
1168:   }
1169:   if (Zdot) {
1170:     if (dm && dm != ts->dm) {
1171:       DMGetNamedGlobalVector(dm,"TSRosW_Zdot",Zdot);
1172:     } else *Zdot = rw->Zdot;
1173:   }
1174:   if (Ystage) {
1175:     if (dm && dm != ts->dm) {
1176:       DMGetNamedGlobalVector(dm,"TSRosW_Ystage",Ystage);
1177:     } else *Ystage = rw->Ystage;
1178:   }
1179:   if (Zstage) {
1180:     if (dm && dm != ts->dm) {
1181:       DMGetNamedGlobalVector(dm,"TSRosW_Zstage",Zstage);
1182:     } else *Zstage = rw->Zstage;
1183:   }
1184:   return(0);
1185: }


1188: static PetscErrorCode TSRosWRestoreVecs(TS ts,DM dm,Vec *Ydot,Vec *Zdot, Vec *Ystage, Vec *Zstage)
1189: {

1193:   if (Ydot) {
1194:     if (dm && dm != ts->dm) {
1195:       DMRestoreNamedGlobalVector(dm,"TSRosW_Ydot",Ydot);
1196:     }
1197:   }
1198:   if (Zdot) {
1199:     if (dm && dm != ts->dm) {
1200:       DMRestoreNamedGlobalVector(dm,"TSRosW_Zdot",Zdot);
1201:     }
1202:   }
1203:   if (Ystage) {
1204:     if (dm && dm != ts->dm) {
1205:       DMRestoreNamedGlobalVector(dm,"TSRosW_Ystage",Ystage);
1206:     }
1207:   }
1208:   if (Zstage) {
1209:     if (dm && dm != ts->dm) {
1210:       DMRestoreNamedGlobalVector(dm,"TSRosW_Zstage",Zstage);
1211:     }
1212:   }
1213:   return(0);
1214: }

1216: static PetscErrorCode DMCoarsenHook_TSRosW(DM fine,DM coarse,void *ctx)
1217: {
1219:   return(0);
1220: }

1222: static PetscErrorCode DMRestrictHook_TSRosW(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1223: {
1224:   TS             ts = (TS)ctx;
1226:   Vec            Ydot,Zdot,Ystage,Zstage;
1227:   Vec            Ydotc,Zdotc,Ystagec,Zstagec;

1230:   TSRosWGetVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage);
1231:   TSRosWGetVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec);
1232:   MatRestrict(restrct,Ydot,Ydotc);
1233:   VecPointwiseMult(Ydotc,rscale,Ydotc);
1234:   MatRestrict(restrct,Ystage,Ystagec);
1235:   VecPointwiseMult(Ystagec,rscale,Ystagec);
1236:   MatRestrict(restrct,Zdot,Zdotc);
1237:   VecPointwiseMult(Zdotc,rscale,Zdotc);
1238:   MatRestrict(restrct,Zstage,Zstagec);
1239:   VecPointwiseMult(Zstagec,rscale,Zstagec);
1240:   TSRosWRestoreVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage);
1241:   TSRosWRestoreVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec);
1242:   return(0);
1243: }


1246: static PetscErrorCode DMSubDomainHook_TSRosW(DM fine,DM coarse,void *ctx)
1247: {
1249:   return(0);
1250: }

1252: static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1253: {
1254:   TS             ts = (TS)ctx;
1256:   Vec            Ydot,Zdot,Ystage,Zstage;
1257:   Vec            Ydots,Zdots,Ystages,Zstages;

1260:   TSRosWGetVecs(ts,dm,&Ydot,&Ystage,&Zdot,&Zstage);
1261:   TSRosWGetVecs(ts,subdm,&Ydots,&Ystages,&Zdots,&Zstages);

1263:   VecScatterBegin(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD);
1264:   VecScatterEnd(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD);

1266:   VecScatterBegin(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD);
1267:   VecScatterEnd(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD);

1269:   VecScatterBegin(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD);
1270:   VecScatterEnd(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD);

1272:   VecScatterBegin(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD);
1273:   VecScatterEnd(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD);

1275:   TSRosWRestoreVecs(ts,dm,&Ydot,&Ystage,&Zdot,&Zstage);
1276:   TSRosWRestoreVecs(ts,subdm,&Ydots,&Ystages,&Zdots,&Zstages);
1277:   return(0);
1278: }

1280: /*
1281:   This defines the nonlinear equation that is to be solved with SNES
1282:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1283: */
1284: static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec U,Vec F,TS ts)
1285: {
1286:   TS_RosW        *ros = (TS_RosW*)ts->data;
1288:   Vec            Ydot,Zdot,Ystage,Zstage;
1289:   PetscReal      shift = ros->scoeff / ts->time_step;
1290:   DM             dm,dmsave;

1293:   SNESGetDM(snes,&dm);
1294:   TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1295:   VecWAXPY(Ydot,shift,U,Zdot);    /* Ydot = shift*U + Zdot */
1296:   VecWAXPY(Ystage,1.0,U,Zstage);  /* Ystage = U + Zstage */
1297:   dmsave = ts->dm;
1298:   ts->dm = dm;
1299:   TSComputeIFunction(ts,ros->stage_time,Ystage,Ydot,F,PETSC_FALSE);
1300:   ts->dm = dmsave;
1301:   TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1302:   return(0);
1303: }

1305: static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec U,Mat A,Mat B,TS ts)
1306: {
1307:   TS_RosW        *ros = (TS_RosW*)ts->data;
1308:   Vec            Ydot,Zdot,Ystage,Zstage;
1309:   PetscReal      shift = ros->scoeff / ts->time_step;
1311:   DM             dm,dmsave;

1314:   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
1315:   SNESGetDM(snes,&dm);
1316:   TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1317:   dmsave = ts->dm;
1318:   ts->dm = dm;
1319:   TSComputeIJacobian(ts,ros->stage_time,Ystage,Ydot,shift,A,B,PETSC_TRUE);
1320:   ts->dm = dmsave;
1321:   TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1322:   return(0);
1323: }

1325: static PetscErrorCode TSRosWTableauSetUp(TS ts)
1326: {
1327:   TS_RosW        *ros = (TS_RosW*)ts->data;
1328:   RosWTableau    tab  = ros->tableau;

1332:   VecDuplicateVecs(ts->vec_sol,tab->s,&ros->Y);
1333:   PetscMalloc1(tab->s,&ros->work);
1334:   return(0);
1335: }

1337: static PetscErrorCode TSSetUp_RosW(TS ts)
1338: {
1339:   TS_RosW        *ros = (TS_RosW*)ts->data;
1341:   DM             dm;
1342:   SNES           snes;

1345:   TSRosWTableauSetUp(ts);
1346:   VecDuplicate(ts->vec_sol,&ros->Ydot);
1347:   VecDuplicate(ts->vec_sol,&ros->Ystage);
1348:   VecDuplicate(ts->vec_sol,&ros->Zdot);
1349:   VecDuplicate(ts->vec_sol,&ros->Zstage);
1350:   VecDuplicate(ts->vec_sol,&ros->vec_sol_prev);
1351:   TSGetDM(ts,&dm);
1352:   DMCoarsenHookAdd(dm,DMCoarsenHook_TSRosW,DMRestrictHook_TSRosW,ts);
1353:   DMSubDomainHookAdd(dm,DMSubDomainHook_TSRosW,DMSubDomainRestrictHook_TSRosW,ts);
1354:   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1355:   TSGetSNES(ts,&snes);
1356:   if (!((PetscObject)snes)->type_name) {
1357:     SNESSetType(snes,SNESKSPONLY);
1358:   }
1359:   return(0);
1360: }
1361: /*------------------------------------------------------------*/

1363: static PetscErrorCode TSSetFromOptions_RosW(PetscOptionItems *PetscOptionsObject,TS ts)
1364: {
1365:   TS_RosW        *ros = (TS_RosW*)ts->data;
1367:   SNES           snes;

1370:   PetscOptionsHead(PetscOptionsObject,"RosW ODE solver options");
1371:   {
1372:     RosWTableauLink link;
1373:     PetscInt        count,choice;
1374:     PetscBool       flg;
1375:     const char      **namelist;

1377:     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
1378:     PetscMalloc1(count,(char***)&namelist);
1379:     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1380:     PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,ros->tableau->name,&choice,&flg);
1381:     if (flg) {TSRosWSetType(ts,namelist[choice]);}
1382:     PetscFree(namelist);

1384:     PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,NULL);
1385:   }
1386:   PetscOptionsTail();
1387:   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1388:   TSGetSNES(ts,&snes);
1389:   if (!((PetscObject)snes)->type_name) {
1390:     SNESSetType(snes,SNESKSPONLY);
1391:   }
1392:   return(0);
1393: }

1395: static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
1396: {
1398:   PetscInt       i;
1399:   size_t         left,count;
1400:   char           *p;

1403:   for (i=0,p=buf,left=len; i<n; i++) {
1404:     PetscSNPrintfCount(p,left,fmt,&count,(double)x[i]);
1405:     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
1406:     left -= count;
1407:     p    += count;
1408:     *p++  = ' ';
1409:   }
1410:   p[i ? 0 : -1] = 0;
1411:   return(0);
1412: }

1414: static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
1415: {
1416:   TS_RosW        *ros = (TS_RosW*)ts->data;
1417:   PetscBool      iascii;

1421:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1422:   if (iascii) {
1423:     RosWTableau tab  = ros->tableau;
1424:     TSRosWType  rostype;
1425:     char        buf[512];
1426:     PetscInt    i;
1427:     PetscReal   abscissa[512];
1428:     TSRosWGetType(ts,&rostype);
1429:     PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);
1430:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ASum);
1431:     PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);
1432:     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1433:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,abscissa);
1434:     PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);
1435:   }
1436:   return(0);
1437: }

1439: static PetscErrorCode TSLoad_RosW(TS ts,PetscViewer viewer)
1440: {
1442:   SNES           snes;
1443:   TSAdapt        adapt;

1446:   TSGetAdapt(ts,&adapt);
1447:   TSAdaptLoad(adapt,viewer);
1448:   TSGetSNES(ts,&snes);
1449:   SNESLoad(snes,viewer);
1450:   /* function and Jacobian context for SNES when used with TS is always ts object */
1451:   SNESSetFunction(snes,NULL,NULL,ts);
1452:   SNESSetJacobian(snes,NULL,NULL,NULL,ts);
1453:   return(0);
1454: }

1456: /*@C
1457:   TSRosWSetType - Set the type of Rosenbrock-W scheme

1459:   Logically collective

1461:   Input Parameter:
1462: +  ts - timestepping context
1463: -  roswtype - type of Rosenbrock-W scheme

1465:   Level: beginner

1467: .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
1468: @*/
1469: PetscErrorCode TSRosWSetType(TS ts,TSRosWType roswtype)
1470: {

1476:   PetscTryMethod(ts,"TSRosWSetType_C",(TS,TSRosWType),(ts,roswtype));
1477:   return(0);
1478: }

1480: /*@C
1481:   TSRosWGetType - Get the type of Rosenbrock-W scheme

1483:   Logically collective

1485:   Input Parameter:
1486: .  ts - timestepping context

1488:   Output Parameter:
1489: .  rostype - type of Rosenbrock-W scheme

1491:   Level: intermediate

1493: .seealso: TSRosWGetType()
1494: @*/
1495: PetscErrorCode TSRosWGetType(TS ts,TSRosWType *rostype)
1496: {

1501:   PetscUseMethod(ts,"TSRosWGetType_C",(TS,TSRosWType*),(ts,rostype));
1502:   return(0);
1503: }

1505: /*@C
1506:   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.

1508:   Logically collective

1510:   Input Parameter:
1511: +  ts - timestepping context
1512: -  flg - PETSC_TRUE to recompute the Jacobian at each stage

1514:   Level: intermediate

1516: .seealso: TSRosWGetType()
1517: @*/
1518: PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1519: {

1524:   PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));
1525:   return(0);
1526: }

1528: static PetscErrorCode  TSRosWGetType_RosW(TS ts,TSRosWType *rostype)
1529: {
1530:   TS_RosW        *ros = (TS_RosW*)ts->data;

1533:   *rostype = ros->tableau->name;
1534:   return(0);
1535: }

1537: static PetscErrorCode  TSRosWSetType_RosW(TS ts,TSRosWType rostype)
1538: {
1539:   TS_RosW         *ros = (TS_RosW*)ts->data;
1540:   PetscErrorCode  ierr;
1541:   PetscBool       match;
1542:   RosWTableauLink link;

1545:   if (ros->tableau) {
1546:     PetscStrcmp(ros->tableau->name,rostype,&match);
1547:     if (match) return(0);
1548:   }
1549:   for (link = RosWTableauList; link; link=link->next) {
1550:     PetscStrcmp(link->tab.name,rostype,&match);
1551:     if (match) {
1552:       if (ts->setupcalled) {TSRosWTableauReset(ts);}
1553:       ros->tableau = &link->tab;
1554:       if (ts->setupcalled) {TSRosWTableauSetUp(ts);}
1555:       ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1556:       return(0);
1557:     }
1558:   }
1559:   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1560:   return(0);
1561: }

1563: static PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1564: {
1565:   TS_RosW *ros = (TS_RosW*)ts->data;

1568:   ros->recompute_jacobian = flg;
1569:   return(0);
1570: }

1572: static PetscErrorCode TSDestroy_RosW(TS ts)
1573: {

1577:   TSReset_RosW(ts);
1578:   if (ts->dm) {
1579:     DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRosW,DMRestrictHook_TSRosW,ts);
1580:     DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRosW,DMSubDomainRestrictHook_TSRosW,ts);
1581:   }
1582:   PetscFree(ts->data);
1583:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",NULL);
1584:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",NULL);
1585:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",NULL);
1586:   return(0);
1587: }

1589: /* ------------------------------------------------------------ */
1590: /*MC
1591:       TSROSW - ODE solver using Rosenbrock-W schemes

1593:   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1594:   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1595:   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().

1597:   Notes:
1598:   This method currently only works with autonomous ODE and DAE.

1600:   Consider trying TSARKIMEX if the stiff part is strongly nonlinear.

1602:   Developer notes:
1603:   Rosenbrock-W methods are typically specified for autonomous ODE

1605: $  udot = f(u)

1607:   by the stage equations

1609: $  k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j

1611:   and step completion formula

1613: $  u_1 = u_0 + sum_j b_j k_j

1615:   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u)
1616:   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1617:   we define new variables for the stage equations

1619: $  y_i = gamma_ij k_j

1621:   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define

1623: $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1}

1625:   to rewrite the method as

1627: $  [M/(h gamma_ii) - J] y_i = f(u_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1628: $  u_1 = u_0 + sum_j bt_j y_j

1630:    where we have introduced the mass matrix M. Continue by defining

1632: $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j

1634:    or, more compactly in tensor notation

1636: $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .

1638:    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1639:    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1640:    equation

1642: $  g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0

1644:    with initial guess y_i = 0.

1646:   Level: beginner

1648: .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister(), TSROSWTHETA1, TSROSWTHETA2, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3,
1649:            TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWGRK4T, TSROSWSHAMP4, TSROSWVELDD4, TSROSW4L
1650: M*/
1651: PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1652: {
1653:   TS_RosW        *ros;

1657:   TSRosWInitializePackage();

1659:   ts->ops->reset          = TSReset_RosW;
1660:   ts->ops->destroy        = TSDestroy_RosW;
1661:   ts->ops->view           = TSView_RosW;
1662:   ts->ops->load           = TSLoad_RosW;
1663:   ts->ops->setup          = TSSetUp_RosW;
1664:   ts->ops->step           = TSStep_RosW;
1665:   ts->ops->interpolate    = TSInterpolate_RosW;
1666:   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1667:   ts->ops->rollback       = TSRollBack_RosW;
1668:   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1669:   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1670:   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;

1672:   ts->usessnes = PETSC_TRUE;

1674:   PetscNewLog(ts,&ros);
1675:   ts->data = (void*)ros;

1677:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",TSRosWGetType_RosW);
1678:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",TSRosWSetType_RosW);
1679:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",TSRosWSetRecomputeJacobian_RosW);

1681:   TSRosWSetType(ts,TSRosWDefault);
1682:   return(0);
1683: }