Actual source code: rosw.c

  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: .  * - 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: .  * - 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: .  * - 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: .  * - 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: +   * -  Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979.
225: -   * -  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: +   * -  Shampine, Implementation of Rosenbrock methods, 1982.
245: -   * -  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: +   * -  van Veldhuizen, D stability and Kaps Rentrop methods, 1984.
265: -   * -  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: .  * -   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: .seealso:  TSRosWRegisterDestroy()
301: @*/
302: PetscErrorCode TSRosWRegisterAll(void)
303: {
304:   if (TSRosWRegisterAllCalled) return 0;
305:   TSRosWRegisterAllCalled = PETSC_TRUE;

307:   {
308:     const PetscReal A = 0;
309:     const PetscReal Gamma = 1;
310:     const PetscReal b = 1;
311:     const PetscReal binterpt=1;

313:     TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,NULL,1,&binterpt);
314:   }

316:   {
317:     const PetscReal A = 0;
318:     const PetscReal Gamma = 0.5;
319:     const PetscReal b = 1;
320:     const PetscReal binterpt=1;

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

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

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

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

368:       binterpt[0][0] = -0.8094010767585034;
369:       binterpt[1][0] = -0.5;
370:       binterpt[2][0] = 2.3094010767585034;
371:       binterpt[0][1] = 0.9641016151377548;
372:       binterpt[1][1] = 0.5;
373:       binterpt[2][1] = -1.4641016151377548;

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

392:     binterpt[0][0]=1.0564298455794094;
393:     binterpt[1][0]=2.296429974281067;
394:     binterpt[2][0]=-1.307599564525376;
395:     binterpt[3][0]=-1.045260255335102;
396:     binterpt[0][1]=-1.3864882699759573;
397:     binterpt[1][1]=-8.262611700275677;
398:     binterpt[2][1]=7.250979895056055;
399:     binterpt[3][1]=2.398120075195581;
400:     binterpt[0][2]=0.5721822314575016;
401:     binterpt[1][2]=4.742931142090097;
402:     binterpt[2][2]=-4.398120075195578;
403:     binterpt[3][2]=-0.9169932983520199;

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

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

435:     PetscReal binterpt[3][2];
436:     binterpt[0][0] = 3.793692883777660870425141387941;
437:     binterpt[1][0] = -2.918692883777660870425141387941;
438:     binterpt[2][0] = 0.125;
439:     binterpt[0][1] = -0.725741064379812106687651020584;
440:     binterpt[1][1] = 0.559074397713145440020984353917;
441:     binterpt[2][1] = 0.16666666666666666666666666666667;

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

461:     binterpt[0][0]=0.089316397477040902157517886164709;
462:     binterpt[1][0]=-0.91068360252295909784248211383529;
463:     binterpt[2][0]=1.8213672050459181956849642276706;
464:     binterpt[0][1]=0.077350269189625764509148780501957;
465:     binterpt[1][1]=1.077350269189625764509148780502;
466:     binterpt[2][1]=-1.1547005383792515290182975610039;

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

471:   {
472:     const PetscReal
473:       A[4][4] = {{0,0,0,0},
474:                  {1./2.,0,0,0},
475:                  {1./2.,1./2.,0,0},
476:                  {1./6.,1./6.,1./6.,0}},
477:       Gamma[4][4] = {{1./2.,0,0,0},
478:                      {0.0,1./4.,0,0},
479:                      {-2.,-2./3.,2./3.,0},
480:                      {1./2.,5./36.,-2./9,0}},
481:       b[4]  = {1./6.,1./6.,1./6.,1./2.},
482:       b2[4] = {1./8.,3./4.,1./8.,0};
483:     PetscReal binterpt[4][3];

485:     binterpt[0][0]=6.25;
486:     binterpt[1][0]=-30.25;
487:     binterpt[2][0]=1.75;
488:     binterpt[3][0]=23.25;
489:     binterpt[0][1]=-9.75;
490:     binterpt[1][1]=58.75;
491:     binterpt[2][1]=-3.25;
492:     binterpt[3][1]=-45.75;
493:     binterpt[0][2]=3.6666666666666666666666666666667;
494:     binterpt[1][2]=-28.333333333333333333333333333333;
495:     binterpt[2][2]=1.6666666666666666666666666666667;
496:     binterpt[3][2]=23.;

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

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

515:     binterpt[0][0]=1.6911764705882352941176470588235;
516:     binterpt[1][0]=3.6813725490196078431372549019608;
517:     binterpt[2][0]=0.23039215686274509803921568627451;
518:     binterpt[3][0]=-4.6029411764705882352941176470588;
519:     binterpt[0][1]=-0.95588235294117647058823529411765;
520:     binterpt[1][1]=-6.2401960784313725490196078431373;
521:     binterpt[2][1]=-0.31862745098039215686274509803922;
522:     binterpt[3][1]=7.5147058823529411764705882352941;
523:     binterpt[0][2]=-0.56862745098039215686274509803922;
524:     binterpt[1][2]=2.7254901960784313725490196078431;
525:     binterpt[2][2]=0.25490196078431372549019607843137;
526:     binterpt[3][2]=-2.4117647058823529411764705882353;

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

531:   {
532:     PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
533:     PetscReal binterpt[4][3];

535:     Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
536:     Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
537:     Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
538:     Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
539:     Gamma[1][2]=0; Gamma[1][3]=0;
540:     Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
541:     Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
542:     Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
543:     Gamma[2][3]=0;
544:     Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
545:     Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
546:     Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
547:     Gamma[3][3]=0;

549:     A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
550:     A[1][0]=0.8717330430169179988320388950590125027645343373957631;
551:     A[1][1]=0; A[1][2]=0; A[1][3]=0;
552:     A[2][0]=0.5275890119763004115618079766722914408876108660811028;
553:     A[2][1]=0.07241098802369958843819203208518599088698057726988732;
554:     A[2][2]=0; A[2][3]=0;
555:     A[3][0]=0.3990960076760701320627260685975778145384666450351314;
556:     A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
557:     A[3][2]=1.038461646937449311660120300601880176655352737312713;
558:     A[3][3]=0;

560:     b[0]=0.1876410243467238251612921333138006734899663569186926;
561:     b[1]=-0.5952974735769549480478230473706443582188442040780541;
562:     b[2]=0.9717899277217721234705114616271378792182450260943198;
563:     b[3]=0.4358665215084589994160194475295062513822671686978816;

565:     b2[0]=0.2147402862233891404862383521089097657790734483804460;
566:     b2[1]=-0.4851622638849390928209050538171743017757490232519684;
567:     b2[2]=0.8687250025203875511662123688667549217531982787600080;
568:     b2[3]=0.4016969751411624011684543450940068201770721128357014;

570:     binterpt[0][0]=2.2565812720167954547104627844105;
571:     binterpt[1][0]=1.349166413351089573796243820819;
572:     binterpt[2][0]=-2.4695174540533503758652847586647;
573:     binterpt[3][0]=-0.13623023131453465264142184656474;
574:     binterpt[0][1]=-3.0826699111559187902922463354557;
575:     binterpt[1][1]=-2.4689115685996042534544925650515;
576:     binterpt[2][1]=5.7428279814696677152129332773553;
577:     binterpt[3][1]=-0.19124650171414467146619437684812;
578:     binterpt[0][2]=1.0137296634858471607430756831148;
579:     binterpt[1][2]=0.52444768167155973161042570784064;
580:     binterpt[2][2]=-2.3015205996945452158771370439586;
581:     binterpt[3][2]=0.76334325453713832352363565300308;

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

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

595:    Not Collective

597:    Level: advanced

599: .seealso: TSRosWRegister(), TSRosWRegisterAll()
600: @*/
601: PetscErrorCode TSRosWRegisterDestroy(void)
602: {
603:   RosWTableauLink link;

605:   while ((link = RosWTableauList)) {
606:     RosWTableau t = &link->tab;
607:     RosWTableauList = link->next;
608:     PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);
609:     PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);
610:     PetscFree2(t->bembed,t->bembedt);
611:     PetscFree(t->binterpt);
612:     PetscFree(t->name);
613:     PetscFree(link);
614:   }
615:   TSRosWRegisterAllCalled = PETSC_FALSE;
616:   return 0;
617: }

619: /*@C
620:   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
621:   from TSInitializePackage().

623:   Level: developer

625: .seealso: PetscInitialize()
626: @*/
627: PetscErrorCode TSRosWInitializePackage(void)
628: {
629:   if (TSRosWPackageInitialized) return 0;
630:   TSRosWPackageInitialized = PETSC_TRUE;
631:   TSRosWRegisterAll();
632:   PetscRegisterFinalize(TSRosWFinalizePackage);
633:   return 0;
634: }

636: /*@C
637:   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
638:   called from PetscFinalize().

640:   Level: developer

642: .seealso: PetscFinalize()
643: @*/
644: PetscErrorCode TSRosWFinalizePackage(void)
645: {
646:   TSRosWPackageInitialized = PETSC_FALSE;
647:   TSRosWRegisterDestroy();
648:   return 0;
649: }

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

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

656:    Input Parameters:
657: +  name - identifier for method
658: .  order - approximation order of method
659: .  s - number of stages, this is the dimension of the matrices below
660: .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
661: .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
662: .  b - Step completion table (dimension s)
663: .  bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
664: .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
665: -  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)

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

