Actual source code: rosw.c

petsc-3.4.5 2014-06-29
  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>                /*I   "petscts.h"   I*/
 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          VecSolPrev;       /* Work vector holding the solution from the previous step (used for interpolation)*/
 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.

226:      Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

228:      Hairer's code ros4.f

230:      Level: intermediate

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

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

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

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

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

244:      References:
245:      Shampine, Implementation of Rosenbrock methods, 1982.

247:      Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

249:      Hairer's code ros4.f

251:      Level: intermediate

253: .seealso: TSROSW, TSROSWGRK4T, TSROSWVELDD4, TSROSW4L
254: M*/

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

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

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

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

265:      References:
266:      van Veldhuizen, D-stability and Kaps-Rentrop methods, 1984.

268:      Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

270:      Hairer's code ros4.f

272:      Level: intermediate

274: .seealso: TSROSW, TSROSWGRK4T, TSROSWSHAMP4, TSROSW4L
275: M*/

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

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

282:      A-stable and L-stable

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

286:      References:
287:      Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

289:      Hairer's code ros4.f

291:      Level: intermediate

293: .seealso: TSROSW, TSROSWGRK4T, TSROSWSHAMP4, TSROSW4L
294: M*/

298: /*@C
299:   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW

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

303:   Level: advanced

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

307: .seealso:  TSRosWRegisterDestroy()
308: @*/
309: PetscErrorCode TSRosWRegisterAll(void)
310: {

314:   if (TSRosWRegisterAllCalled) return(0);
315:   TSRosWRegisterAllCalled = PETSC_TRUE;

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

323:     TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,NULL,1,&binterpt);
324:   }

