1e27a552bSJed Brown /* 261692a83SJed Brown Code for timestepping with Rosenbrock W methods 3e27a552bSJed Brown 4e27a552bSJed Brown Notes: 5e27a552bSJed Brown The general system is written as 6e27a552bSJed Brown 761692a83SJed Brown G(t,X,Xdot) = F(t,X) 8e27a552bSJed Brown 961692a83SJed Brown where G represents the stiff part of the physics and F represents the non-stiff part. 1061692a83SJed Brown This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian. 11e27a552bSJed Brown 12e27a552bSJed Brown */ 13e27a552bSJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 14e27a552bSJed Brown 1561692a83SJed Brown #include <../src/mat/blockinvert.h> 1661692a83SJed Brown 17ae6f9490SJed Brown static const TSRosWType TSRosWDefault = TSROSWRA34PW2; 18e27a552bSJed Brown static PetscBool TSRosWRegisterAllCalled; 19e27a552bSJed Brown static PetscBool TSRosWPackageInitialized; 20e27a552bSJed Brown 2161692a83SJed Brown typedef struct _RosWTableau *RosWTableau; 2261692a83SJed Brown struct _RosWTableau { 23e27a552bSJed Brown char *name; 24e27a552bSJed Brown PetscInt order; /* Classical approximation order of the method */ 25e27a552bSJed Brown PetscInt s; /* Number of stages */ 26f4aed992SEmil Constantinescu PetscInt pinterp; /* Interpolation order */ 2761692a83SJed Brown PetscReal *A; /* Propagation table, strictly lower triangular */ 2861692a83SJed Brown PetscReal *Gamma; /* Stage table, lower triangular with nonzero diagonal */ 29c17803e7SJed Brown PetscBool *GammaZeroDiag; /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */ 3043b21953SEmil Constantinescu PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/ 3161692a83SJed Brown PetscReal *b; /* Step completion table */ 32fe7e6d57SJed Brown PetscReal *bembed; /* Step completion table for embedded method of order one less */ 3361692a83SJed Brown PetscReal *ASum; /* Row sum of A */ 3461692a83SJed Brown PetscReal *GammaSum; /* Row sum of Gamma, only needed for non-autonomous systems */ 3561692a83SJed Brown PetscReal *At; /* Propagation table in transformed variables */ 3661692a83SJed Brown PetscReal *bt; /* Step completion table in transformed variables */ 37fe7e6d57SJed Brown PetscReal *bembedt; /* Step completion table of order one less in transformed variables */ 3861692a83SJed Brown PetscReal *GammaInv; /* Inverse of Gamma, used for transformed variables */ 398d59e960SJed Brown PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 40*3ca35412SEmil Constantinescu PetscReal *binterpt; /* Dense output formula */ 41e27a552bSJed Brown }; 4261692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink; 4361692a83SJed Brown struct _RosWTableauLink { 4461692a83SJed Brown struct _RosWTableau tab; 4561692a83SJed Brown RosWTableauLink next; 46e27a552bSJed Brown }; 4761692a83SJed Brown static RosWTableauLink RosWTableauList; 48e27a552bSJed Brown 49e27a552bSJed Brown typedef struct { 5061692a83SJed Brown RosWTableau tableau; 5161692a83SJed Brown Vec *Y; /* States computed during the step, used to complete the step */ 52e27a552bSJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 5361692a83SJed Brown Vec Ystage; /* Work vector for the state value at each stage */ 5461692a83SJed Brown Vec Zdot; /* Ydot = Zdot + shift*Y */ 5561692a83SJed Brown Vec Zstage; /* Y = Zstage + Y */ 56*3ca35412SEmil Constantinescu Vec VecSolPrev; /* Work vector holding the solution from the previous step (used for interpolation)*/ 571c3436cfSJed Brown PetscScalar *work; /* Scalar work space of length number of stages, used to prepare VecMAXPY() */ 58e27a552bSJed Brown PetscReal shift; 59e27a552bSJed Brown PetscReal stage_time; 60c17803e7SJed Brown PetscReal stage_explicit; /* Flag indicates that the current stage is explicit */ 6161692a83SJed Brown PetscBool recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */ 62108c343cSJed Brown TSStepStatus status; 63e27a552bSJed Brown } TS_RosW; 64e27a552bSJed Brown 65fe7e6d57SJed Brown /*MC 663606a31eSEmil Constantinescu TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method). 673606a31eSEmil Constantinescu 683606a31eSEmil Constantinescu Only an approximate Jacobian is needed. 693606a31eSEmil Constantinescu 703606a31eSEmil Constantinescu Level: intermediate 713606a31eSEmil Constantinescu 723606a31eSEmil Constantinescu .seealso: TSROSW 733606a31eSEmil Constantinescu M*/ 743606a31eSEmil Constantinescu 753606a31eSEmil Constantinescu /*MC 763606a31eSEmil Constantinescu TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method). 773606a31eSEmil Constantinescu 783606a31eSEmil Constantinescu Only an approximate Jacobian is needed. 793606a31eSEmil Constantinescu 803606a31eSEmil Constantinescu Level: intermediate 813606a31eSEmil Constantinescu 823606a31eSEmil Constantinescu .seealso: TSROSW 833606a31eSEmil Constantinescu M*/ 843606a31eSEmil Constantinescu 853606a31eSEmil Constantinescu /*MC 86fe7e6d57SJed Brown TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme. 87fe7e6d57SJed Brown 88fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P. 89fe7e6d57SJed Brown 90fe7e6d57SJed Brown Level: intermediate 91fe7e6d57SJed Brown 92fe7e6d57SJed Brown .seealso: TSROSW 93fe7e6d57SJed Brown M*/ 94fe7e6d57SJed Brown 95fe7e6d57SJed Brown /*MC 96fe7e6d57SJed Brown TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme. 97fe7e6d57SJed Brown 98fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M. 99fe7e6d57SJed Brown 100fe7e6d57SJed Brown Level: intermediate 101fe7e6d57SJed Brown 102fe7e6d57SJed Brown .seealso: TSROSW 103fe7e6d57SJed Brown M*/ 104fe7e6d57SJed Brown 105fe7e6d57SJed Brown /*MC 106fe7e6d57SJed Brown TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1. 107fe7e6d57SJed Brown 108fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 109fe7e6d57SJed Brown 110fe7e6d57SJed Brown 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. 111fe7e6d57SJed Brown 112fe7e6d57SJed Brown References: 113fe7e6d57SJed Brown Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 114fe7e6d57SJed Brown 115fe7e6d57SJed Brown Level: intermediate 116fe7e6d57SJed Brown 117fe7e6d57SJed Brown .seealso: TSROSW 118fe7e6d57SJed Brown M*/ 119fe7e6d57SJed Brown 120fe7e6d57SJed Brown /*MC 121fe7e6d57SJed Brown TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1. 122fe7e6d57SJed Brown 123fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 124fe7e6d57SJed Brown 125fe7e6d57SJed Brown This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48. 126fe7e6d57SJed Brown 127fe7e6d57SJed Brown References: 128fe7e6d57SJed Brown Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 129fe7e6d57SJed Brown 130fe7e6d57SJed Brown Level: intermediate 131fe7e6d57SJed Brown 132fe7e6d57SJed Brown .seealso: TSROSW 133fe7e6d57SJed Brown M*/ 134fe7e6d57SJed Brown 135ef3c5b88SJed Brown /*MC 136ef3c5b88SJed Brown TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme 137ef3c5b88SJed Brown 138ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 139ef3c5b88SJed Brown 140ef3c5b88SJed Brown Both the third order and embedded second order methods are stiffly accurate and L-stable. 141ef3c5b88SJed Brown 142ef3c5b88SJed Brown References: 143ef3c5b88SJed Brown Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 144ef3c5b88SJed Brown 145ef3c5b88SJed Brown Level: intermediate 146ef3c5b88SJed Brown 147ef3c5b88SJed Brown .seealso: TSROSW, TSROSWSANDU3 148ef3c5b88SJed Brown M*/ 149ef3c5b88SJed Brown 150ef3c5b88SJed Brown /*MC 151ef3c5b88SJed Brown TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme 152ef3c5b88SJed Brown 153ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 154ef3c5b88SJed Brown 155ef3c5b88SJed Brown The third order method is L-stable, but not stiffly accurate. 156ef3c5b88SJed Brown The second order embedded method is strongly A-stable with R(infty) = 0.5. 157ef3c5b88SJed Brown The internal stages are L-stable. 158ef3c5b88SJed Brown This method is called ROS3 in the paper. 159ef3c5b88SJed Brown 160ef3c5b88SJed Brown References: 161ef3c5b88SJed Brown Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 162ef3c5b88SJed Brown 163ef3c5b88SJed Brown Level: intermediate 164ef3c5b88SJed Brown 165ef3c5b88SJed Brown .seealso: TSROSW, TSROSWRODAS3 166ef3c5b88SJed Brown M*/ 167ef3c5b88SJed Brown 168961f28d0SJed Brown /*MC 169961f28d0SJed Brown TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages 170961f28d0SJed Brown 171961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 172961f28d0SJed Brown 173961f28d0SJed Brown A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3) 174961f28d0SJed Brown 175961f28d0SJed Brown References: 176961f28d0SJed Brown Emil Constantinescu 177961f28d0SJed Brown 178961f28d0SJed Brown Level: intermediate 179961f28d0SJed Brown 18043b21953SEmil Constantinescu .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP 181961f28d0SJed Brown M*/ 182961f28d0SJed Brown 183961f28d0SJed Brown /*MC 184961f28d0SJed Brown TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 185961f28d0SJed Brown 186961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 187961f28d0SJed Brown 188961f28d0SJed Brown L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 189961f28d0SJed Brown 190961f28d0SJed Brown References: 191961f28d0SJed Brown Emil Constantinescu 192961f28d0SJed Brown 193961f28d0SJed Brown Level: intermediate 194961f28d0SJed Brown 19543b21953SEmil Constantinescu .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP 196961f28d0SJed Brown M*/ 197961f28d0SJed Brown 198961f28d0SJed Brown /*MC 19943b21953SEmil Constantinescu TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 200961f28d0SJed Brown 201961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 202961f28d0SJed Brown 203961f28d0SJed Brown L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 204961f28d0SJed Brown 205961f28d0SJed Brown References: 206961f28d0SJed Brown Emil Constantinescu 207961f28d0SJed Brown 208961f28d0SJed Brown Level: intermediate 209961f28d0SJed Brown 210961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP 211961f28d0SJed Brown M*/ 212961f28d0SJed Brown 213e27a552bSJed Brown #undef __FUNCT__ 214e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll" 215e27a552bSJed Brown /*@C 216e27a552bSJed Brown TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW 217e27a552bSJed Brown 218e27a552bSJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 219e27a552bSJed Brown 220e27a552bSJed Brown Level: advanced 221e27a552bSJed Brown 222e27a552bSJed Brown .keywords: TS, TSRosW, register, all 223e27a552bSJed Brown 224e27a552bSJed Brown .seealso: TSRosWRegisterDestroy() 225e27a552bSJed Brown @*/ 226e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void) 227e27a552bSJed Brown { 228e27a552bSJed Brown PetscErrorCode ierr; 229e27a552bSJed Brown 230e27a552bSJed Brown PetscFunctionBegin; 231e27a552bSJed Brown if (TSRosWRegisterAllCalled) PetscFunctionReturn(0); 232e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_TRUE; 233e27a552bSJed Brown 234e27a552bSJed Brown { 2353606a31eSEmil Constantinescu const PetscReal 2363606a31eSEmil Constantinescu A = 0, 2373606a31eSEmil Constantinescu Gamma = 1, 2383606a31eSEmil Constantinescu b = 1; 239f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 2403606a31eSEmil Constantinescu } 2413606a31eSEmil Constantinescu 2423606a31eSEmil Constantinescu { 2433606a31eSEmil Constantinescu const PetscReal 2443606a31eSEmil Constantinescu A= 0, 2453606a31eSEmil Constantinescu Gamma = 0.5, 2463606a31eSEmil Constantinescu b = 1; 247f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 2483606a31eSEmil Constantinescu } 2493606a31eSEmil Constantinescu 2503606a31eSEmil Constantinescu { 25161692a83SJed Brown const PetscReal g = 1. + 1./PetscSqrtReal(2.0); 252e27a552bSJed Brown const PetscReal 25361692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 25461692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2551c3436cfSJed Brown b[2] = {0.5,0.5}, 2561c3436cfSJed Brown b1[2] = {1.0,0.0}; 257f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1,0,PETSC_NULL);CHKERRQ(ierr); 258e27a552bSJed Brown } 259e27a552bSJed Brown { 26061692a83SJed Brown const PetscReal g = 1. - 1./PetscSqrtReal(2.0); 261e27a552bSJed Brown const PetscReal 26261692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 26361692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2641c3436cfSJed Brown b[2] = {0.5,0.5}, 2651c3436cfSJed Brown b1[2] = {1.0,0.0}; 266f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1,0,PETSC_NULL);CHKERRQ(ierr); 267fe7e6d57SJed Brown } 268fe7e6d57SJed Brown { 269fe7e6d57SJed Brown const PetscReal g = 7.8867513459481287e-01; 270fe7e6d57SJed Brown const PetscReal 271fe7e6d57SJed Brown A[3][3] = {{0,0,0}, 272fe7e6d57SJed Brown {1.5773502691896257e+00,0,0}, 273fe7e6d57SJed Brown {0.5,0,0}}, 274fe7e6d57SJed Brown Gamma[3][3] = {{g,0,0}, 275fe7e6d57SJed Brown {-1.5773502691896257e+00,g,0}, 27625833a80SEmil Constantinescu {-6.7075317547305480e-01,-1.7075317547305482e-01,g}}, 277fe7e6d57SJed Brown b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01}, 278fe7e6d57SJed Brown b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01}; 279f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr); 280fe7e6d57SJed Brown } 281fe7e6d57SJed Brown { 282*3ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 283fe7e6d57SJed Brown const PetscReal g = 4.3586652150845900e-01; 284fe7e6d57SJed Brown const PetscReal 285fe7e6d57SJed Brown A[4][4] = {{0,0,0,0}, 286fe7e6d57SJed Brown {8.7173304301691801e-01,0,0,0}, 287fe7e6d57SJed Brown {8.4457060015369423e-01,-1.1299064236484185e-01,0,0}, 288fe7e6d57SJed Brown {0,0,1.,0}}, 289fe7e6d57SJed Brown Gamma[4][4] = {{g,0,0,0}, 290fe7e6d57SJed Brown {-8.7173304301691801e-01,g,0,0}, 291fe7e6d57SJed Brown {-9.0338057013044082e-01,5.4180672388095326e-02,g,0}, 292fe7e6d57SJed Brown {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}}, 293fe7e6d57SJed Brown b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01}, 294*3ca35412SEmil Constantinescu b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01}; 295*3ca35412SEmil Constantinescu 296*3ca35412SEmil Constantinescu binterpt[0][0]=1.0564298455794094; 297*3ca35412SEmil Constantinescu binterpt[1][0]=2.296429974281067; 298*3ca35412SEmil Constantinescu binterpt[2][0]=-1.307599564525376; 299*3ca35412SEmil Constantinescu binterpt[3][0]=-1.045260255335102; 300*3ca35412SEmil Constantinescu binterpt[0][1]=-1.3864882699759573; 301*3ca35412SEmil Constantinescu binterpt[1][1]=-8.262611700275677; 302*3ca35412SEmil Constantinescu binterpt[2][1]=7.250979895056055; 303*3ca35412SEmil Constantinescu binterpt[3][1]=2.398120075195581; 304*3ca35412SEmil Constantinescu binterpt[0][2]=0.5721822314575016; 305*3ca35412SEmil Constantinescu binterpt[1][2]=4.742931142090097; 306*3ca35412SEmil Constantinescu binterpt[2][2]=-4.398120075195578; 307*3ca35412SEmil Constantinescu binterpt[3][2]=-0.9169932983520199; 308*3ca35412SEmil Constantinescu 309*3ca35412SEmil Constantinescu ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr); 310e27a552bSJed Brown } 311ef3c5b88SJed Brown { 312ef3c5b88SJed Brown const PetscReal g = 0.5; 313ef3c5b88SJed Brown const PetscReal 314ef3c5b88SJed Brown A[4][4] = {{0,0,0,0}, 315ef3c5b88SJed Brown {0,0,0,0}, 316ef3c5b88SJed Brown {1.,0,0,0}, 317ef3c5b88SJed Brown {0.75,-0.25,0.5,0}}, 318ef3c5b88SJed Brown Gamma[4][4] = {{g,0,0,0}, 319ef3c5b88SJed Brown {1.,g,0,0}, 320ef3c5b88SJed Brown {-0.25,-0.25,g,0}, 321ef3c5b88SJed Brown {1./12,1./12,-2./3,g}}, 322ef3c5b88SJed Brown b[4] = {5./6,-1./6,-1./6,0.5}, 323ef3c5b88SJed Brown b2[4] = {0.75,-0.25,0.5,0}; 324f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr); 325ef3c5b88SJed Brown } 326ef3c5b88SJed Brown { 327ef3c5b88SJed Brown const PetscReal g = 0.43586652150845899941601945119356; 328ef3c5b88SJed Brown const PetscReal 329ef3c5b88SJed Brown A[3][3] = {{0,0,0}, 330ef3c5b88SJed Brown {g,0,0}, 331ef3c5b88SJed Brown {g,0,0}}, 332ef3c5b88SJed Brown Gamma[3][3] = {{g,0,0}, 333ef3c5b88SJed Brown {-0.19294655696029095575009695436041,g,0}, 334ef3c5b88SJed Brown {0,1.74927148125794685173529749738960,g}}, 335ef3c5b88SJed Brown b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829}, 336ef3c5b88SJed Brown b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619}; 337f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr); 338ef3c5b88SJed Brown } 339b1c69cc3SEmil Constantinescu { 340aaf9cf16SJed Brown const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0; 341b1c69cc3SEmil Constantinescu const PetscReal 342b1c69cc3SEmil Constantinescu A[3][3] = {{0,0,0}, 343b1c69cc3SEmil Constantinescu {1,0,0}, 344b1c69cc3SEmil Constantinescu {0.25,0.25,0}}, 345b1c69cc3SEmil Constantinescu Gamma[3][3] = {{0,0,0}, 346aaf9cf16SJed Brown {(-3.0-s3)/6.0,g,0}, 347aaf9cf16SJed Brown {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}}, 348b1c69cc3SEmil Constantinescu b[3] = {1./6.,1./6.,2./3.}, 349b1c69cc3SEmil Constantinescu b2[3] = {1./4.,1./4.,1./2.}; 350f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr); 351b1c69cc3SEmil Constantinescu } 352b1c69cc3SEmil Constantinescu 353b1c69cc3SEmil Constantinescu { 354b1c69cc3SEmil Constantinescu const PetscReal 355b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 356b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 357b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 358b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 359b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 360b1c69cc3SEmil Constantinescu {0.0,1./4.,0,0}, 361b1c69cc3SEmil Constantinescu {-2.,-2./3.,2./3.,0}, 362b1c69cc3SEmil Constantinescu {1./2.,5./36.,-2./9,0}}, 363b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 364b1c69cc3SEmil Constantinescu b2[4] = {1./8.,3./4.,1./8.,0}; 365f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr); 366b1c69cc3SEmil Constantinescu } 367b1c69cc3SEmil Constantinescu 368b1c69cc3SEmil Constantinescu { 369b1c69cc3SEmil Constantinescu const PetscReal 370b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 371b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 372b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 373b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 374b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 375b1c69cc3SEmil Constantinescu {0.0,3./4.,0,0}, 376b1c69cc3SEmil Constantinescu {-2./3.,-23./9.,2./9.,0}, 377b1c69cc3SEmil Constantinescu {1./18.,65./108.,-2./27,0}}, 378b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 379b1c69cc3SEmil Constantinescu b2[4] = {3./16.,10./16.,3./16.,0}; 380f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr); 381b1c69cc3SEmil Constantinescu } 382753f8adbSEmil Constantinescu 383753f8adbSEmil Constantinescu { 384753f8adbSEmil Constantinescu PetscReal A[4][4],Gamma[4][4],b[4],b2[4]; 385*3ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 386753f8adbSEmil Constantinescu 387753f8adbSEmil Constantinescu Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816; 38805e8e825SJed Brown Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0; 389753f8adbSEmil Constantinescu Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476; 390753f8adbSEmil Constantinescu Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816; 39105e8e825SJed Brown Gamma[1][2]=0; Gamma[1][3]=0; 392753f8adbSEmil Constantinescu Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903; 393753f8adbSEmil Constantinescu Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131; 394753f8adbSEmil Constantinescu Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816; 39505e8e825SJed Brown Gamma[2][3]=0; 396753f8adbSEmil Constantinescu Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783; 397753f8adbSEmil Constantinescu Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984; 398753f8adbSEmil Constantinescu Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198; 399753f8adbSEmil Constantinescu Gamma[3][3]=0; 400753f8adbSEmil Constantinescu 40105e8e825SJed Brown A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0; 402753f8adbSEmil Constantinescu A[1][0]=0.8717330430169179988320388950590125027645343373957631; 40305e8e825SJed Brown A[1][1]=0; A[1][2]=0; A[1][3]=0; 404753f8adbSEmil Constantinescu A[2][0]=0.5275890119763004115618079766722914408876108660811028; 405753f8adbSEmil Constantinescu A[2][1]=0.07241098802369958843819203208518599088698057726988732; 40605e8e825SJed Brown A[2][2]=0; A[2][3]=0; 407753f8adbSEmil Constantinescu A[3][0]=0.3990960076760701320627260685975778145384666450351314; 408753f8adbSEmil Constantinescu A[3][1]=-0.4375576546135194437228463747348862825846903771419953; 409753f8adbSEmil Constantinescu A[3][2]=1.038461646937449311660120300601880176655352737312713; 41005e8e825SJed Brown A[3][3]=0; 411753f8adbSEmil Constantinescu 412753f8adbSEmil Constantinescu b[0]=0.1876410243467238251612921333138006734899663569186926; 413753f8adbSEmil Constantinescu b[1]=-0.5952974735769549480478230473706443582188442040780541; 414753f8adbSEmil Constantinescu b[2]=0.9717899277217721234705114616271378792182450260943198; 415753f8adbSEmil Constantinescu b[3]=0.4358665215084589994160194475295062513822671686978816; 416753f8adbSEmil Constantinescu 417753f8adbSEmil Constantinescu b2[0]=0.2147402862233891404862383521089097657790734483804460; 418753f8adbSEmil Constantinescu b2[1]=-0.4851622638849390928209050538171743017757490232519684; 419753f8adbSEmil Constantinescu b2[2]=0.8687250025203875511662123688667549217531982787600080; 420753f8adbSEmil Constantinescu b2[3]=0.4016969751411624011684543450940068201770721128357014; 421753f8adbSEmil Constantinescu 422*3ca35412SEmil Constantinescu binterpt[0][0]=2.2565812720167954547104627844105; 423*3ca35412SEmil Constantinescu binterpt[1][0]=1.349166413351089573796243820819; 424*3ca35412SEmil Constantinescu binterpt[2][0]=-2.4695174540533503758652847586647; 425*3ca35412SEmil Constantinescu binterpt[3][0]=-0.13623023131453465264142184656474; 426*3ca35412SEmil Constantinescu binterpt[0][1]=-3.0826699111559187902922463354557; 427*3ca35412SEmil Constantinescu binterpt[1][1]=-2.4689115685996042534544925650515; 428*3ca35412SEmil Constantinescu binterpt[2][1]=5.7428279814696677152129332773553; 429*3ca35412SEmil Constantinescu binterpt[3][1]=-0.19124650171414467146619437684812; 430*3ca35412SEmil Constantinescu binterpt[0][2]=1.0137296634858471607430756831148; 431*3ca35412SEmil Constantinescu binterpt[1][2]=0.52444768167155973161042570784064; 432*3ca35412SEmil Constantinescu binterpt[2][2]=-2.3015205996945452158771370439586; 433*3ca35412SEmil Constantinescu binterpt[3][2]=0.76334325453713832352363565300308; 434f4aed992SEmil Constantinescu 435f73f8d2cSSatish Balay ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr); 436753f8adbSEmil Constantinescu } 437753f8adbSEmil Constantinescu 438e27a552bSJed Brown PetscFunctionReturn(0); 439e27a552bSJed Brown } 440e27a552bSJed Brown 441e27a552bSJed Brown #undef __FUNCT__ 442e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy" 443e27a552bSJed Brown /*@C 444e27a552bSJed Brown TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 445e27a552bSJed Brown 446e27a552bSJed Brown Not Collective 447e27a552bSJed Brown 448e27a552bSJed Brown Level: advanced 449e27a552bSJed Brown 450e27a552bSJed Brown .keywords: TSRosW, register, destroy 451e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 452e27a552bSJed Brown @*/ 453e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void) 454e27a552bSJed Brown { 455e27a552bSJed Brown PetscErrorCode ierr; 45661692a83SJed Brown RosWTableauLink link; 457e27a552bSJed Brown 458e27a552bSJed Brown PetscFunctionBegin; 45961692a83SJed Brown while ((link = RosWTableauList)) { 46061692a83SJed Brown RosWTableau t = &link->tab; 46161692a83SJed Brown RosWTableauList = link->next; 46261692a83SJed Brown ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr); 46343b21953SEmil Constantinescu ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr); 464fe7e6d57SJed Brown ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr); 465f4aed992SEmil Constantinescu ierr = PetscFree(t->binterpt);CHKERRQ(ierr); 466e27a552bSJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 467e27a552bSJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 468e27a552bSJed Brown } 469e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_FALSE; 470e27a552bSJed Brown PetscFunctionReturn(0); 471e27a552bSJed Brown } 472e27a552bSJed Brown 473e27a552bSJed Brown #undef __FUNCT__ 474e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage" 475e27a552bSJed Brown /*@C 476e27a552bSJed Brown TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called 477e27a552bSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW() 478e27a552bSJed Brown when using static libraries. 479e27a552bSJed Brown 480e27a552bSJed Brown Input Parameter: 481e27a552bSJed Brown path - The dynamic library path, or PETSC_NULL 482e27a552bSJed Brown 483e27a552bSJed Brown Level: developer 484e27a552bSJed Brown 485e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package 486e27a552bSJed Brown .seealso: PetscInitialize() 487e27a552bSJed Brown @*/ 488e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[]) 489e27a552bSJed Brown { 490e27a552bSJed Brown PetscErrorCode ierr; 491e27a552bSJed Brown 492e27a552bSJed Brown PetscFunctionBegin; 493e27a552bSJed Brown if (TSRosWPackageInitialized) PetscFunctionReturn(0); 494e27a552bSJed Brown TSRosWPackageInitialized = PETSC_TRUE; 495e27a552bSJed Brown ierr = TSRosWRegisterAll();CHKERRQ(ierr); 496e27a552bSJed Brown ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr); 497e27a552bSJed Brown PetscFunctionReturn(0); 498e27a552bSJed Brown } 499e27a552bSJed Brown 500e27a552bSJed Brown #undef __FUNCT__ 501e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage" 502e27a552bSJed Brown /*@C 503e27a552bSJed Brown TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is 504e27a552bSJed Brown called from PetscFinalize(). 505e27a552bSJed Brown 506e27a552bSJed Brown Level: developer 507e27a552bSJed Brown 508e27a552bSJed Brown .keywords: Petsc, destroy, package 509e27a552bSJed Brown .seealso: PetscFinalize() 510e27a552bSJed Brown @*/ 511e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void) 512e27a552bSJed Brown { 513e27a552bSJed Brown PetscErrorCode ierr; 514e27a552bSJed Brown 515e27a552bSJed Brown PetscFunctionBegin; 516e27a552bSJed Brown TSRosWPackageInitialized = PETSC_FALSE; 517e27a552bSJed Brown ierr = TSRosWRegisterDestroy();CHKERRQ(ierr); 518e27a552bSJed Brown PetscFunctionReturn(0); 519e27a552bSJed Brown } 520e27a552bSJed Brown 521e27a552bSJed Brown #undef __FUNCT__ 522e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister" 523e27a552bSJed Brown /*@C 52461692a83SJed Brown TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 525e27a552bSJed Brown 526e27a552bSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 527e27a552bSJed Brown 528e27a552bSJed Brown Input Parameters: 529e27a552bSJed Brown + name - identifier for method 530e27a552bSJed Brown . order - approximation order of method 531e27a552bSJed Brown . s - number of stages, this is the dimension of the matrices below 53261692a83SJed Brown . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 53361692a83SJed Brown . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 534fe7e6d57SJed Brown . b - Step completion table (dimension s) 535fe7e6d57SJed Brown - bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available) 536f4aed992SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt 537f4aed992SEmil Constantinescu . binterpt - Coefficients of the interpolation formula (dimension s*pinterp) 538e27a552bSJed Brown 539e27a552bSJed Brown Notes: 54061692a83SJed Brown Several Rosenbrock W methods are provided, this function is only needed to create new methods. 541e27a552bSJed Brown 542e27a552bSJed Brown Level: advanced 543e27a552bSJed Brown 544e27a552bSJed Brown .keywords: TS, register 545e27a552bSJed Brown 546e27a552bSJed Brown .seealso: TSRosW 547e27a552bSJed Brown @*/ 548e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s, 549f4aed992SEmil Constantinescu const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[], 550f4aed992SEmil Constantinescu PetscInt pinterp,const PetscReal binterpt[]) 551e27a552bSJed Brown { 552e27a552bSJed Brown PetscErrorCode ierr; 55361692a83SJed Brown RosWTableauLink link; 55461692a83SJed Brown RosWTableau t; 55561692a83SJed Brown PetscInt i,j,k; 55661692a83SJed Brown PetscScalar *GammaInv; 557e27a552bSJed Brown 558e27a552bSJed Brown PetscFunctionBegin; 559fe7e6d57SJed Brown PetscValidCharPointer(name,1); 560fe7e6d57SJed Brown PetscValidPointer(A,4); 561fe7e6d57SJed Brown PetscValidPointer(Gamma,5); 562fe7e6d57SJed Brown PetscValidPointer(b,6); 563fe7e6d57SJed Brown if (bembed) PetscValidPointer(bembed,7); 564fe7e6d57SJed Brown 565e27a552bSJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 566e27a552bSJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 567e27a552bSJed Brown t = &link->tab; 568e27a552bSJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 569e27a552bSJed Brown t->order = order; 570e27a552bSJed Brown t->s = s; 57161692a83SJed Brown ierr = PetscMalloc5(s*s,PetscReal,&t->A,s*s,PetscReal,&t->Gamma,s,PetscReal,&t->b,s,PetscReal,&t->ASum,s,PetscReal,&t->GammaSum);CHKERRQ(ierr); 57243b21953SEmil Constantinescu ierr = PetscMalloc5(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag,s*s,PetscReal,&t->GammaExplicitCorr);CHKERRQ(ierr); 573e27a552bSJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 57461692a83SJed Brown ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 57543b21953SEmil Constantinescu ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 57661692a83SJed Brown ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 577fe7e6d57SJed Brown if (bembed) { 578fe7e6d57SJed Brown ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr); 579fe7e6d57SJed Brown ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 580fe7e6d57SJed Brown } 58161692a83SJed Brown for (i=0; i<s; i++) { 58261692a83SJed Brown t->ASum[i] = 0; 58361692a83SJed Brown t->GammaSum[i] = 0; 58461692a83SJed Brown for (j=0; j<s; j++) { 58561692a83SJed Brown t->ASum[i] += A[i*s+j]; 586fe7e6d57SJed Brown t->GammaSum[i] += Gamma[i*s+j]; 58761692a83SJed Brown } 58861692a83SJed Brown } 58961692a83SJed Brown ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */ 59061692a83SJed Brown for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i]; 591fd96d5b0SEmil Constantinescu for (i=0; i<s; i++) { 592fd96d5b0SEmil Constantinescu if (Gamma[i*s+i] == 0.0) { 593fd96d5b0SEmil Constantinescu GammaInv[i*s+i] = 1.0; 594c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_TRUE; 595fd96d5b0SEmil Constantinescu } else { 596c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_FALSE; 597fd96d5b0SEmil Constantinescu } 598fd96d5b0SEmil Constantinescu } 599fd96d5b0SEmil Constantinescu 60061692a83SJed Brown switch (s) { 60161692a83SJed Brown case 1: GammaInv[0] = 1./GammaInv[0]; break; 60261692a83SJed Brown case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break; 60361692a83SJed Brown case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break; 60461692a83SJed Brown case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break; 60561692a83SJed Brown case 5: { 60661692a83SJed Brown PetscInt ipvt5[5]; 60761692a83SJed Brown MatScalar work5[5*5]; 60861692a83SJed Brown ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break; 60961692a83SJed Brown } 61061692a83SJed Brown case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break; 61161692a83SJed Brown case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break; 61261692a83SJed Brown default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s); 61361692a83SJed Brown } 61461692a83SJed Brown for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 61561692a83SJed Brown ierr = PetscFree(GammaInv);CHKERRQ(ierr); 61643b21953SEmil Constantinescu 61743b21953SEmil Constantinescu for (i=0; i<s; i++) { 61843b21953SEmil Constantinescu for (k=0; k<i+1; k++) { 61943b21953SEmil Constantinescu t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]); 62043b21953SEmil Constantinescu for (j=k+1; j<i+1; j++) { 62143b21953SEmil Constantinescu t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]); 62243b21953SEmil Constantinescu } 62343b21953SEmil Constantinescu } 62443b21953SEmil Constantinescu } 62543b21953SEmil Constantinescu 62661692a83SJed Brown for (i=0; i<s; i++) { 62761692a83SJed Brown for (j=0; j<s; j++) { 62861692a83SJed Brown t->At[i*s+j] = 0; 62961692a83SJed Brown for (k=0; k<s; k++) { 63061692a83SJed Brown t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j]; 63161692a83SJed Brown } 63261692a83SJed Brown } 63361692a83SJed Brown t->bt[i] = 0; 63461692a83SJed Brown for (j=0; j<s; j++) { 63561692a83SJed Brown t->bt[i] += t->b[j] * t->GammaInv[j*s+i]; 63661692a83SJed Brown } 637fe7e6d57SJed Brown if (bembed) { 638fe7e6d57SJed Brown t->bembedt[i] = 0; 639fe7e6d57SJed Brown for (j=0; j<s; j++) { 640fe7e6d57SJed Brown t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i]; 641fe7e6d57SJed Brown } 642fe7e6d57SJed Brown } 64361692a83SJed Brown } 6448d59e960SJed Brown t->ccfl = 1.0; /* Fix this */ 6458d59e960SJed Brown 646f4aed992SEmil Constantinescu t->pinterp = pinterp; 647*3ca35412SEmil Constantinescu ierr = PetscMalloc(s*pinterp*sizeof(binterpt[0]),&t->binterpt);CHKERRQ(ierr); 648*3ca35412SEmil Constantinescu ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 64961692a83SJed Brown link->next = RosWTableauList; 65061692a83SJed Brown RosWTableauList = link; 651e27a552bSJed Brown PetscFunctionReturn(0); 652e27a552bSJed Brown } 653e27a552bSJed Brown 654e27a552bSJed Brown #undef __FUNCT__ 6551c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW" 6561c3436cfSJed Brown /* 6571c3436cfSJed Brown The step completion formula is 6581c3436cfSJed Brown 6591c3436cfSJed Brown x1 = x0 + b^T Y 6601c3436cfSJed Brown 6611c3436cfSJed Brown where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 6621c3436cfSJed Brown updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 6631c3436cfSJed Brown 6641c3436cfSJed Brown x1e = x0 + be^T Y 6651c3436cfSJed Brown = x1 - b^T Y + be^T Y 6661c3436cfSJed Brown = x1 + (be - b)^T Y 6671c3436cfSJed Brown 6681c3436cfSJed Brown so we can evaluate the method of different order even after the step has been optimistically completed. 6691c3436cfSJed Brown */ 6701c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done) 6711c3436cfSJed Brown { 6721c3436cfSJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 6731c3436cfSJed Brown RosWTableau tab = ros->tableau; 6741c3436cfSJed Brown PetscScalar *w = ros->work; 6751c3436cfSJed Brown PetscInt i; 6761c3436cfSJed Brown PetscErrorCode ierr; 6771c3436cfSJed Brown 6781c3436cfSJed Brown PetscFunctionBegin; 6791c3436cfSJed Brown if (order == tab->order) { 680108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */ 6811c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 682de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bt[i]; 683de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 684108c343cSJed Brown } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 6851c3436cfSJed Brown if (done) *done = PETSC_TRUE; 6861c3436cfSJed Brown PetscFunctionReturn(0); 6871c3436cfSJed Brown } else if (order == tab->order-1) { 6881c3436cfSJed Brown if (!tab->bembedt) goto unavailable; 689108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */ 6901c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 691de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i]; 692de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 693108c343cSJed Brown } else { /* Use rollback-and-recomplete formula (bembedt - bt) */ 694108c343cSJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 695108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 696108c343cSJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 6971c3436cfSJed Brown } 6981c3436cfSJed Brown if (done) *done = PETSC_TRUE; 6991c3436cfSJed Brown PetscFunctionReturn(0); 7001c3436cfSJed Brown } 7011c3436cfSJed Brown unavailable: 7021c3436cfSJed Brown if (done) *done = PETSC_FALSE; 7031c3436cfSJed Brown else SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_SUP,"Rosenbrock-W '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 7041c3436cfSJed Brown PetscFunctionReturn(0); 7051c3436cfSJed Brown } 7061c3436cfSJed Brown 7071c3436cfSJed Brown #undef __FUNCT__ 708e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW" 709e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts) 710e27a552bSJed Brown { 71161692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 71261692a83SJed Brown RosWTableau tab = ros->tableau; 713e27a552bSJed Brown const PetscInt s = tab->s; 7141c3436cfSJed Brown const PetscReal *At = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv; 7150feba352SEmil Constantinescu const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr; 716c17803e7SJed Brown const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 71761692a83SJed Brown PetscScalar *w = ros->work; 7187d4bf2deSEmil Constantinescu Vec *Y = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage; 719e27a552bSJed Brown SNES snes; 7201c3436cfSJed Brown TSAdapt adapt; 7211c3436cfSJed Brown PetscInt i,j,its,lits,reject,next_scheme; 722cdbf8f93SLisandro Dalcin PetscReal next_time_step; 7231c3436cfSJed Brown PetscBool accept; 724e27a552bSJed Brown PetscErrorCode ierr; 7250feba352SEmil Constantinescu MatStructure str; 726e27a552bSJed Brown 727e27a552bSJed Brown PetscFunctionBegin; 728e27a552bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 729cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 7301c3436cfSJed Brown accept = PETSC_TRUE; 731108c343cSJed Brown ros->status = TS_STEP_INCOMPLETE; 732e27a552bSJed Brown 73397335746SJed Brown for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 7341c3436cfSJed Brown const PetscReal h = ts->time_step; 735*3ca35412SEmil Constantinescu ierr = VecCopy(ts->vec_sol,ros->VecSolPrev);CHKERRQ(ierr);/*move this at the end*/ 736e27a552bSJed Brown for (i=0; i<s; i++) { 7371c3436cfSJed Brown ros->stage_time = ts->ptime + h*ASum[i]; 738c17803e7SJed Brown if (GammaZeroDiag[i]) { 739c17803e7SJed Brown ros->stage_explicit = PETSC_TRUE; 740fd96d5b0SEmil Constantinescu ros->shift = 1./h; 741c17803e7SJed Brown } else { 742c17803e7SJed Brown ros->stage_explicit = PETSC_FALSE; 74361692a83SJed Brown ros->shift = 1./(h*Gamma[i*s+i]); 744fd96d5b0SEmil Constantinescu } 74561692a83SJed Brown 74661692a83SJed Brown ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr); 747de19f811SJed Brown for (j=0; j<i; j++) w[j] = At[i*s+j]; 748de19f811SJed Brown ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 74961692a83SJed Brown 75061692a83SJed Brown for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j]; 75161692a83SJed Brown ierr = VecZeroEntries(Zdot);CHKERRQ(ierr); 75261692a83SJed Brown ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr); 75361692a83SJed Brown 754e27a552bSJed Brown /* Initial guess taken from last stage */ 75561692a83SJed Brown ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr); 75661692a83SJed Brown 7577d4bf2deSEmil Constantinescu if (!ros->stage_explicit) { 75861692a83SJed Brown if (!ros->recompute_jacobian && !i) { 75961692a83SJed Brown ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 76061692a83SJed Brown } 76161692a83SJed Brown ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr); 762e27a552bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 763e27a552bSJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 764e27a552bSJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 76597335746SJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 76697335746SJed Brown ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 76797335746SJed Brown if (!accept) goto reject_step; 7687d4bf2deSEmil Constantinescu } else { 7691ce71dffSSatish Balay Mat J,Jp; 7700feba352SEmil Constantinescu ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */ 7710feba352SEmil Constantinescu ierr = TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);CHKERRQ(ierr); 7720feba352SEmil Constantinescu ierr = VecScale(Y[i],-1.0); 7730feba352SEmil Constantinescu ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/ 7740feba352SEmil Constantinescu 7750feba352SEmil Constantinescu ierr = VecZeroEntries(Zstage);CHKERRQ(ierr); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */ 7760feba352SEmil Constantinescu for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j]; 7770feba352SEmil Constantinescu ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 7780feba352SEmil Constantinescu /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */ 7790feba352SEmil Constantinescu str = SAME_NONZERO_PATTERN; 7800feba352SEmil Constantinescu ierr = TSGetIJacobian(ts,&J,&Jp,PETSC_NULL,PETSC_NULL); 7810feba352SEmil Constantinescu ierr = TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);CHKERRQ(ierr); 7820feba352SEmil Constantinescu ierr = MatMult(J,Zstage,Zdot); 7830feba352SEmil Constantinescu 7840feba352SEmil Constantinescu ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); 7850feba352SEmil Constantinescu ierr = VecScale(Y[i],h); 7867d4bf2deSEmil Constantinescu ts->linear_its += 1; 7877d4bf2deSEmil Constantinescu } 788e27a552bSJed Brown } 7891c3436cfSJed Brown ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 790108c343cSJed Brown ros->status = TS_STEP_PENDING; 791e27a552bSJed Brown 7921c3436cfSJed Brown /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 7931c3436cfSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 7941c3436cfSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 7958d59e960SJed Brown ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 7961c3436cfSJed Brown ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 7971c3436cfSJed Brown if (accept) { 7981c3436cfSJed Brown /* ignore next_scheme for now */ 799e27a552bSJed Brown ts->ptime += ts->time_step; 800cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 801e27a552bSJed Brown ts->steps++; 802108c343cSJed Brown ros->status = TS_STEP_COMPLETE; 8031c3436cfSJed Brown break; 8041c3436cfSJed Brown } else { /* Roll back the current step */ 8051c3436cfSJed Brown for (i=0; i<s; i++) w[i] = -tab->bt[i]; 8061c3436cfSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr); 8071c3436cfSJed Brown ts->time_step = next_time_step; 808108c343cSJed Brown ros->status = TS_STEP_INCOMPLETE; 8091c3436cfSJed Brown } 810476b6736SJed Brown reject_step: continue; 8111c3436cfSJed Brown } 812b2ce242eSJed Brown if (ros->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 813e27a552bSJed Brown PetscFunctionReturn(0); 814e27a552bSJed Brown } 815e27a552bSJed Brown 816e27a552bSJed Brown #undef __FUNCT__ 817e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW" 818e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X) 819e27a552bSJed Brown { 82061692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 821f4aed992SEmil Constantinescu PetscInt s = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j; 822f4aed992SEmil Constantinescu PetscReal h; 823f4aed992SEmil Constantinescu PetscReal tt,t; 824f4aed992SEmil Constantinescu PetscScalar *bt; 825f4aed992SEmil Constantinescu const PetscReal *Bt = ros->tableau->binterpt; 826f4aed992SEmil Constantinescu PetscErrorCode ierr; 827f4aed992SEmil Constantinescu const PetscReal *GammaInv = ros->tableau->GammaInv; 828f4aed992SEmil Constantinescu PetscScalar *w = ros->work; 829f4aed992SEmil Constantinescu Vec *Y = ros->Y; 830e27a552bSJed Brown 831e27a552bSJed Brown PetscFunctionBegin; 832f4aed992SEmil Constantinescu if (!Bt) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name); 833f4aed992SEmil Constantinescu 834f4aed992SEmil Constantinescu switch (ros->status) { 835f4aed992SEmil Constantinescu case TS_STEP_INCOMPLETE: 836f4aed992SEmil Constantinescu case TS_STEP_PENDING: 837f4aed992SEmil Constantinescu h = ts->time_step; 838f4aed992SEmil Constantinescu t = (itime - ts->ptime)/h; 839f4aed992SEmil Constantinescu break; 840f4aed992SEmil Constantinescu case TS_STEP_COMPLETE: 841f4aed992SEmil Constantinescu h = ts->time_step_prev; 842f4aed992SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 843f4aed992SEmil Constantinescu break; 844f4aed992SEmil Constantinescu default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus"); 845f4aed992SEmil Constantinescu } 846*3ca35412SEmil Constantinescu ierr = PetscMalloc(s*sizeof(bt[0]),&bt);CHKERRQ(ierr); 847f4aed992SEmil Constantinescu for (i=0; i<s; i++) bt[i] = 0; 848f4aed992SEmil Constantinescu for (j=0,tt=t; j<pinterp; j++,tt*=t) { 849f4aed992SEmil Constantinescu for (i=0; i<s; i++) { 850*3ca35412SEmil Constantinescu bt[i] += Bt[i*pinterp+j] * tt; 851f4aed992SEmil Constantinescu } 852f4aed992SEmil Constantinescu } 853f4aed992SEmil Constantinescu 854f4aed992SEmil Constantinescu /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */ 855f4aed992SEmil Constantinescu /*X<-0*/ 856f4aed992SEmil Constantinescu ierr = VecZeroEntries(X);CHKERRQ(ierr); 857f4aed992SEmil Constantinescu 858f4aed992SEmil Constantinescu /*X<- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */ 859*3ca35412SEmil Constantinescu for (j=0; j<s; j++) w[j]=0; 860*3ca35412SEmil Constantinescu for (j=0; j<s; j++) { 861*3ca35412SEmil Constantinescu for (i=j; i<s; i++) { 862*3ca35412SEmil Constantinescu w[j] += bt[i]*GammaInv[i*s+j]; 863f4aed992SEmil Constantinescu } 864*3ca35412SEmil Constantinescu } 865*3ca35412SEmil Constantinescu ierr = VecMAXPY(X,i,w,Y);CHKERRQ(ierr); 866f4aed992SEmil Constantinescu 867f4aed992SEmil Constantinescu /*X<-y(t) + X*/ 868*3ca35412SEmil Constantinescu ierr = VecAXPY(X,1.0,ros->VecSolPrev);CHKERRQ(ierr); 869f4aed992SEmil Constantinescu 870f4aed992SEmil Constantinescu ierr = PetscFree(bt);CHKERRQ(ierr); 871f4aed992SEmil Constantinescu 872e27a552bSJed Brown PetscFunctionReturn(0); 873e27a552bSJed Brown } 874e27a552bSJed Brown 875e27a552bSJed Brown /*------------------------------------------------------------*/ 876e27a552bSJed Brown #undef __FUNCT__ 877e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW" 878e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts) 879e27a552bSJed Brown { 88061692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 881e27a552bSJed Brown PetscInt s; 882e27a552bSJed Brown PetscErrorCode ierr; 883e27a552bSJed Brown 884e27a552bSJed Brown PetscFunctionBegin; 88561692a83SJed Brown if (!ros->tableau) PetscFunctionReturn(0); 88661692a83SJed Brown s = ros->tableau->s; 88761692a83SJed Brown ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr); 88861692a83SJed Brown ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr); 88961692a83SJed Brown ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr); 89061692a83SJed Brown ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr); 89161692a83SJed Brown ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr); 892*3ca35412SEmil Constantinescu ierr = VecDestroy(&ros->VecSolPrev);CHKERRQ(ierr); 89361692a83SJed Brown ierr = PetscFree(ros->work);CHKERRQ(ierr); 894e27a552bSJed Brown PetscFunctionReturn(0); 895e27a552bSJed Brown } 896e27a552bSJed Brown 897e27a552bSJed Brown #undef __FUNCT__ 898e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW" 899e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts) 900e27a552bSJed Brown { 901e27a552bSJed Brown PetscErrorCode ierr; 902e27a552bSJed Brown 903e27a552bSJed Brown PetscFunctionBegin; 904e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 905e27a552bSJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 906e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr); 907e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr); 90861692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr); 909e27a552bSJed Brown PetscFunctionReturn(0); 910e27a552bSJed Brown } 911e27a552bSJed Brown 912e27a552bSJed Brown /* 913e27a552bSJed Brown This defines the nonlinear equation that is to be solved with SNES 914e27a552bSJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 915e27a552bSJed Brown */ 916e27a552bSJed Brown #undef __FUNCT__ 917e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW" 918e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts) 919e27a552bSJed Brown { 92061692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 921e27a552bSJed Brown PetscErrorCode ierr; 922e27a552bSJed Brown 923e27a552bSJed Brown PetscFunctionBegin; 92461692a83SJed Brown ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */ 92561692a83SJed Brown ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr); /* Ystage = X + Zstage */ 92661692a83SJed Brown ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr); 927e27a552bSJed Brown PetscFunctionReturn(0); 928e27a552bSJed Brown } 929e27a552bSJed Brown 930e27a552bSJed Brown #undef __FUNCT__ 931e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW" 932e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 933e27a552bSJed Brown { 93461692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 935e27a552bSJed Brown PetscErrorCode ierr; 936e27a552bSJed Brown 937e27a552bSJed Brown PetscFunctionBegin; 93861692a83SJed Brown /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 93961692a83SJed Brown ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 940e27a552bSJed Brown PetscFunctionReturn(0); 941e27a552bSJed Brown } 942e27a552bSJed Brown 943e27a552bSJed Brown #undef __FUNCT__ 944e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW" 945e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts) 946e27a552bSJed Brown { 94761692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 94861692a83SJed Brown RosWTableau tab = ros->tableau; 949e27a552bSJed Brown PetscInt s = tab->s; 950e27a552bSJed Brown PetscErrorCode ierr; 951e27a552bSJed Brown 952e27a552bSJed Brown PetscFunctionBegin; 95361692a83SJed Brown if (!ros->tableau) { 954e27a552bSJed Brown ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr); 955e27a552bSJed Brown } 95661692a83SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr); 95761692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr); 95861692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr); 95961692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr); 96061692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr); 961*3ca35412SEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&ros->VecSolPrev);CHKERRQ(ierr); 96261692a83SJed Brown ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr); 963e27a552bSJed Brown PetscFunctionReturn(0); 964e27a552bSJed Brown } 965e27a552bSJed Brown /*------------------------------------------------------------*/ 966e27a552bSJed Brown 967e27a552bSJed Brown #undef __FUNCT__ 968e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW" 969e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts) 970e27a552bSJed Brown { 97161692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 972e27a552bSJed Brown PetscErrorCode ierr; 97361692a83SJed Brown char rostype[256]; 974e27a552bSJed Brown 975e27a552bSJed Brown PetscFunctionBegin; 976e27a552bSJed Brown ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr); 977e27a552bSJed Brown { 97861692a83SJed Brown RosWTableauLink link; 979e27a552bSJed Brown PetscInt count,choice; 980e27a552bSJed Brown PetscBool flg; 981e27a552bSJed Brown const char **namelist; 98261692a83SJed Brown SNES snes; 98361692a83SJed Brown 98461692a83SJed Brown ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr); 98561692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) ; 986e27a552bSJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 98761692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 98861692a83SJed Brown ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr); 98961692a83SJed Brown ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr); 990e27a552bSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 99161692a83SJed Brown 99261692a83SJed Brown ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr); 99361692a83SJed Brown 99461692a83SJed Brown /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 99561692a83SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 99661692a83SJed Brown if (!((PetscObject)snes)->type_name) { 99761692a83SJed Brown ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 99861692a83SJed Brown } 99961692a83SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1000e27a552bSJed Brown } 1001e27a552bSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1002e27a552bSJed Brown PetscFunctionReturn(0); 1003e27a552bSJed Brown } 1004e27a552bSJed Brown 1005e27a552bSJed Brown #undef __FUNCT__ 1006e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray" 1007e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 1008e27a552bSJed Brown { 1009e27a552bSJed Brown PetscErrorCode ierr; 1010e408995aSJed Brown PetscInt i; 1011e408995aSJed Brown size_t left,count; 1012e27a552bSJed Brown char *p; 1013e27a552bSJed Brown 1014e27a552bSJed Brown PetscFunctionBegin; 1015e408995aSJed Brown for (i=0,p=buf,left=len; i<n; i++) { 1016e408995aSJed Brown ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 1017e27a552bSJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 1018e27a552bSJed Brown left -= count; 1019e27a552bSJed Brown p += count; 1020e27a552bSJed Brown *p++ = ' '; 1021e27a552bSJed Brown } 1022e27a552bSJed Brown p[i ? 0 : -1] = 0; 1023e27a552bSJed Brown PetscFunctionReturn(0); 1024e27a552bSJed Brown } 1025e27a552bSJed Brown 1026e27a552bSJed Brown #undef __FUNCT__ 1027e27a552bSJed Brown #define __FUNCT__ "TSView_RosW" 1028e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer) 1029e27a552bSJed Brown { 103061692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 103161692a83SJed Brown RosWTableau tab = ros->tableau; 1032e27a552bSJed Brown PetscBool iascii; 1033e27a552bSJed Brown PetscErrorCode ierr; 1034e27a552bSJed Brown 1035e27a552bSJed Brown PetscFunctionBegin; 1036e27a552bSJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1037e27a552bSJed Brown if (iascii) { 103861692a83SJed Brown const TSRosWType rostype; 1039e408995aSJed Brown PetscInt i; 1040e408995aSJed Brown PetscReal abscissa[512]; 1041e27a552bSJed Brown char buf[512]; 104261692a83SJed Brown ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr); 104361692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Rosenbrock-W %s\n",rostype);CHKERRQ(ierr); 1044e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr); 104561692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A = %s\n",buf);CHKERRQ(ierr); 1046e408995aSJed Brown for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 1047e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr); 1048e408995aSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr); 1049e27a552bSJed Brown } 1050e27a552bSJed Brown ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 1051e27a552bSJed Brown PetscFunctionReturn(0); 1052e27a552bSJed Brown } 1053e27a552bSJed Brown 1054e27a552bSJed Brown #undef __FUNCT__ 1055e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType" 1056e27a552bSJed Brown /*@C 105761692a83SJed Brown TSRosWSetType - Set the type of Rosenbrock-W scheme 1058e27a552bSJed Brown 1059e27a552bSJed Brown Logically collective 1060e27a552bSJed Brown 1061e27a552bSJed Brown Input Parameter: 1062e27a552bSJed Brown + ts - timestepping context 106361692a83SJed Brown - rostype - type of Rosenbrock-W scheme 1064e27a552bSJed Brown 1065020d8f30SJed Brown Level: beginner 1066e27a552bSJed Brown 1067020d8f30SJed Brown .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3 1068e27a552bSJed Brown @*/ 106961692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype) 1070e27a552bSJed Brown { 1071e27a552bSJed Brown PetscErrorCode ierr; 1072e27a552bSJed Brown 1073e27a552bSJed Brown PetscFunctionBegin; 1074e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 107561692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr); 1076e27a552bSJed Brown PetscFunctionReturn(0); 1077e27a552bSJed Brown } 1078e27a552bSJed Brown 1079e27a552bSJed Brown #undef __FUNCT__ 1080e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType" 1081e27a552bSJed Brown /*@C 108261692a83SJed Brown TSRosWGetType - Get the type of Rosenbrock-W scheme 1083e27a552bSJed Brown 1084e27a552bSJed Brown Logically collective 1085e27a552bSJed Brown 1086e27a552bSJed Brown Input Parameter: 1087e27a552bSJed Brown . ts - timestepping context 1088e27a552bSJed Brown 1089e27a552bSJed Brown Output Parameter: 109061692a83SJed Brown . rostype - type of Rosenbrock-W scheme 1091e27a552bSJed Brown 1092e27a552bSJed Brown Level: intermediate 1093e27a552bSJed Brown 1094e27a552bSJed Brown .seealso: TSRosWGetType() 1095e27a552bSJed Brown @*/ 109661692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype) 1097e27a552bSJed Brown { 1098e27a552bSJed Brown PetscErrorCode ierr; 1099e27a552bSJed Brown 1100e27a552bSJed Brown PetscFunctionBegin; 1101e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 110261692a83SJed Brown ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr); 1103e27a552bSJed Brown PetscFunctionReturn(0); 1104e27a552bSJed Brown } 1105e27a552bSJed Brown 1106e27a552bSJed Brown #undef __FUNCT__ 110761692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian" 1108e27a552bSJed Brown /*@C 110961692a83SJed Brown TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 1110e27a552bSJed Brown 1111e27a552bSJed Brown Logically collective 1112e27a552bSJed Brown 1113e27a552bSJed Brown Input Parameter: 1114e27a552bSJed Brown + ts - timestepping context 111561692a83SJed Brown - flg - PETSC_TRUE to recompute the Jacobian at each stage 1116e27a552bSJed Brown 1117e27a552bSJed Brown Level: intermediate 1118e27a552bSJed Brown 1119e27a552bSJed Brown .seealso: TSRosWGetType() 1120e27a552bSJed Brown @*/ 112161692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg) 1122e27a552bSJed Brown { 1123e27a552bSJed Brown PetscErrorCode ierr; 1124e27a552bSJed Brown 1125e27a552bSJed Brown PetscFunctionBegin; 1126e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 112761692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1128e27a552bSJed Brown PetscFunctionReturn(0); 1129e27a552bSJed Brown } 1130e27a552bSJed Brown 1131e27a552bSJed Brown EXTERN_C_BEGIN 1132e27a552bSJed Brown #undef __FUNCT__ 1133e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW" 113461692a83SJed Brown PetscErrorCode TSRosWGetType_RosW(TS ts,const TSRosWType *rostype) 1135e27a552bSJed Brown { 113661692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1137e27a552bSJed Brown PetscErrorCode ierr; 1138e27a552bSJed Brown 1139e27a552bSJed Brown PetscFunctionBegin; 114061692a83SJed Brown if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);} 114161692a83SJed Brown *rostype = ros->tableau->name; 1142e27a552bSJed Brown PetscFunctionReturn(0); 1143e27a552bSJed Brown } 1144e27a552bSJed Brown #undef __FUNCT__ 1145e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW" 114661692a83SJed Brown PetscErrorCode TSRosWSetType_RosW(TS ts,const TSRosWType rostype) 1147e27a552bSJed Brown { 114861692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1149e27a552bSJed Brown PetscErrorCode ierr; 1150e27a552bSJed Brown PetscBool match; 115161692a83SJed Brown RosWTableauLink link; 1152e27a552bSJed Brown 1153e27a552bSJed Brown PetscFunctionBegin; 115461692a83SJed Brown if (ros->tableau) { 115561692a83SJed Brown ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr); 1156e27a552bSJed Brown if (match) PetscFunctionReturn(0); 1157e27a552bSJed Brown } 115861692a83SJed Brown for (link = RosWTableauList; link; link=link->next) { 115961692a83SJed Brown ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr); 1160e27a552bSJed Brown if (match) { 1161e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 116261692a83SJed Brown ros->tableau = &link->tab; 1163e27a552bSJed Brown PetscFunctionReturn(0); 1164e27a552bSJed Brown } 1165e27a552bSJed Brown } 116661692a83SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype); 1167e27a552bSJed Brown PetscFunctionReturn(0); 1168e27a552bSJed Brown } 116961692a83SJed Brown 1170e27a552bSJed Brown #undef __FUNCT__ 117161692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW" 117261692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg) 1173e27a552bSJed Brown { 117461692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1175e27a552bSJed Brown 1176e27a552bSJed Brown PetscFunctionBegin; 117761692a83SJed Brown ros->recompute_jacobian = flg; 1178e27a552bSJed Brown PetscFunctionReturn(0); 1179e27a552bSJed Brown } 1180e27a552bSJed Brown EXTERN_C_END 1181e27a552bSJed Brown 1182e27a552bSJed Brown /* ------------------------------------------------------------ */ 1183e27a552bSJed Brown /*MC 1184020d8f30SJed Brown TSROSW - ODE solver using Rosenbrock-W schemes 1185e27a552bSJed Brown 1186e27a552bSJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1187e27a552bSJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1188e27a552bSJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1189e27a552bSJed Brown 1190e27a552bSJed Brown Notes: 119161692a83SJed Brown This method currently only works with autonomous ODE and DAE. 119261692a83SJed Brown 119361692a83SJed Brown Developer notes: 119461692a83SJed Brown Rosenbrock-W methods are typically specified for autonomous ODE 119561692a83SJed Brown 119661692a83SJed Brown $ xdot = f(x) 119761692a83SJed Brown 119861692a83SJed Brown by the stage equations 119961692a83SJed Brown 120061692a83SJed Brown $ k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 120161692a83SJed Brown 120261692a83SJed Brown and step completion formula 120361692a83SJed Brown 120461692a83SJed Brown $ x_1 = x_0 + sum_j b_j k_j 120561692a83SJed Brown 120661692a83SJed Brown with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x) 120761692a83SJed Brown and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 120861692a83SJed Brown we define new variables for the stage equations 120961692a83SJed Brown 121061692a83SJed Brown $ y_i = gamma_ij k_j 121161692a83SJed Brown 121261692a83SJed Brown The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 121361692a83SJed Brown 121461692a83SJed Brown $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i} 121561692a83SJed Brown 121661692a83SJed Brown to rewrite the method as 121761692a83SJed Brown 121861692a83SJed Brown $ [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j 121961692a83SJed Brown $ x_1 = x_0 + sum_j bt_j y_j 122061692a83SJed Brown 122161692a83SJed Brown where we have introduced the mass matrix M. Continue by defining 122261692a83SJed Brown 122361692a83SJed Brown $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 122461692a83SJed Brown 122561692a83SJed Brown or, more compactly in tensor notation 122661692a83SJed Brown 122761692a83SJed Brown $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 122861692a83SJed Brown 122961692a83SJed Brown Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 123061692a83SJed Brown stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 123161692a83SJed Brown equation 123261692a83SJed Brown 123361692a83SJed Brown $ g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 123461692a83SJed Brown 123561692a83SJed Brown with initial guess y_i = 0. 1236e27a552bSJed Brown 1237e27a552bSJed Brown Level: beginner 1238e27a552bSJed Brown 1239020d8f30SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister() 1240e27a552bSJed Brown 1241e27a552bSJed Brown M*/ 1242e27a552bSJed Brown EXTERN_C_BEGIN 1243e27a552bSJed Brown #undef __FUNCT__ 1244e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW" 1245e27a552bSJed Brown PetscErrorCode TSCreate_RosW(TS ts) 1246e27a552bSJed Brown { 124761692a83SJed Brown TS_RosW *ros; 1248e27a552bSJed Brown PetscErrorCode ierr; 1249e27a552bSJed Brown 1250e27a552bSJed Brown PetscFunctionBegin; 1251e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1252e27a552bSJed Brown ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1253e27a552bSJed Brown #endif 1254e27a552bSJed Brown 1255e27a552bSJed Brown ts->ops->reset = TSReset_RosW; 1256e27a552bSJed Brown ts->ops->destroy = TSDestroy_RosW; 1257e27a552bSJed Brown ts->ops->view = TSView_RosW; 1258e27a552bSJed Brown ts->ops->setup = TSSetUp_RosW; 1259e27a552bSJed Brown ts->ops->step = TSStep_RosW; 1260e27a552bSJed Brown ts->ops->interpolate = TSInterpolate_RosW; 12611c3436cfSJed Brown ts->ops->evaluatestep = TSEvaluateStep_RosW; 1262e27a552bSJed Brown ts->ops->setfromoptions = TSSetFromOptions_RosW; 1263e27a552bSJed Brown ts->ops->snesfunction = SNESTSFormFunction_RosW; 1264e27a552bSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1265e27a552bSJed Brown 126661692a83SJed Brown ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr); 126761692a83SJed Brown ts->data = (void*)ros; 1268e27a552bSJed Brown 1269e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr); 1270e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr); 127161692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr); 1272e27a552bSJed Brown PetscFunctionReturn(0); 1273e27a552bSJed Brown } 1274e27a552bSJed Brown EXTERN_C_END 1275