670:    Level: advanced

672: .seealso: TSRosW
673: @*/
674: PetscErrorCode TSRosWRegister(TSRosWType name,PetscInt order,PetscInt s,const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[],
675:                               PetscInt pinterp,const PetscReal binterpt[])
676: {
677:   RosWTableauLink link;
678:   RosWTableau     t;
679:   PetscInt        i,j,k;
680:   PetscScalar     *GammaInv;


688:   TSRosWInitializePackage();
689:   PetscNew(&link);
690:   t        = &link->tab;
691:   PetscStrallocpy(name,&t->name);
692:   t->order = order;
693:   t->s     = s;
694:   PetscMalloc5(s*s,&t->A,s*s,&t->Gamma,s,&t->b,s,&t->ASum,s,&t->GammaSum);
695:   PetscMalloc5(s*s,&t->At,s,&t->bt,s*s,&t->GammaInv,s,&t->GammaZeroDiag,s*s,&t->GammaExplicitCorr);
696:   PetscArraycpy(t->A,A,s*s);
697:   PetscArraycpy(t->Gamma,Gamma,s*s);
698:   PetscArraycpy(t->GammaExplicitCorr,Gamma,s*s);
699:   PetscArraycpy(t->b,b,s);
700:   if (bembed) {
701:     PetscMalloc2(s,&t->bembed,s,&t->bembedt);
702:     PetscArraycpy(t->bembed,bembed,s);
703:   }
704:   for (i=0; i<s; i++) {
705:     t->ASum[i]     = 0;
706:     t->GammaSum[i] = 0;
707:     for (j=0; j<s; j++) {
708:       t->ASum[i]     += A[i*s+j];
709:       t->GammaSum[i] += Gamma[i*s+j];
710:     }
711:   }
712:   PetscMalloc1(s*s,&GammaInv); /* Need to use Scalar for inverse, then convert back to Real */
713:   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
714:   for (i=0; i<s; i++) {
715:     if (Gamma[i*s+i] == 0.0) {
716:       GammaInv[i*s+i] = 1.0;
717:       t->GammaZeroDiag[i] = PETSC_TRUE;
718:     } else {
719:       t->GammaZeroDiag[i] = PETSC_FALSE;
720:     }
721:   }

723:   switch (s) {
724:   case 1: GammaInv[0] = 1./GammaInv[0]; break;
725:   case 2: PetscKernel_A_gets_inverse_A_2(GammaInv,0,PETSC_FALSE,NULL); break;
726:   case 3: PetscKernel_A_gets_inverse_A_3(GammaInv,0,PETSC_FALSE,NULL); break;
727:   case 4: PetscKernel_A_gets_inverse_A_4(GammaInv,0,PETSC_FALSE,NULL); break;
728:   case 5: {
729:     PetscInt  ipvt5[5];
730:     MatScalar work5[5*5];
731:     PetscKernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0,PETSC_FALSE,NULL); break;
732:   }
733:   case 6: PetscKernel_A_gets_inverse_A_6(GammaInv,0,PETSC_FALSE,NULL); break;
734:   case 7: PetscKernel_A_gets_inverse_A_7(GammaInv,0,PETSC_FALSE,NULL); break;
735:   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
736:   }
737:   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
738:   PetscFree(GammaInv);