326:   {
327:     const PetscReal A = 0;
328:     const PetscReal Gamma = 0.5;
329:     const PetscReal b = 1;
330:     const PetscReal binterpt=1;

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

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

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

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

378:       binterpt[0][0] = -0.8094010767585034;
379:       binterpt[1][0] = -0.5;
380:       binterpt[2][0] = 2.3094010767585034;
381:       binterpt[0][1] = 0.9641016151377548;
382:       binterpt[1][1] = 0.5;
383:       binterpt[2][1] = -1.4641016151377548;

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

402:     binterpt[0][0]=1.0564298455794094;
403:     binterpt[1][0]=2.296429974281067;
404:     binterpt[2][0]=-1.307599564525376;
405:     binterpt[3][0]=-1.045260255335102;
406:     binterpt[0][1]=-1.3864882699759573;
407:     binterpt[1][1]=-8.262611700275677;
408:     binterpt[2][1]=7.250979895056055;
409:     binterpt[3][1]=2.398120075195581;
410:     binterpt[0][2]=0.5721822314575016;
411:     binterpt[1][2]=4.742931142090097;
412:     binterpt[2][2]=-4.398120075195578;
413:     binterpt[3][2]=-0.9169932983520199;

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

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

445:     PetscReal binterpt[3][2];
446:     binterpt[0][0] = 3.793692883777660870425141387941;
447:     binterpt[1][0] = -2.918692883777660870425141387941;
448:     binterpt[2][0] = 0.125;
449:     binterpt[0][1] = -0.725741064379812106687651020584;
450:     binterpt[1][1] = 0.559074397713145440020984353917;
451:     binterpt[2][1] = 0.16666666666666666666666666666667;

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

471:     binterpt[0][0]=0.089316397477040902157517886164709;
472:     binterpt[1][0]=-0.91068360252295909784248211383529;
473:     binterpt[2][0]=1.8213672050459181956849642276706;
474:     binterpt[0][1]=0.077350269189625764509148780501957;
475:     binterpt[1][1]=1.077350269189625764509148780502;
476:     binterpt[2][1]=-1.1547005383792515290182975610039;

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

481:   {
482:     const PetscReal
483:       A[4][4] = {{0,0,0,0},
484:                  {1./2.,0,0,0},
485:                  {1./2.,1./2.,0,0},
486:                  {1./6.,1./6.,1./6.,0}},
487:       Gamma[4][4] = {{1./2.,0,0,0},
488:                      {0.0,1./4.,0,0},
489:                      {-2.,-2./3.,2./3.,0},
490:                      {1./2.,5./36.,-2./9,0}},
491:       b[4]  = {1./6.,1./6.,1./6.,1./2.},
492:       b2[4] = {1./8.,3./4.,1./8.,0};
493:     PetscReal binterpt[4][3];

495:     binterpt[0][0]=6.25;
496:     binterpt[1][0]=-30.25;
497:     binterpt[2][0]=1.75;
498:     binterpt[3][0]=23.25;
499:     binterpt[0][1]=-9.75;
500:     binterpt[1][1]=58.75;
501:     binterpt[2][1]=-3.25;
502:     binterpt[3][1]=-45.75;
503:     binterpt[0][2]=3.6666666666666666666666666666667;
504:     binterpt[1][2]=-28.333333333333333333333333333333;
505:     binterpt[2][2]=1.6666666666666666666666666666667;
506:     binterpt[3][2]=23.;

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

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

525:     binterpt[0][0]=1.6911764705882352941176470588235;
526:     binterpt[1][0]=3.6813725490196078431372549019608;
527:     binterpt[2][0]=0.23039215686274509803921568627451;
528:     binterpt[3][0]=-4.6029411764705882352941176470588;
529:     binterpt[0][1]=-0.95588235294117647058823529411765;
530:     binterpt[1][1]=-6.2401960784313725490196078431373;
531:     binterpt[2][1]=-0.31862745098039215686274509803922;
532:     binterpt[3][1]=7.5147058823529411764705882352941;
533:     binterpt[0][2]=-0.56862745098039215686274509803922;
534:     binterpt[1][2]=2.7254901960784313725490196078431;
535:     binterpt[2][2]=0.25490196078431372549019607843137;
536:     binterpt[3][2]=-2.4117647058823529411764705882353;

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

541:   {
542:     PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
543:     PetscReal binterpt[4][3];

545:     Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
546:     Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
547:     Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
548:     Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
549:     Gamma[1][2]=0; Gamma[1][3]=0;
550:     Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
551:     Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
552:     Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
553:     Gamma[2][3]=0;
554:     Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
555:     Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
556:     Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
557:     Gamma[3][3]=0;

559:     A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
560:     A[1][0]=0.8717330430169179988320388950590125027645343373957631;
561:     A[1][1]=0; A[1][2]=0; A[1][3]=0;
562:     A[2][0]=0.5275890119763004115618079766722914408876108660811028;
563:     A[2][1]=0.07241098802369958843819203208518599088698057726988732;
564:     A[2][2]=0; A[2][3]=0;
565:     A[3][0]=0.3990960076760701320627260685975778145384666450351314;
566:     A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
567:     A[3][2]=1.038461646937449311660120300601880176655352737312713;
568:     A[3][3]=0;

570:     b[0]=0.1876410243467238251612921333138006734899663569186926;
571:     b[1]=-0.5952974735769549480478230473706443582188442040780541;
572:     b[2]=0.9717899277217721234705114616271378792182450260943198;
573:     b[3]=0.4358665215084589994160194475295062513822671686978816;

575:     b2[0]=0.2147402862233891404862383521089097657790734483804460;
576:     b2[1]=-0.4851622638849390928209050538171743017757490232519684;
577:     b2[2]=0.8687250025203875511662123688667549217531982787600080;
578:     b2[3]=0.4016969751411624011684543450940068201770721128357014;

580:     binterpt[0][0]=2.2565812720167954547104627844105;
581:     binterpt[1][0]=1.349166413351089573796243820819;
582:     binterpt[2][0]=-2.4695174540533503758652847586647;
583:     binterpt[3][0]=-0.13623023131453465264142184656474;
584:     binterpt[0][1]=-3.0826699111559187902922463354557;
585:     binterpt[1][1]=-2.4689115685996042534544925650515;
586:     binterpt[2][1]=5.7428279814696677152129332773553;
587:     binterpt[3][1]=-0.19124650171414467146619437684812;
588:     binterpt[0][2]=1.0137296634858471607430756831148;
589:     binterpt[1][2]=0.52444768167155973161042570784064;
590:     binterpt[2][2]=-2.3015205996945452158771370439586;
591:     binterpt[3][2]=0.76334325453713832352363565300308;

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



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

609:    Not Collective

611:    Level: advanced

613: .keywords: TSRosW, register, destroy
614: .seealso: TSRosWRegister(), TSRosWRegisterAll()
615: @*/
616: PetscErrorCode TSRosWRegisterDestroy(void)
617: {
618:   PetscErrorCode  ierr;
619:   RosWTableauLink link;

622:   while ((link = RosWTableauList)) {
623:     RosWTableau t = &link->tab;
624:     RosWTableauList = link->next;
625:     PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);
626:     PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);
627:     PetscFree2(t->bembed,t->bembedt);
628:     PetscFree(t->binterpt);
629:     PetscFree(t->name);
630:     PetscFree(link);
631:   }
632:   TSRosWRegisterAllCalled = PETSC_FALSE;
633:   return(0);
634: }

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

643:   Level: developer

645: .keywords: TS, TSRosW, initialize, package
646: .seealso: PetscInitialize()
647: @*/
648: PetscErrorCode TSRosWInitializePackage(void)
649: {

653:   if (TSRosWPackageInitialized) return(0);
654:   TSRosWPackageInitialized = PETSC_TRUE;
655:   TSRosWRegisterAll();
656:   PetscRegisterFinalize(TSRosWFinalizePackage);
657:   return(0);
658: }

662: /*@C
663:   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
664:   called from PetscFinalize().

666:   Level: developer

668: .keywords: Petsc, destroy, package
669: .seealso: PetscFinalize()
670: @*/
671: PetscErrorCode TSRosWFinalizePackage(void)
672: {

676:   TSRosWPackageInitialized = PETSC_FALSE;
677:   TSRosWRegisterDestroy();
678:   return(0);
679: }

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

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

688:    Input Parameters:
689: +  name - identifier for method
690: .  order - approximation order of method
691: .  s - number of stages, this is the dimension of the matrices below
692: .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
693: .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
694: .  b - Step completion table (dimension s)
695: .  bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
696: .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
697: -  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)

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

702:    Level: advanced

704: .keywords: TS, register