740:   for (i=0; i<s; i++) {
741:     for (k=0; k<i+1; k++) {
742:       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
743:       for (j=k+1; j<i+1; j++) {
744:         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
745:       }
746:     }
747:   }

749:   for (i=0; i<s; i++) {
750:     for (j=0; j<s; j++) {
751:       t->At[i*s+j] = 0;
752:       for (k=0; k<s; k++) {
753:         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
754:       }
755:     }
756:     t->bt[i] = 0;
757:     for (j=0; j<s; j++) {
758:       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
759:     }
760:     if (bembed) {
761:       t->bembedt[i] = 0;
762:       for (j=0; j<s; j++) {
763:         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
764:       }
765:     }
766:   }
767:   t->ccfl = 1.0;                /* Fix this */

769:   t->pinterp = pinterp;
770:   PetscMalloc1(s*pinterp,&t->binterpt);
771:   PetscArraycpy(t->binterpt,binterpt,s*pinterp);
772:   link->next = RosWTableauList;
773:   RosWTableauList = link;
774:   return 0;
775: }

777: /*@C
778:    TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices

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

782:    Input Parameters:
783: +  name - identifier for method
784: .  gamma - leading coefficient (diagonal entry)
785: .  a2 - design parameter, see Table 7.2 of Hairer&Wanner
786: .  a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
787: .  b3 - design parameter, see Table 7.2 of Hairer&Wanner
788: .  beta43 - design parameter or PETSC_DEFAULT to use Equation 7.21 of Hairer&Wanner
789: -  e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer

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

796:    Level: developer

798: .seealso: TSRosW, TSRosWRegister()
799: @*/
800: PetscErrorCode TSRosWRegisterRos4(TSRosWType name,PetscReal gamma,PetscReal a2,PetscReal a3,PetscReal b3,PetscReal e4)
801: {
802:   /* Declare numeric constants so they can be quad precision without being truncated at double */
803:   const PetscReal one = 1,two = 2,three = 3,four = 4,five = 5,six = 6,eight = 8,twelve = 12,twenty = 20,twentyfour = 24,
804:     p32 = one/six - gamma + gamma*gamma,
805:     p42 = one/eight - gamma/three,
806:     p43 = one/twelve - gamma/three,
807:     p44 = one/twentyfour - gamma/two + three/two*gamma*gamma - gamma*gamma*gamma,
808:     p56 = one/twenty - gamma/four;
809:   PetscReal   a4,a32,a42,a43,b1,b2,b4,beta2p,beta3p,beta4p,beta32,beta42,beta43,beta32beta2p,beta4jbetajp;
810:   PetscReal   A[4][4],Gamma[4][4],b[4],bm[4];
811:   PetscScalar M[3][3],rhs[3];

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

818:   /* Solve order conditions 7.15a, 7.15c, 7.15e */
819:   M[0][0] = one; M[0][1] = one;      M[0][2] = one;      /* 7.15a */
820:   M[1][0] = 0.0; M[1][1] = a2*a2;    M[1][2] = a4*a4;    /* 7.15c */
821:   M[2][0] = 0.0; M[2][1] = a2*a2*a2; M[2][2] = a4*a4*a4; /* 7.15e */
822:   rhs[0]  = one - b3;
823:   rhs[1]  = one/three - a3*a3*b3;
824:   rhs[2]  = one/four - a3*a3*a3*b3;
825:   PetscKernel_A_gets_inverse_A_3(&M[0][0],0,PETSC_FALSE,NULL);
826:   b1      = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
827:   b2      = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
828:   b4      = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);

830:   /* Step 3 */
831:   beta43       = (p56 - a2*p43) / (b4*a3*a3*(a3 - a2)); /* 7.21 */
832:   beta32beta2p =  p44 / (b4*beta43);                    /* 7.15h */
833:   beta4jbetajp = (p32 - b3*beta32beta2p) / b4;
834:   M[0][0]      = b2;                                    M[0][1] = b3;                 M[0][2] = b4;
835:   M[1][0]      = a4*a4*beta32beta2p-a3*a3*beta4jbetajp; M[1][1] = a2*a2*beta4jbetajp; M[1][2] = -a2*a2*beta32beta2p;
836:   M[2][0]      = b4*beta43*a3*a3-p43;                   M[2][1] = -b4*beta43*a2*a2;   M[2][2] = 0;
837:   rhs[0]       = one/two - gamma; rhs[1] = 0; rhs[2] = -a2*a2*p32;
838:   PetscKernel_A_gets_inverse_A_3(&M[0][0],0,PETSC_FALSE,NULL);
839:   beta2p       = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
840:   beta3p       = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
841:   beta4p       = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);

843:   /* Step 4: back-substitute */
844:   beta32 = beta32beta2p / beta2p;
845:   beta42 = (beta4jbetajp - beta43*beta3p) / beta2p;

847:   /* Step 5: 7.15f and 7.20, then 7.16 */
848:   a43 = 0;
849:   a32 = p42 / (b3*a3*beta2p + b4*a4*beta2p);
850:   a42 = a32;

852:   A[0][0]     = 0;          A[0][1] = 0;   A[0][2] = 0;   A[0][3] = 0;
853:   A[1][0]     = a2;         A[1][1] = 0;   A[1][2] = 0;   A[1][3] = 0;
854:   A[2][0]     = a3-a32;     A[2][1] = a32; A[2][2] = 0;   A[2][3] = 0;
855:   A[3][0]     = a4-a43-a42; A[3][1] = a42; A[3][2] = a43; A[3][3] = 0;
856:   Gamma[0][0] = gamma;                        Gamma[0][1] = 0;              Gamma[0][2] = 0;              Gamma[0][3] = 0;
857:   Gamma[1][0] = beta2p-A[1][0];               Gamma[1][1] = gamma;          Gamma[1][2] = 0;              Gamma[1][3] = 0;
858:   Gamma[2][0] = beta3p-beta32-A[2][0];        Gamma[2][1] = beta32-A[2][1]; Gamma[2][2] = gamma;          Gamma[2][3] = 0;
859:   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;
860:   b[0] = b1; b[1] = b2; b[2] = b3; b[3] = b4;

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

868:   {
869:     const PetscReal misfit = a2*a2*bm[1] + a3*a3*bm[2] + a4*a4*bm[3] - one/three;
871:   }
872:   TSRosWRegister(name,4,4,&A[0][0],&Gamma[0][0],b,bm,0,NULL);
873:   return 0;
874: }

876: /*
877:  The step completion formula is

879:  x1 = x0 + b^T Y

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

884:  x1e = x0 + be^T Y
885:      = x1 - b^T Y + be^T Y
886:      = x1 + (be - b)^T Y

888:  so we can evaluate the method of different order even after the step has been optimistically completed.
889: */
890: static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec U,PetscBool *done)
891: {
892:   TS_RosW        *ros = (TS_RosW*)ts->data;
893:   RosWTableau    tab  = ros->tableau;
894:   PetscScalar    *w   = ros->work;
895:   PetscInt       i;

897:   if (order == tab->order) {
898:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
899:       VecCopy(ts->vec_sol,U);
900:       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
901:       VecMAXPY(U,tab->s,w,ros->Y);
902:     } else VecCopy(ts->vec_sol,U);
903:     if (done) *done = PETSC_TRUE;
904:     return 0;
905:   } else if (order == tab->order-1) {
906:     if (!tab->bembedt) goto unavailable;
907:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
908:       VecCopy(ts->vec_sol,U);
909:       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
910:       VecMAXPY(U,tab->s,w,ros->Y);
911:     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
912:       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
913:       VecCopy(ts->vec_sol,U);
914:       VecMAXPY(U,tab->s,w,ros->Y);
915:     }
916:     if (done) *done = PETSC_TRUE;
917:     return 0;
918:   }
919:   unavailable:
920:   if (done) *done = PETSC_FALSE;
921:   else SETERRQ(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);
922:   return 0;
923: }

925: static PetscErrorCode TSRollBack_RosW(TS ts)
926: {
927:   TS_RosW        *ros = (TS_RosW*)ts->data;

929:   VecCopy(ros->vec_sol_prev,ts->vec_sol);
930:   return 0;
931: }