706: .seealso: TSRosW
707: @*/
708: PetscErrorCode TSRosWRegister(TSRosWType name,PetscInt order,PetscInt s,const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[],
709:                               PetscInt pinterp,const PetscReal binterpt[])
710: {
711:   PetscErrorCode  ierr;
712:   RosWTableauLink link;
713:   RosWTableau     t;
714:   PetscInt        i,j,k;
715:   PetscScalar     *GammaInv;


724:   PetscMalloc(sizeof(*link),&link);
725:   PetscMemzero(link,sizeof(*link));
726:   t        = &link->tab;
727:   PetscStrallocpy(name,&t->name);
728:   t->order = order;
729:   t->s     = s;
730:   PetscMalloc5(s*s,PetscReal,&t->A,s*s,PetscReal,&t->Gamma,s,PetscReal,&t->b,s,PetscReal,&t->ASum,s,PetscReal,&t->GammaSum);
731:   PetscMalloc5(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag,s*s,PetscReal,&t->GammaExplicitCorr);
732:   PetscMemcpy(t->A,A,s*s*sizeof(A[0]));
733:   PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));
734:   PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));
735:   PetscMemcpy(t->b,b,s*sizeof(b[0]));
736:   if (bembed) {
737:     PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);
738:     PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));
739:   }
740:   for (i=0; i<s; i++) {
741:     t->ASum[i]     = 0;
742:     t->GammaSum[i] = 0;
743:     for (j=0; j<s; j++) {
744:       t->ASum[i]     += A[i*s+j];
745:       t->GammaSum[i] += Gamma[i*s+j];
746:     }
747:   }
748:   PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv); /* Need to use Scalar for inverse, then convert back to Real */
749:   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
750:   for (i=0; i<s; i++) {
751:     if (Gamma[i*s+i] == 0.0) {
752:       GammaInv[i*s+i] = 1.0;
753:       t->GammaZeroDiag[i] = PETSC_TRUE;
754:     } else {
755:       t->GammaZeroDiag[i] = PETSC_FALSE;
756:     }
757:   }

759:   switch (s) {
760:   case 1: GammaInv[0] = 1./GammaInv[0]; break;
761:   case 2: PetscKernel_A_gets_inverse_A_2(GammaInv,0); break;
762:   case 3: PetscKernel_A_gets_inverse_A_3(GammaInv,0); break;
763:   case 4: PetscKernel_A_gets_inverse_A_4(GammaInv,0); break;
764:   case 5: {
765:     PetscInt  ipvt5[5];
766:     MatScalar work5[5*5];
767:     PetscKernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0); break;
768:   }
769:   case 6: PetscKernel_A_gets_inverse_A_6(GammaInv,0); break;
770:   case 7: PetscKernel_A_gets_inverse_A_7(GammaInv,0); break;
771:   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
772:   }
773:   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
774:   PetscFree(GammaInv);

776:   for (i=0; i<s; i++) {
777:     for (k=0; k<i+1; k++) {
778:       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
779:       for (j=k+1; j<i+1; j++) {
780:         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
781:       }
782:     }
783:   }

785:   for (i=0; i<s; i++) {
786:     for (j=0; j<s; j++) {
787:       t->At[i*s+j] = 0;
788:       for (k=0; k<s; k++) {
789:         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
790:       }
791:     }
792:     t->bt[i] = 0;
793:     for (j=0; j<s; j++) {
794:       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
795:     }
796:     if (bembed) {
797:       t->bembedt[i] = 0;
798:       for (j=0; j<s; j++) {
799:         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
800:       }
801:     }
802:   }
803:   t->ccfl = 1.0;                /* Fix this */

805:   t->pinterp = pinterp;
806:   PetscMalloc(s*pinterp*sizeof(binterpt[0]),&t->binterpt);
807:   PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));
808:   link->next = RosWTableauList;
809:   RosWTableauList = link;
810:   return(0);
811: }

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

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

820:    Input Parameters:
821: +  name - identifier for method
822: .  gamma - leading coefficient (diagonal entry)
823: .  a2 - design parameter, see Table 7.2 of Hairer&Wanner
824: .  a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
825: .  b3 - design parameter, see Table 7.2 of Hairer&Wanner
826: .  beta43 - design parameter or PETSC_DEFAULT to use Equation 7.21 of Hairer&Wanner
827: .  e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer

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

834:    Level: developer

836: .keywords: TS, register

838: .seealso: TSRosW, TSRosWRegister()
839: @*/
840: PetscErrorCode TSRosWRegisterRos4(TSRosWType name,PetscReal gamma,PetscReal a2,PetscReal a3,PetscReal b3,PetscReal e4)
841: {
843:   /* Declare numeric constants so they can be quad precision without being truncated at double */
844:   const PetscReal one = 1,two = 2,three = 3,four = 4,five = 5,six = 6,eight = 8,twelve = 12,twenty = 20,twentyfour = 24,
845:     p32 = one/six - gamma + gamma*gamma,
846:     p42 = one/eight - gamma/three,
847:     p43 = one/twelve - gamma/three,
848:     p44 = one/twentyfour - gamma/two + three/two*gamma*gamma - gamma*gamma*gamma,
849:     p56 = one/twenty - gamma/four;
850:   PetscReal   a4,a32,a42,a43,b1,b2,b4,beta2p,beta3p,beta4p,beta32,beta42,beta43,beta32beta2p,beta4jbetajp;
851:   PetscReal   A[4][4],Gamma[4][4],b[4],bm[4];
852:   PetscScalar M[3][3],rhs[3];

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

860:   /* Solve order conditions 7.15a, 7.15c, 7.15e */
861:   M[0][0] = one; M[0][1] = one;      M[0][2] = one;      /* 7.15a */
862:   M[1][0] = 0.0; M[1][1] = a2*a2;    M[1][2] = a4*a4;    /* 7.15c */
863:   M[2][0] = 0.0; M[2][1] = a2*a2*a2; M[2][2] = a4*a4*a4; /* 7.15e */
864:   rhs[0]  = one - b3;
865:   rhs[1]  = one/three - a3*a3*b3;
866:   rhs[2]  = one/four - a3*a3*a3*b3;
867:   PetscKernel_A_gets_inverse_A_3(&M[0][0],0);
868:   b1      = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
869:   b2      = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
870:   b4      = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);