933: static PetscErrorCode TSStep_RosW(TS ts)
934: {
935:   TS_RosW         *ros = (TS_RosW*)ts->data;
936:   RosWTableau     tab  = ros->tableau;
937:   const PetscInt  s    = tab->s;
938:   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
939:   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
940:   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
941:   PetscScalar     *w   = ros->work;
942:   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
943:   SNES            snes;
944:   TSAdapt         adapt;
945:   PetscInt        i,j,its,lits;
946:   PetscInt        rejections = 0;
947:   PetscBool       stageok,accept = PETSC_TRUE;
948:   PetscReal       next_time_step = ts->time_step;
949:   PetscInt        lag;

951:   if (!ts->steprollback) {
952:     VecCopy(ts->vec_sol,ros->vec_sol_prev);
953:   }

955:   ros->status = TS_STEP_INCOMPLETE;
956:   while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
957:     const PetscReal h = ts->time_step;
958:     for (i=0; i<s; i++) {
959:       ros->stage_time = ts->ptime + h*ASum[i];
960:       TSPreStage(ts,ros->stage_time);
961:       if (GammaZeroDiag[i]) {
962:         ros->stage_explicit = PETSC_TRUE;
963:         ros->scoeff         = 1.;
964:       } else {
965:         ros->stage_explicit = PETSC_FALSE;
966:         ros->scoeff         = 1./Gamma[i*s+i];
967:       }

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

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

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

980:       if (!ros->stage_explicit) {
981:         TSGetSNES(ts,&snes);
982:         if (!ros->recompute_jacobian && !i) {
983:           SNESGetLagJacobian(snes,&lag);
984:           if (lag == 1) {  /* use did not set a nontrivial lag, so lag over all stages */
985:             SNESSetLagJacobian(snes,-2); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */
986:           }
987:         }
988:         SNESSolve(snes,NULL,Y[i]);
989:         if (!ros->recompute_jacobian && i == s-1 && lag == 1) {
990:           SNESSetLagJacobian(snes,lag); /* Set lag back to 1 so we know user did not set it */
991:         }
992:         SNESGetIterationNumber(snes,&its);
993:         SNESGetLinearSolveIterations(snes,&lits);
994:         ts->snes_its += its; ts->ksp_its += lits;
995:       } else {
996:         Mat J,Jp;
997:         VecZeroEntries(Ydot); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
998:         TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);
999:         VecScale(Y[i],-1.0);
1000:         VecAXPY(Y[i],-1.0,Zdot); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/

1002:         VecZeroEntries(Zstage); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
1003:         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
1004:         VecMAXPY(Zstage,i,w,Y);

1006:         /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
1007:         TSGetIJacobian(ts,&J,&Jp,NULL,NULL);
1008:         TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,J,Jp,PETSC_FALSE);
1009:         MatMult(J,Zstage,Zdot);
1010:         VecAXPY(Y[i],-1.0,Zdot);
1011:         ts->ksp_its += 1;

1013:         VecScale(Y[i],h);
1014:       }
1015:       TSPostStage(ts,ros->stage_time,i,Y);
1016:       TSGetAdapt(ts,&adapt);
1017:       TSAdaptCheckStage(adapt,ts,ros->stage_time,Y[i],&stageok);
1018:       if (!stageok) goto reject_step;
1019:     }

1021:     ros->status = TS_STEP_INCOMPLETE;
1022:     TSEvaluateStep_RosW(ts,tab->order,ts->vec_sol,NULL);
1023:     ros->status = TS_STEP_PENDING;
1024:     TSGetAdapt(ts,&adapt);
1025:     TSAdaptCandidatesClear(adapt);
1026:     TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);
1027:     TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);
1028:     ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1029:     if (!accept) { /* Roll back the current step */
1030:       TSRollBack_RosW(ts);
1031:       ts->time_step = next_time_step;
1032:       goto reject_step;
1033:     }

1035:     ts->ptime += ts->time_step;
1036:     ts->time_step = next_time_step;
1037:     break;

1039:   reject_step:
1040:     ts->reject++; accept = PETSC_FALSE;
1041:     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1042:       ts->reason = TS_DIVERGED_STEP_REJECTED;
1043:       PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);
1044:     }
1045:   }
1046:   return 0;
1047: }

1049: static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec U)
1050: {
1051:   TS_RosW         *ros = (TS_RosW*)ts->data;
1052:   PetscInt        s    = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j;
1053:   PetscReal       h;
1054:   PetscReal       tt,t;
1055:   PetscScalar     *bt;
1056:   const PetscReal *Bt = ros->tableau->binterpt;
1057:   const PetscReal *GammaInv = ros->tableau->GammaInv;
1058:   PetscScalar     *w        = ros->work;
1059:   Vec             *Y        = ros->Y;


1063:   switch (ros->status) {
1064:   case TS_STEP_INCOMPLETE:
1065:   case TS_STEP_PENDING:
1066:     h = ts->time_step;
1067:     t = (itime - ts->ptime)/h;
1068:     break;
1069:   case TS_STEP_COMPLETE:
1070:     h = ts->ptime - ts->ptime_prev;
1071:     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
1072:     break;
1073:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
1074:   }
1075:   PetscMalloc1(s,&bt);
1076:   for (i=0; i<s; i++) bt[i] = 0;
1077:   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
1078:     for (i=0; i<s; i++) {
1079:       bt[i] += Bt[i*pinterp+j] * tt;
1080:     }
1081:   }

1083:   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1084:   /* U <- 0*/
1085:   VecZeroEntries(U);
1086:   /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
1087:   for (j=0; j<s; j++) w[j] = 0;
1088:   for (j=0; j<s; j++) {
1089:     for (i=j; i<s; i++) {
1090:       w[j] +=  bt[i]*GammaInv[i*s+j];
1091:     }
1092:   }
1093:   VecMAXPY(U,i,w,Y);
1094:   /* U <- y(t) + U */
1095:   VecAXPY(U,1,ros->vec_sol_prev);

1097:   PetscFree(bt);
1098:   return 0;
1099: }

1101: /*------------------------------------------------------------*/

1103: static PetscErrorCode TSRosWTableauReset(TS ts)
1104: {
1105:   TS_RosW        *ros = (TS_RosW*)ts->data;
1106:   RosWTableau    tab  = ros->tableau;

1108:   if (!tab) return 0;
1109:   VecDestroyVecs(tab->s,&ros->Y);
1110:   PetscFree(ros->work);
1111:   return 0;
1112: }

1114: static PetscErrorCode TSReset_RosW(TS ts)
1115: {
1116:   TS_RosW        *ros = (TS_RosW*)ts->data;

1118:   TSRosWTableauReset(ts);
1119:   VecDestroy(&ros->Ydot);
1120:   VecDestroy(&ros->Ystage);
1121:   VecDestroy(&ros->Zdot);
1122:   VecDestroy(&ros->Zstage);
1123:   VecDestroy(&ros->vec_sol_prev);
1124:   return 0;
1125: }

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

1131:   if (Ydot) {
1132:     if (dm && dm != ts->dm) {
1133:       DMGetNamedGlobalVector(dm,"TSRosW_Ydot",Ydot);
1134:     } else *Ydot = rw->Ydot;
1135:   }
1136:   if (Zdot) {
1137:     if (dm && dm != ts->dm) {
1138:       DMGetNamedGlobalVector(dm,"TSRosW_Zdot",Zdot);
1139:     } else *Zdot = rw->Zdot;
1140:   }
1141:   if (Ystage) {
1142:     if (dm && dm != ts->dm) {
1143:       DMGetNamedGlobalVector(dm,"TSRosW_Ystage",Ystage);
1144:     } else *Ystage = rw->Ystage;
1145:   }
1146:   if (Zstage) {
1147:     if (dm && dm != ts->dm) {
1148:       DMGetNamedGlobalVector(dm,"TSRosW_Zstage",Zstage);
1149:     } else *Zstage = rw->Zstage;
1150:   }
1151:   return 0;
1152: }

1154: static PetscErrorCode TSRosWRestoreVecs(TS ts,DM dm,Vec *Ydot,Vec *Zdot, Vec *Ystage, Vec *Zstage)
1155: {
1156:   if (Ydot) {
1157:     if (dm && dm != ts->dm) {
1158:       DMRestoreNamedGlobalVector(dm,"TSRosW_Ydot",Ydot);
1159:     }
1160:   }
1161:   if (Zdot) {
1162:     if (dm && dm != ts->dm) {
1163:       DMRestoreNamedGlobalVector(dm,"TSRosW_Zdot",Zdot);
1164:     }
1165:   }
1166:   if (Ystage) {
1167:     if (dm && dm != ts->dm) {
1168:       DMRestoreNamedGlobalVector(dm,"TSRosW_Ystage",Ystage);
1169:     }
1170:   }
1171:   if (Zstage) {
1172:     if (dm && dm != ts->dm) {
1173:       DMRestoreNamedGlobalVector(dm,"TSRosW_Zstage",Zstage);
1174:     }
1175:   }
1176:   return 0;
1177: }

1179: static PetscErrorCode DMCoarsenHook_TSRosW(DM fine,DM coarse,void *ctx)
1180: {
1181:   return 0;
1182: }

1184: static PetscErrorCode DMRestrictHook_TSRosW(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1185: {
1186:   TS             ts = (TS)ctx;
1187:   Vec            Ydot,Zdot,Ystage,Zstage;
1188:   Vec            Ydotc,Zdotc,Ystagec,Zstagec;

1190:   TSRosWGetVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage);
1191:   TSRosWGetVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec);
1192:   MatRestrict(restrct,Ydot,Ydotc);
1193:   VecPointwiseMult(Ydotc,rscale,Ydotc);
1194:   MatRestrict(restrct,Ystage,Ystagec);
1195:   VecPointwiseMult(Ystagec,rscale,Ystagec);
1196:   MatRestrict(restrct,Zdot,Zdotc);
1197:   VecPointwiseMult(Zdotc,rscale,Zdotc);
1198:   MatRestrict(restrct,Zstage,Zstagec);
1199:   VecPointwiseMult(Zstagec,rscale,Zstagec);
1200:   TSRosWRestoreVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage);
1201:   TSRosWRestoreVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec);
1202:   return 0;
1203: }

1205: static PetscErrorCode DMSubDomainHook_TSRosW(DM fine,DM coarse,void *ctx)
1206: {
1207:   return 0;
1208: }

1210: static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1211: {
1212:   TS             ts = (TS)ctx;
1213:   Vec            Ydot,Zdot,Ystage,Zstage;
1214:   Vec            Ydots,Zdots,Ystages,Zstages;

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

1219:   VecScatterBegin(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD);
1220:   VecScatterEnd(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD);

1222:   VecScatterBegin(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD);
1223:   VecScatterEnd(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD);

1225:   VecScatterBegin(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD);
1226:   VecScatterEnd(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD);

1228:   VecScatterBegin(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD);
1229:   VecScatterEnd(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD);

1231:   TSRosWRestoreVecs(ts,dm,&Ydot,&Ystage,&Zdot,&Zstage);
1232:   TSRosWRestoreVecs(ts,subdm,&Ydots,&Ystages,&Zdots,&Zstages);
1233:   return 0;
1234: }

1236: /*
1237:   This defines the nonlinear equation that is to be solved with SNES
1238:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1239: */
1240: static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec U,Vec F,TS ts)
1241: {
1242:   TS_RosW        *ros = (TS_RosW*)ts->data;
1243:   Vec            Ydot,Zdot,Ystage,Zstage;
1244:   PetscReal      shift = ros->scoeff / ts->time_step;
1245:   DM             dm,dmsave;

1247:   SNESGetDM(snes,&dm);
1248:   TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1249:   VecWAXPY(Ydot,shift,U,Zdot);    /* Ydot = shift*U + Zdot */
1250:   VecWAXPY(Ystage,1.0,U,Zstage);  /* Ystage = U + Zstage */
1251:   dmsave = ts->dm;
1252:   ts->dm = dm;
1253:   TSComputeIFunction(ts,ros->stage_time,Ystage,Ydot,F,PETSC_FALSE);
1254:   ts->dm = dmsave;
1255:   TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1256:   return 0;
1257: }