872:   /* Step 3 */
873:   beta43       = (p56 - a2*p43) / (b4*a3*a3*(a3 - a2)); /* 7.21 */
874:   beta32beta2p =  p44 / (b4*beta43);                    /* 7.15h */
875:   beta4jbetajp = (p32 - b3*beta32beta2p) / b4;
876:   M[0][0]      = b2;                                    M[0][1] = b3;                 M[0][2] = b4;
877:   M[1][0]      = a4*a4*beta32beta2p-a3*a3*beta4jbetajp; M[1][1] = a2*a2*beta4jbetajp; M[1][2] = -a2*a2*beta32beta2p;
878:   M[2][0]      = b4*beta43*a3*a3-p43;                   M[2][1] = -b4*beta43*a2*a2;   M[2][2] = 0;
879:   rhs[0]       = one/two - gamma; rhs[1] = 0; rhs[2] = -a2*a2*p32;
880:   PetscKernel_A_gets_inverse_A_3(&M[0][0],0);
881:   beta2p       = PetscRealPart(M[0][0]*rhs[0] + M[0][1]*rhs[1] + M[0][2]*rhs[2]);
882:   beta3p       = PetscRealPart(M[1][0]*rhs[0] + M[1][1]*rhs[1] + M[1][2]*rhs[2]);
883:   beta4p       = PetscRealPart(M[2][0]*rhs[0] + M[2][1]*rhs[1] + M[2][2]*rhs[2]);

885:   /* Step 4: back-substitute */
886:   beta32 = beta32beta2p / beta2p;
887:   beta42 = (beta4jbetajp - beta43*beta3p) / beta2p;

889:   /* Step 5: 7.15f and 7.20, then 7.16 */
890:   a43 = 0;
891:   a32 = p42 / (b3*a3*beta2p + b4*a4*beta2p);
892:   a42 = a32;

894:   A[0][0]     = 0;          A[0][1] = 0;   A[0][2] = 0;   A[0][3] = 0;
895:   A[1][0]     = a2;         A[1][1] = 0;   A[1][2] = 0;   A[1][3] = 0;
896:   A[2][0]     = a3-a32;     A[2][1] = a32; A[2][2] = 0;   A[2][3] = 0;
897:   A[3][0]     = a4-a43-a42; A[3][1] = a42; A[3][2] = a43; A[3][3] = 0;
898:   Gamma[0][0] = gamma;                        Gamma[0][1] = 0;              Gamma[0][2] = 0;              Gamma[0][3] = 0;
899:   Gamma[1][0] = beta2p-A[1][0];               Gamma[1][1] = gamma;          Gamma[1][2] = 0;              Gamma[1][3] = 0;
900:   Gamma[2][0] = beta3p-beta32-A[2][0];        Gamma[2][1] = beta32-A[2][1]; Gamma[2][2] = gamma;          Gamma[2][3] = 0;
901:   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;
902:   b[0] = b1; b[1] = b2; b[2] = b3; b[3] = b4;

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

910:   {
911:     const PetscReal misfit = a2*a2*bm[1] + a3*a3*bm[2] + a4*a4*bm[3] - one/three;
912:     if (PetscAbs(misfit) > PETSC_SMALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Assumptions violated, could not construct a third order embedded method");
913:   }
914:   TSRosWRegister(name,4,4,&A[0][0],&Gamma[0][0],b,bm,0,NULL);
915:   return(0);
916: }

920: /*
921:  The step completion formula is

923:  x1 = x0 + b^T Y

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

928:  x1e = x0 + be^T Y
929:      = x1 - b^T Y + be^T Y
930:      = x1 + (be - b)^T Y

932:  so we can evaluate the method of different order even after the step has been optimistically completed.
933: */
934: static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec U,PetscBool *done)
935: {
936:   TS_RosW        *ros = (TS_RosW*)ts->data;
937:   RosWTableau    tab  = ros->tableau;
938:   PetscScalar    *w   = ros->work;
939:   PetscInt       i;

943:   if (order == tab->order) {
944:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
945:       VecCopy(ts->vec_sol,U);
946:       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
947:       VecMAXPY(U,tab->s,w,ros->Y);
948:     } else {VecCopy(ts->vec_sol,U);}
949:     if (done) *done = PETSC_TRUE;
950:     return(0);
951:   } else if (order == tab->order-1) {
952:     if (!tab->bembedt) goto unavailable;
953:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
954:       VecCopy(ts->vec_sol,U);
955:       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
956:       VecMAXPY(U,tab->s,w,ros->Y);
957:     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
958:       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
959:       VecCopy(ts->vec_sol,U);
960:       VecMAXPY(U,tab->s,w,ros->Y);
961:     }
962:     if (done) *done = PETSC_TRUE;
963:     return(0);
964:   }
965:   unavailable:
966:   if (done) *done = PETSC_FALSE;
967:   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Rosenbrock-W '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order);
968:   return(0);
969: }