1259: static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec U,Mat A,Mat B,TS ts)
1260: {
1261:   TS_RosW        *ros = (TS_RosW*)ts->data;
1262:   Vec            Ydot,Zdot,Ystage,Zstage;
1263:   PetscReal      shift = ros->scoeff / ts->time_step;
1264:   DM             dm,dmsave;

1266:   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
1267:   SNESGetDM(snes,&dm);
1268:   TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1269:   dmsave = ts->dm;
1270:   ts->dm = dm;
1271:   TSComputeIJacobian(ts,ros->stage_time,Ystage,Ydot,shift,A,B,PETSC_TRUE);
1272:   ts->dm = dmsave;
1273:   TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1274:   return 0;
1275: }

1277: static PetscErrorCode TSRosWTableauSetUp(TS ts)
1278: {
1279:   TS_RosW        *ros = (TS_RosW*)ts->data;
1280:   RosWTableau    tab  = ros->tableau;

1282:   VecDuplicateVecs(ts->vec_sol,tab->s,&ros->Y);
1283:   PetscMalloc1(tab->s,&ros->work);
1284:   return 0;
1285: }

1287: static PetscErrorCode TSSetUp_RosW(TS ts)
1288: {
1289:   TS_RosW        *ros = (TS_RosW*)ts->data;
1290:   DM             dm;
1291:   SNES           snes;
1292:   TSRHSJacobian  rhsjacobian;

1294:   TSRosWTableauSetUp(ts);
1295:   VecDuplicate(ts->vec_sol,&ros->Ydot);
1296:   VecDuplicate(ts->vec_sol,&ros->Ystage);
1297:   VecDuplicate(ts->vec_sol,&ros->Zdot);
1298:   VecDuplicate(ts->vec_sol,&ros->Zstage);
1299:   VecDuplicate(ts->vec_sol,&ros->vec_sol_prev);
1300:   TSGetDM(ts,&dm);
1301:   DMCoarsenHookAdd(dm,DMCoarsenHook_TSRosW,DMRestrictHook_TSRosW,ts);
1302:   DMSubDomainHookAdd(dm,DMSubDomainHook_TSRosW,DMSubDomainRestrictHook_TSRosW,ts);
1303:   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1304:   TSGetSNES(ts,&snes);
1305:   if (!((PetscObject)snes)->type_name) {
1306:     SNESSetType(snes,SNESKSPONLY);
1307:   }
1308:   DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);
1309:   if (rhsjacobian == TSComputeRHSJacobianConstant) {
1310:     Mat Amat,Pmat;

1312:     /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */
1313:     SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL);
1314:     if (Amat && Amat == ts->Arhs) {
1315:       if (Amat == Pmat) {
1316:         MatDuplicate(ts->Arhs,MAT_COPY_VALUES,&Amat);
1317:         SNESSetJacobian(snes,Amat,Amat,NULL,NULL);
1318:       } else {
1319:         MatDuplicate(ts->Arhs,MAT_COPY_VALUES,&Amat);
1320:         SNESSetJacobian(snes,Amat,NULL,NULL,NULL);
1321:         if (Pmat && Pmat == ts->Brhs) {
1322:           MatDuplicate(ts->Brhs,MAT_COPY_VALUES,&Pmat);
1323:           SNESSetJacobian(snes,NULL,Pmat,NULL,NULL);
1324:           MatDestroy(&Pmat);
1325:         }
1326:       }
1327:       MatDestroy(&Amat);
1328:     }
1329:   }
1330:   return 0;
1331: }
1332: /*------------------------------------------------------------*/

1334: static PetscErrorCode TSSetFromOptions_RosW(PetscOptionItems *PetscOptionsObject,TS ts)
1335: {
1336:   TS_RosW        *ros = (TS_RosW*)ts->data;
1337:   SNES           snes;

1339:   PetscOptionsHead(PetscOptionsObject,"RosW ODE solver options");
1340:   {
1341:     RosWTableauLink link;
1342:     PetscInt        count,choice;
1343:     PetscBool       flg;
1344:     const char      **namelist;

1346:     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
1347:     PetscMalloc1(count,(char***)&namelist);
1348:     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1349:     PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,ros->tableau->name,&choice,&flg);
1350:     if (flg) TSRosWSetType(ts,namelist[choice]);
1351:     PetscFree(namelist);

1353:     PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,NULL);
1354:   }
1355:   PetscOptionsTail();
1356:   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1357:   TSGetSNES(ts,&snes);
1358:   if (!((PetscObject)snes)->type_name) {
1359:     SNESSetType(snes,SNESKSPONLY);
1360:   }
1361:   return 0;
1362: }

1364: static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
1365: {
1366:   TS_RosW        *ros = (TS_RosW*)ts->data;
1367:   PetscBool      iascii;

1369:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1370:   if (iascii) {
1371:     RosWTableau tab  = ros->tableau;
1372:     TSRosWType  rostype;
1373:     char        buf[512];
1374:     PetscInt    i;
1375:     PetscReal   abscissa[512];
1376:     TSRosWGetType(ts,&rostype);
1377:     PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);
1378:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ASum);
1379:     PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);
1380:     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1381:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,abscissa);
1382:     PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);
1383:   }
1384:   return 0;
1385: }

1387: static PetscErrorCode TSLoad_RosW(TS ts,PetscViewer viewer)
1388: {
1389:   SNES           snes;
1390:   TSAdapt        adapt;

1392:   TSGetAdapt(ts,&adapt);
1393:   TSAdaptLoad(adapt,viewer);
1394:   TSGetSNES(ts,&snes);
1395:   SNESLoad(snes,viewer);
1396:   /* function and Jacobian context for SNES when used with TS is always ts object */
1397:   SNESSetFunction(snes,NULL,NULL,ts);
1398:   SNESSetJacobian(snes,NULL,NULL,NULL,ts);
1399:   return 0;
1400: }