973: static PetscErrorCode TSStep_RosW(TS ts)
974: {
975:   TS_RosW         *ros = (TS_RosW*)ts->data;
976:   RosWTableau     tab  = ros->tableau;
977:   const PetscInt  s    = tab->s;
978:   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
979:   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
980:   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
981:   PetscScalar     *w   = ros->work;
982:   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
983:   SNES            snes;
984:   TSAdapt         adapt;
985:   PetscInt        i,j,its,lits,reject,next_scheme;
986:   PetscReal       next_time_step;
987:   PetscBool       accept;
988:   PetscErrorCode  ierr;
989:   MatStructure    str;

992:   TSGetSNES(ts,&snes);
993:   next_time_step = ts->time_step;
994:   accept         = PETSC_TRUE;
995:   ros->status    = TS_STEP_INCOMPLETE;

997:   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
998:     const PetscReal h = ts->time_step;
999:     TSPreStep(ts);
1000:     VecCopy(ts->vec_sol,ros->VecSolPrev); /*move this at the end*/
1001:     for (i=0; i<s; i++) {
1002:       ros->stage_time = ts->ptime + h*ASum[i];
1003:       TSPreStage(ts,ros->stage_time);
1004:       if (GammaZeroDiag[i]) {
1005:         ros->stage_explicit = PETSC_TRUE;
1006:         ros->scoeff         = 1.;
1007:       } else {
1008:         ros->stage_explicit = PETSC_FALSE;
1009:         ros->scoeff         = 1./Gamma[i*s+i];
1010:       }

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

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

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

1023:       if (!ros->stage_explicit) {
1024:         if (!ros->recompute_jacobian && !i) {
1025:           SNESSetLagJacobian(snes,-2); /* Recompute the Jacobian on this solve, but not again */
1026:         }
1027:         SNESSolve(snes,NULL,Y[i]);
1028:         SNESGetIterationNumber(snes,&its);
1029:         SNESGetLinearSolveIterations(snes,&lits);
1030:         ts->snes_its += its; ts->ksp_its += lits;
1031:         TSGetAdapt(ts,&adapt);
1032:         TSAdaptCheckStage(adapt,ts,&accept);
1033:         if (!accept) goto reject_step;
1034:       } else {
1035:         Mat J,Jp;
1036:         VecZeroEntries(Ydot); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
1037:         TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);
1038:         VecScale(Y[i],-1.0);
1039:         VecAXPY(Y[i],-1.0,Zdot); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/

1041:         VecZeroEntries(Zstage); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
1042:         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
1043:         VecMAXPY(Zstage,i,w,Y);
1044:         /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
1045:         str  = SAME_NONZERO_PATTERN;
1046:         TSGetIJacobian(ts,&J,&Jp,NULL,NULL);
1047:         TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);
1048:         MatMult(J,Zstage,Zdot);

1050:         VecAXPY(Y[i],-1.0,Zdot);
1051:         VecScale(Y[i],h);
1052:         ts->ksp_its += 1;
1053:       }
1054:     }
1055:     TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);
1056:     ros->status = TS_STEP_PENDING;

1058:     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
1059:     TSGetAdapt(ts,&adapt);
1060:     TSAdaptCandidatesClear(adapt);
1061:     TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);
1062:     TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);
1063:     if (accept) {
1064:       /* ignore next_scheme for now */
1065:       ts->ptime    += ts->time_step;
1066:       ts->time_step = next_time_step;
1067:       ts->steps++;
1068:       ros->status = TS_STEP_COMPLETE;
1069:       break;
1070:     } else {                    /* Roll back the current step */
1071:       for (i=0; i<s; i++) w[i] = -tab->bt[i];
1072:       VecMAXPY(ts->vec_sol,s,w,Y);
1073:       ts->time_step = next_time_step;
1074:       ros->status   = TS_STEP_INCOMPLETE;
1075:     }
1076: reject_step: continue;
1077:   }
1078:   if (ros->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
1079:   return(0);
1080: }

1084: static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec U)
1085: {
1086:   TS_RosW         *ros = (TS_RosW*)ts->data;
1087:   PetscInt        s    = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j;
1088:   PetscReal       h;
1089:   PetscReal       tt,t;
1090:   PetscScalar     *bt;
1091:   const PetscReal *Bt = ros->tableau->binterpt;
1092:   PetscErrorCode  ierr;
1093:   const PetscReal *GammaInv = ros->tableau->GammaInv;
1094:   PetscScalar     *w        = ros->work;
1095:   Vec             *Y        = ros->Y;

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

1100:   switch (ros->status) {
1101:   case TS_STEP_INCOMPLETE:
1102:   case TS_STEP_PENDING:
1103:     h = ts->time_step;
1104:     t = (itime - ts->ptime)/h;
1105:     break;
1106:   case TS_STEP_COMPLETE:
1107:     h = ts->time_step_prev;
1108:     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
1109:     break;
1110:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
1111:   }
1112:   PetscMalloc(s*sizeof(bt[0]),&bt);
1113:   for (i=0; i<s; i++) bt[i] = 0;
1114:   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
1115:     for (i=0; i<s; i++) {
1116:       bt[i] += Bt[i*pinterp+j] * tt;
1117:     }
1118:   }

1120:   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1121:   /*U<-0*/
1122:   VecZeroEntries(U);

1124:   /*U<- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
1125:   for (j=0; j<s; j++) w[j]=0;
1126:   for (j=0; j<s; j++) {
1127:     for (i=j; i<s; i++) {
1128:       w[j] +=  bt[i]*GammaInv[i*s+j];
1129:     }
1130:   }
1131:   VecMAXPY(U,i,w,Y);

1133:   /*X<-y(t) + X*/
1134:   VecAXPY(U,1.0,ros->VecSolPrev);

1136:   PetscFree(bt);
1137:   return(0);
1138: }

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

1150:   if (!ros->tableau) return(0);
1151:   s    = ros->tableau->s;
1152:   VecDestroyVecs(s,&ros->Y);
1153:   VecDestroy(&ros->Ydot);
1154:   VecDestroy(&ros->Ystage);
1155:   VecDestroy(&ros->Zdot);
1156:   VecDestroy(&ros->Zstage);
1157:   VecDestroy(&ros->VecSolPrev);
1158:   PetscFree(ros->work);
1159:   return(0);
1160: }

1164: static PetscErrorCode TSDestroy_RosW(TS ts)
1165: {

1169:   TSReset_RosW(ts);
1170:   PetscFree(ts->data);
1171:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",NULL);
1172:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",NULL);
1173:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",NULL);
1174:   return(0);
1175: }


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

1186:   if (Ydot) {
1187:     if (dm && dm != ts->dm) {
1188:       DMGetNamedGlobalVector(dm,"TSRosW_Ydot",Ydot);
1189:     } else *Ydot = rw->Ydot;
1190:   }
1191:   if (Zdot) {
1192:     if (dm && dm != ts->dm) {
1193:       DMGetNamedGlobalVector(dm,"TSRosW_Zdot",Zdot);
1194:     } else *Zdot = rw->Zdot;
1195:   }
1196:   if (Ystage) {
1197:     if (dm && dm != ts->dm) {
1198:       DMGetNamedGlobalVector(dm,"TSRosW_Ystage",Ystage);
1199:     } else *Ystage = rw->Ystage;
1200:   }
1201:   if (Zstage) {
1202:     if (dm && dm != ts->dm) {
1203:       DMGetNamedGlobalVector(dm,"TSRosW_Zstage",Zstage);
1204:     } else *Zstage = rw->Zstage;
1205:   }
1206:   return(0);
1207: }


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

1217:   if (Ydot) {
1218:     if (dm && dm != ts->dm) {
1219:       DMRestoreNamedGlobalVector(dm,"TSRosW_Ydot",Ydot);
1220:     }
1221:   }
1222:   if (Zdot) {
1223:     if (dm && dm != ts->dm) {
1224:       DMRestoreNamedGlobalVector(dm,"TSRosW_Zdot",Zdot);
1225:     }
1226:   }
1227:   if (Ystage) {
1228:     if (dm && dm != ts->dm) {
1229:       DMRestoreNamedGlobalVector(dm,"TSRosW_Ystage",Ystage);
1230:     }
1231:   }
1232:   if (Zstage) {
1233:     if (dm && dm != ts->dm) {
1234:       DMRestoreNamedGlobalVector(dm,"TSRosW_Zstage",Zstage);
1235:     }
1236:   }
1237:   return(0);
1238: }

1242: static PetscErrorCode DMCoarsenHook_TSRosW(DM fine,DM coarse,void *ctx)
1243: {
1245:   return(0);
1246: }

1250: static PetscErrorCode DMRestrictHook_TSRosW(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1251: {
1252:   TS             ts = (TS)ctx;
1254:   Vec            Ydot,Zdot,Ystage,Zstage;
1255:   Vec            Ydotc,Zdotc,Ystagec,Zstagec;

1258:   TSRosWGetVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage);
1259:   TSRosWGetVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec);
1260:   MatRestrict(restrct,Ydot,Ydotc);
1261:   VecPointwiseMult(Ydotc,rscale,Ydotc);
1262:   MatRestrict(restrct,Ystage,Ystagec);
1263:   VecPointwiseMult(Ystagec,rscale,Ystagec);
1264:   MatRestrict(restrct,Zdot,Zdotc);
1265:   VecPointwiseMult(Zdotc,rscale,Zdotc);
1266:   MatRestrict(restrct,Zstage,Zstagec);
1267:   VecPointwiseMult(Zstagec,rscale,Zstagec);
1268:   TSRosWRestoreVecs(ts,fine,&Ydot,&Ystage,&Zdot,&Zstage);
1269:   TSRosWRestoreVecs(ts,coarse,&Ydotc,&Ystagec,&Zdotc,&Zstagec);
1270:   return(0);
1271: }


1276: static PetscErrorCode DMSubDomainHook_TSRosW(DM fine,DM coarse,void *ctx)
1277: {
1279:   return(0);
1280: }

1284: static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1285: {
1286:   TS             ts = (TS)ctx;
1288:   Vec            Ydot,Zdot,Ystage,Zstage;
1289:   Vec            Ydots,Zdots,Ystages,Zstages;

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

1295:   VecScatterBegin(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD);
1296:   VecScatterEnd(gscat,Ydot,Ydots,INSERT_VALUES,SCATTER_FORWARD);

1298:   VecScatterBegin(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD);
1299:   VecScatterEnd(gscat,Ystage,Ystages,INSERT_VALUES,SCATTER_FORWARD);

1301:   VecScatterBegin(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD);
1302:   VecScatterEnd(gscat,Zdot,Zdots,INSERT_VALUES,SCATTER_FORWARD);

1304:   VecScatterBegin(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD);
1305:   VecScatterEnd(gscat,Zstage,Zstages,INSERT_VALUES,SCATTER_FORWARD);

1307:   TSRosWRestoreVecs(ts,dm,&Ydot,&Ystage,&Zdot,&Zstage);
1308:   TSRosWRestoreVecs(ts,subdm,&Ydots,&Ystages,&Zdots,&Zstages);
1309:   return(0);
1310: }