1402: /*@C
1403:   TSRosWSetType - Set the type of Rosenbrock-W scheme

1405:   Logically collective

1407:   Input Parameters:
1408: +  ts - timestepping context
1409: -  roswtype - type of Rosenbrock-W scheme

1411:   Level: beginner

1413: .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
1414: @*/
1415: PetscErrorCode TSRosWSetType(TS ts,TSRosWType roswtype)
1416: {
1419:   PetscTryMethod(ts,"TSRosWSetType_C",(TS,TSRosWType),(ts,roswtype));
1420:   return 0;
1421: }

1423: /*@C
1424:   TSRosWGetType - Get the type of Rosenbrock-W scheme

1426:   Logically collective

1428:   Input Parameter:
1429: .  ts - timestepping context

1431:   Output Parameter:
1432: .  rostype - type of Rosenbrock-W scheme

1434:   Level: intermediate

1436: .seealso: TSRosWGetType()
1437: @*/
1438: PetscErrorCode TSRosWGetType(TS ts,TSRosWType *rostype)
1439: {
1441:   PetscUseMethod(ts,"TSRosWGetType_C",(TS,TSRosWType*),(ts,rostype));
1442:   return 0;
1443: }

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

1448:   Logically collective

1450:   Input Parameters:
1451: +  ts - timestepping context
1452: -  flg - PETSC_TRUE to recompute the Jacobian at each stage

1454:   Level: intermediate

1456: .seealso: TSRosWGetType()
1457: @*/
1458: PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1459: {
1461:   PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));
1462:   return 0;
1463: }

1465: static PetscErrorCode  TSRosWGetType_RosW(TS ts,TSRosWType *rostype)
1466: {
1467:   TS_RosW        *ros = (TS_RosW*)ts->data;

1469:   *rostype = ros->tableau->name;
1470:   return 0;
1471: }

1473: static PetscErrorCode  TSRosWSetType_RosW(TS ts,TSRosWType rostype)
1474: {
1475:   TS_RosW         *ros = (TS_RosW*)ts->data;
1476:   PetscBool       match;
1477:   RosWTableauLink link;

1479:   if (ros->tableau) {
1480:     PetscStrcmp(ros->tableau->name,rostype,&match);
1481:     if (match) return 0;
1482:   }
1483:   for (link = RosWTableauList; link; link=link->next) {
1484:     PetscStrcmp(link->tab.name,rostype,&match);
1485:     if (match) {
1486:       if (ts->setupcalled) TSRosWTableauReset(ts);
1487:       ros->tableau = &link->tab;
1488:       if (ts->setupcalled) TSRosWTableauSetUp(ts);
1489:       ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1490:       return 0;
1491:     }
1492:   }
1493:   SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1494: }

1496: static PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1497: {
1498:   TS_RosW *ros = (TS_RosW*)ts->data;

1500:   ros->recompute_jacobian = flg;
1501:   return 0;
1502: }

1504: static PetscErrorCode TSDestroy_RosW(TS ts)
1505: {
1506:   TSReset_RosW(ts);
1507:   if (ts->dm) {
1508:     DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRosW,DMRestrictHook_TSRosW,ts);
1509:     DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRosW,DMSubDomainRestrictHook_TSRosW,ts);
1510:   }
1511:   PetscFree(ts->data);
1512:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",NULL);
1513:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",NULL);
1514:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",NULL);
1515:   return 0;
1516: }

1518: /* ------------------------------------------------------------ */
1519: /*MC
1520:       TSROSW - ODE solver using Rosenbrock-W schemes

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

1526:   Notes:
1527:   This method currently only works with autonomous ODE and DAE.

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

1531:   Since this uses a single linear solve per time-step if you wish to lag the jacobian or preconditioner computation you must use also -snes_lag_jacobian_persists true or -snes_lag_jacobian_preconditioner true

1533:   Developer Notes:
1534:   Rosenbrock-W methods are typically specified for autonomous ODE

1536: $  udot = f(u)

1538:   by the stage equations

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

1542:   and step completion formula

1544: $  u_1 = u_0 + sum_j b_j k_j

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

1550: $  y_i = gamma_ij k_j

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

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

1556:   to rewrite the method as

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

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

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

1565:    or, more compactly in tensor notation

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

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

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

1575:    with initial guess y_i = 0.

1577:   Level: beginner

1579: .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister(), TSROSWTHETA1, TSROSWTHETA2, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3,
1580:            TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWGRK4T, TSROSWSHAMP4, TSROSWVELDD4, TSROSW4L
1581: M*/
1582: PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1583: {
1584:   TS_RosW        *ros;

1586:   TSRosWInitializePackage();

1588:   ts->ops->reset          = TSReset_RosW;
1589:   ts->ops->destroy        = TSDestroy_RosW;
1590:   ts->ops->view           = TSView_RosW;
1591:   ts->ops->load           = TSLoad_RosW;
1592:   ts->ops->setup          = TSSetUp_RosW;
1593:   ts->ops->step           = TSStep_RosW;
1594:   ts->ops->interpolate    = TSInterpolate_RosW;
1595:   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1596:   ts->ops->rollback       = TSRollBack_RosW;
1597:   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1598:   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1599:   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;

1601:   ts->usessnes = PETSC_TRUE;

1603:   PetscNewLog(ts,&ros);
1604:   ts->data = (void*)ros;

1606:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",TSRosWGetType_RosW);
1607:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",TSRosWSetType_RosW);
1608:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",TSRosWSetRecomputeJacobian_RosW);

1610:   TSRosWSetType(ts,TSRosWDefault);
1611:   return 0;
1612: }