1312: /*
1313:   This defines the nonlinear equation that is to be solved with SNES
1314:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1315: */
1318: static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec U,Vec F,TS ts)
1319: {
1320:   TS_RosW        *ros = (TS_RosW*)ts->data;
1322:   Vec            Ydot,Zdot,Ystage,Zstage;
1323:   PetscReal      shift = ros->scoeff / ts->time_step;
1324:   DM             dm,dmsave;

1327:   SNESGetDM(snes,&dm);
1328:   TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1329:   VecWAXPY(Ydot,shift,U,Zdot);    /* Ydot = shift*U + Zdot */
1330:   VecWAXPY(Ystage,1.0,U,Zstage);  /* Ystage = U + Zstage */
1331:   dmsave = ts->dm;
1332:   ts->dm = dm;
1333:   TSComputeIFunction(ts,ros->stage_time,Ystage,Ydot,F,PETSC_FALSE);
1334:   ts->dm = dmsave;
1335:   TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1336:   return(0);
1337: }

1341: static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec U,Mat *A,Mat *B,MatStructure *str,TS ts)
1342: {
1343:   TS_RosW        *ros = (TS_RosW*)ts->data;
1344:   Vec            Ydot,Zdot,Ystage,Zstage;
1345:   PetscReal      shift = ros->scoeff / ts->time_step;
1347:   DM             dm,dmsave;

1350:   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
1351:   SNESGetDM(snes,&dm);
1352:   TSRosWGetVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1353:   dmsave = ts->dm;
1354:   ts->dm = dm;
1355:   TSComputeIJacobian(ts,ros->stage_time,Ystage,Ydot,shift,A,B,str,PETSC_TRUE);
1356:   ts->dm = dmsave;
1357:   TSRosWRestoreVecs(ts,dm,&Ydot,&Zdot,&Ystage,&Zstage);
1358:   return(0);
1359: }

1363: static PetscErrorCode TSSetUp_RosW(TS ts)
1364: {
1365:   TS_RosW        *ros = (TS_RosW*)ts->data;
1366:   RosWTableau    tab  = ros->tableau;
1367:   PetscInt       s    = tab->s;
1369:   DM             dm;

1372:   if (!ros->tableau) {
1373:     TSRosWSetType(ts,TSRosWDefault);
1374:   }
1375:   VecDuplicateVecs(ts->vec_sol,s,&ros->Y);
1376:   VecDuplicate(ts->vec_sol,&ros->Ydot);
1377:   VecDuplicate(ts->vec_sol,&ros->Ystage);
1378:   VecDuplicate(ts->vec_sol,&ros->Zdot);
1379:   VecDuplicate(ts->vec_sol,&ros->Zstage);
1380:   VecDuplicate(ts->vec_sol,&ros->VecSolPrev);
1381:   PetscMalloc(s*sizeof(ros->work[0]),&ros->work);
1382:   TSGetDM(ts,&dm);
1383:   if (dm) {
1384:     DMCoarsenHookAdd(dm,DMCoarsenHook_TSRosW,DMRestrictHook_TSRosW,ts);
1385:     DMSubDomainHookAdd(dm,DMSubDomainHook_TSRosW,DMSubDomainRestrictHook_TSRosW,ts);
1386:   }
1387:   return(0);
1388: }
1389: /*------------------------------------------------------------*/

1393: static PetscErrorCode TSSetFromOptions_RosW(TS ts)
1394: {
1395:   TS_RosW        *ros = (TS_RosW*)ts->data;
1397:   char           rostype[256];

1400:   PetscOptionsHead("RosW ODE solver options");
1401:   {
1402:     RosWTableauLink link;
1403:     PetscInt        count,choice;
1404:     PetscBool       flg;
1405:     const char      **namelist;
1406:     SNES            snes;

1408:     PetscStrncpy(rostype,TSRosWDefault,sizeof(rostype));
1409:     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
1410:     PetscMalloc(count*sizeof(char*),&namelist);
1411:     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1412:     PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);
1413:     TSRosWSetType(ts,flg ? namelist[choice] : rostype);
1414:     PetscFree(namelist);

1416:     PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,NULL);

1418:     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1419:     TSGetSNES(ts,&snes);
1420:     if (!((PetscObject)snes)->type_name) {
1421:       SNESSetType(snes,SNESKSPONLY);
1422:     }
1423:     SNESSetFromOptions(snes);
1424:   }
1425:   PetscOptionsTail();
1426:   return(0);
1427: }

1431: static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
1432: {
1434:   PetscInt       i;
1435:   size_t         left,count;
1436:   char           *p;

1439:   for (i=0,p=buf,left=len; i<n; i++) {
1440:     PetscSNPrintfCount(p,left,fmt,&count,x[i]);
1441:     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
1442:     left -= count;
1443:     p    += count;
1444:     *p++  = ' ';
1445:   }
1446:   p[i ? 0 : -1] = 0;
1447:   return(0);
1448: }

1452: static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
1453: {
1454:   TS_RosW        *ros = (TS_RosW*)ts->data;
1455:   RosWTableau    tab  = ros->tableau;
1456:   PetscBool      iascii;
1458:   TSAdapt        adapt;

1461:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1462:   if (iascii) {
1463:     TSRosWType rostype;
1464:     PetscInt   i;
1465:     PetscReal  abscissa[512];
1466:     char       buf[512];
1467:     TSRosWGetType(ts,&rostype);
1468:     PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);
1469:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ASum);
1470:     PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);
1471:     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1472:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,abscissa);
1473:     PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);
1474:   }
1475:   TSGetAdapt(ts,&adapt);
1476:   TSAdaptView(adapt,viewer);
1477:   SNESView(ts->snes,viewer);
1478:   return(0);
1479: }

1483: /*@C
1484:   TSRosWSetType - Set the type of Rosenbrock-W scheme

1486:   Logically collective

1488:   Input Parameter:
1489: +  ts - timestepping context
1490: -  rostype - type of Rosenbrock-W scheme

1492:   Level: beginner

1494: .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
1495: @*/
1496: PetscErrorCode TSRosWSetType(TS ts,TSRosWType rostype)
1497: {

1502:   PetscTryMethod(ts,"TSRosWSetType_C",(TS,TSRosWType),(ts,rostype));
1503:   return(0);
1504: }

1508: /*@C
1509:   TSRosWGetType - Get the type of Rosenbrock-W scheme

1511:   Logically collective

1513:   Input Parameter:
1514: .  ts - timestepping context

1516:   Output Parameter:
1517: .  rostype - type of Rosenbrock-W scheme

1519:   Level: intermediate

1521: .seealso: TSRosWGetType()
1522: @*/
1523: PetscErrorCode TSRosWGetType(TS ts,TSRosWType *rostype)
1524: {

1529:   PetscUseMethod(ts,"TSRosWGetType_C",(TS,TSRosWType*),(ts,rostype));
1530:   return(0);
1531: }

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

1538:   Logically collective

1540:   Input Parameter:
1541: +  ts - timestepping context
1542: -  flg - PETSC_TRUE to recompute the Jacobian at each stage

1544:   Level: intermediate

1546: .seealso: TSRosWGetType()
1547: @*/
1548: PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1549: {

1554:   PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));
1555:   return(0);
1556: }

1560: PetscErrorCode  TSRosWGetType_RosW(TS ts,TSRosWType *rostype)
1561: {
1562:   TS_RosW        *ros = (TS_RosW*)ts->data;

1566:   if (!ros->tableau) {TSRosWSetType(ts,TSRosWDefault);}
1567:   *rostype = ros->tableau->name;
1568:   return(0);
1569: }

1573: PetscErrorCode  TSRosWSetType_RosW(TS ts,TSRosWType rostype)
1574: {
1575:   TS_RosW         *ros = (TS_RosW*)ts->data;
1576:   PetscErrorCode  ierr;
1577:   PetscBool       match;
1578:   RosWTableauLink link;

1581:   if (ros->tableau) {
1582:     PetscStrcmp(ros->tableau->name,rostype,&match);
1583:     if (match) return(0);
1584:   }
1585:   for (link = RosWTableauList; link; link=link->next) {
1586:     PetscStrcmp(link->tab.name,rostype,&match);
1587:     if (match) {
1588:       TSReset_RosW(ts);
1589:       ros->tableau = &link->tab;
1590:       return(0);
1591:     }
1592:   }
1593:   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1594:   return(0);
1595: }

1599: PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1600: {
1601:   TS_RosW *ros = (TS_RosW*)ts->data;

1604:   ros->recompute_jacobian = flg;
1605:   return(0);
1606: }


1609: /* ------------------------------------------------------------ */
1610: /*MC
1611:       TSROSW - ODE solver using Rosenbrock-W schemes

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

1617:   Notes:
1618:   This method currently only works with autonomous ODE and DAE.

1620:   Developer notes:
1621:   Rosenbrock-W methods are typically specified for autonomous ODE

1623: $  udot = f(u)

1625:   by the stage equations

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

1629:   and step completion formula

1631: $  u_1 = u_0 + sum_j b_j k_j

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

1637: $  y_i = gamma_ij k_j

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

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

1643:   to rewrite the method as

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

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

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

1652:    or, more compactly in tensor notation

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

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

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

1662:    with initial guess y_i = 0.

1664:   Level: beginner

1666: .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister(), TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3,
1667:            TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWGRK4T, TSROSWSHAMP4, TSROSWVELDD4, TSROSW4L
1668: M*/
1671: PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1672: {
1673:   TS_RosW        *ros;

1677: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1678:   TSRosWInitializePackage();
1679: #endif

1681:   ts->ops->reset          = TSReset_RosW;
1682:   ts->ops->destroy        = TSDestroy_RosW;
1683:   ts->ops->view           = TSView_RosW;
1684:   ts->ops->setup          = TSSetUp_RosW;
1685:   ts->ops->step           = TSStep_RosW;
1686:   ts->ops->interpolate    = TSInterpolate_RosW;
1687:   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1688:   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1689:   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1690:   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;

1692:   PetscNewLog(ts,TS_RosW,&ros);
1693:   ts->data = (void*)ros;

1695:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWGetType_C",TSRosWGetType_RosW);
1696:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetType_C",TSRosWSetType_RosW);
1697:   PetscObjectComposeFunction((PetscObject)ts,"TSRosWSetRecomputeJacobian_C",TSRosWSetRecomputeJacobian_RosW);
1698:   return(0);
1699: }