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 */ 403ca35412SEmil 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 */ 563ca35412SEmil 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 184998eb97aSJed Brown TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four 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 199998eb97aSJed Brown TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four 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, 2381f80e275SEmil Constantinescu b = 1, 2391f80e275SEmil Constantinescu binterpt=1; 2401f80e275SEmil Constantinescu 2411f80e275SEmil Constantinescu ierr = TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,PETSC_NULL,1,&binterpt);CHKERRQ(ierr); 2423606a31eSEmil Constantinescu } 2433606a31eSEmil Constantinescu 2443606a31eSEmil Constantinescu { 2453606a31eSEmil Constantinescu const PetscReal 2463606a31eSEmil Constantinescu A= 0, 2473606a31eSEmil Constantinescu Gamma = 0.5, 2481f80e275SEmil Constantinescu b = 1, 2491f80e275SEmil Constantinescu binterpt=1; 2501f80e275SEmil Constantinescu ierr = TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,PETSC_NULL,1,&binterpt);CHKERRQ(ierr); 2513606a31eSEmil Constantinescu } 2523606a31eSEmil Constantinescu 2533606a31eSEmil Constantinescu { 25461692a83SJed Brown const PetscReal g = 1. + 1./PetscSqrtReal(2.0); 255e27a552bSJed Brown const PetscReal 25661692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 25761692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2581c3436cfSJed Brown b[2] = {0.5,0.5}, 2591c3436cfSJed Brown b1[2] = {1.0,0.0}; 2601f80e275SEmil Constantinescu PetscReal binterpt[2][2]; 2611f80e275SEmil Constantinescu binterpt[0][0]=g-1.0; 2621f80e275SEmil Constantinescu binterpt[1][0]=2.0-g; 2631f80e275SEmil Constantinescu binterpt[0][1]=g-1.5; 2641f80e275SEmil Constantinescu binterpt[1][1]=1.5-g; 2651f80e275SEmil Constantinescu ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]);CHKERRQ(ierr); 266e27a552bSJed Brown } 267e27a552bSJed Brown { 26861692a83SJed Brown const PetscReal g = 1. - 1./PetscSqrtReal(2.0); 269e27a552bSJed Brown const PetscReal 27061692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 27161692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2721c3436cfSJed Brown b[2] = {0.5,0.5}, 2731c3436cfSJed Brown b1[2] = {1.0,0.0}; 2741f80e275SEmil Constantinescu PetscReal binterpt[2][2]; 2751f80e275SEmil Constantinescu binterpt[0][0]=g-1.0; 2761f80e275SEmil Constantinescu binterpt[1][0]=2.0-g; 2771f80e275SEmil Constantinescu binterpt[0][1]=g-1.5; 2781f80e275SEmil Constantinescu binterpt[1][1]=1.5-g; 2791f80e275SEmil Constantinescu ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1,2,&binterpt[0][0]);CHKERRQ(ierr); 280fe7e6d57SJed Brown } 281fe7e6d57SJed Brown { 282fe7e6d57SJed Brown const PetscReal g = 7.8867513459481287e-01; 2831f80e275SEmil Constantinescu PetscReal binterpt[3][2]; 284fe7e6d57SJed Brown const PetscReal 285fe7e6d57SJed Brown A[3][3] = {{0,0,0}, 286fe7e6d57SJed Brown {1.5773502691896257e+00,0,0}, 287fe7e6d57SJed Brown {0.5,0,0}}, 288fe7e6d57SJed Brown Gamma[3][3] = {{g,0,0}, 289fe7e6d57SJed Brown {-1.5773502691896257e+00,g,0}, 29025833a80SEmil Constantinescu {-6.7075317547305480e-01,-1.7075317547305482e-01,g}}, 291fe7e6d57SJed Brown b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01}, 292fe7e6d57SJed Brown b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01}; 2931f80e275SEmil Constantinescu 2941f80e275SEmil Constantinescu binterpt[0][0]=-0.8094010767585034; 2951f80e275SEmil Constantinescu binterpt[1][0]=-0.5; 2961f80e275SEmil Constantinescu binterpt[2][0]=2.3094010767585034; 2971f80e275SEmil Constantinescu binterpt[0][1]=0.9641016151377548; 2981f80e275SEmil Constantinescu binterpt[1][1]=0.5; 2991f80e275SEmil Constantinescu binterpt[2][1]=-1.4641016151377548; 3001f80e275SEmil Constantinescu ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);CHKERRQ(ierr); 301fe7e6d57SJed Brown } 302fe7e6d57SJed Brown { 3033ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 304fe7e6d57SJed Brown const PetscReal g = 4.3586652150845900e-01; 305fe7e6d57SJed Brown const PetscReal 306fe7e6d57SJed Brown A[4][4] = {{0,0,0,0}, 307fe7e6d57SJed Brown {8.7173304301691801e-01,0,0,0}, 308fe7e6d57SJed Brown {8.4457060015369423e-01,-1.1299064236484185e-01,0,0}, 309fe7e6d57SJed Brown {0,0,1.,0}}, 310fe7e6d57SJed Brown Gamma[4][4] = {{g,0,0,0}, 311fe7e6d57SJed Brown {-8.7173304301691801e-01,g,0,0}, 312fe7e6d57SJed Brown {-9.0338057013044082e-01,5.4180672388095326e-02,g,0}, 313fe7e6d57SJed Brown {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}}, 314fe7e6d57SJed Brown b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01}, 3153ca35412SEmil Constantinescu b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01}; 3163ca35412SEmil Constantinescu 3173ca35412SEmil Constantinescu binterpt[0][0]=1.0564298455794094; 3183ca35412SEmil Constantinescu binterpt[1][0]=2.296429974281067; 3193ca35412SEmil Constantinescu binterpt[2][0]=-1.307599564525376; 3203ca35412SEmil Constantinescu binterpt[3][0]=-1.045260255335102; 3213ca35412SEmil Constantinescu binterpt[0][1]=-1.3864882699759573; 3223ca35412SEmil Constantinescu binterpt[1][1]=-8.262611700275677; 3233ca35412SEmil Constantinescu binterpt[2][1]=7.250979895056055; 3243ca35412SEmil Constantinescu binterpt[3][1]=2.398120075195581; 3253ca35412SEmil Constantinescu binterpt[0][2]=0.5721822314575016; 3263ca35412SEmil Constantinescu binterpt[1][2]=4.742931142090097; 3273ca35412SEmil Constantinescu binterpt[2][2]=-4.398120075195578; 3283ca35412SEmil Constantinescu binterpt[3][2]=-0.9169932983520199; 3293ca35412SEmil Constantinescu 3303ca35412SEmil Constantinescu ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr); 331e27a552bSJed Brown } 332ef3c5b88SJed Brown { 333ef3c5b88SJed Brown const PetscReal g = 0.5; 334ef3c5b88SJed Brown const PetscReal 335ef3c5b88SJed Brown A[4][4] = {{0,0,0,0}, 336ef3c5b88SJed Brown {0,0,0,0}, 337ef3c5b88SJed Brown {1.,0,0,0}, 338ef3c5b88SJed Brown {0.75,-0.25,0.5,0}}, 339ef3c5b88SJed Brown Gamma[4][4] = {{g,0,0,0}, 340ef3c5b88SJed Brown {1.,g,0,0}, 341ef3c5b88SJed Brown {-0.25,-0.25,g,0}, 342ef3c5b88SJed Brown {1./12,1./12,-2./3,g}}, 343ef3c5b88SJed Brown b[4] = {5./6,-1./6,-1./6,0.5}, 344ef3c5b88SJed Brown b2[4] = {0.75,-0.25,0.5,0}; 345f4aed992SEmil Constantinescu ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr); 346ef3c5b88SJed Brown } 347ef3c5b88SJed Brown { 348ef3c5b88SJed Brown const PetscReal g = 0.43586652150845899941601945119356; 349ef3c5b88SJed Brown const PetscReal 350ef3c5b88SJed Brown A[3][3] = {{0,0,0}, 351ef3c5b88SJed Brown {g,0,0}, 352ef3c5b88SJed Brown {g,0,0}}, 353ef3c5b88SJed Brown Gamma[3][3] = {{g,0,0}, 354ef3c5b88SJed Brown {-0.19294655696029095575009695436041,g,0}, 355ef3c5b88SJed Brown {0,1.74927148125794685173529749738960,g}}, 356ef3c5b88SJed Brown b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829}, 357ef3c5b88SJed Brown b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619}; 3581f80e275SEmil Constantinescu 3591f80e275SEmil Constantinescu PetscReal binterpt[3][2]; 3601f80e275SEmil Constantinescu binterpt[0][0]=3.793692883777660870425141387941; 3611f80e275SEmil Constantinescu binterpt[1][0]=-2.918692883777660870425141387941; 3621f80e275SEmil Constantinescu binterpt[2][0]=0.125; 3631f80e275SEmil Constantinescu binterpt[0][1]=-0.725741064379812106687651020584; 3641f80e275SEmil Constantinescu binterpt[1][1]=0.559074397713145440020984353917; 3651f80e275SEmil Constantinescu binterpt[2][1]=0.16666666666666666666666666666667; 3661f80e275SEmil Constantinescu 3671f80e275SEmil Constantinescu ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);CHKERRQ(ierr); 368ef3c5b88SJed Brown } 369b1c69cc3SEmil Constantinescu { 370aaf9cf16SJed Brown const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0; 371b1c69cc3SEmil Constantinescu const PetscReal 372b1c69cc3SEmil Constantinescu A[3][3] = {{0,0,0}, 373b1c69cc3SEmil Constantinescu {1,0,0}, 374b1c69cc3SEmil Constantinescu {0.25,0.25,0}}, 375b1c69cc3SEmil Constantinescu Gamma[3][3] = {{0,0,0}, 376aaf9cf16SJed Brown {(-3.0-s3)/6.0,g,0}, 377aaf9cf16SJed Brown {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}}, 378b1c69cc3SEmil Constantinescu b[3] = {1./6.,1./6.,2./3.}, 379b1c69cc3SEmil Constantinescu b2[3] = {1./4.,1./4.,1./2.}; 380c0cb691aSEmil Constantinescu 381c0cb691aSEmil Constantinescu PetscReal binterpt[3][2]; 382c0cb691aSEmil Constantinescu binterpt[0][0]=0.089316397477040902157517886164709; 383c0cb691aSEmil Constantinescu binterpt[1][0]=-0.91068360252295909784248211383529; 384c0cb691aSEmil Constantinescu binterpt[2][0]=1.8213672050459181956849642276706; 385c0cb691aSEmil Constantinescu binterpt[0][1]=0.077350269189625764509148780501957; 386c0cb691aSEmil Constantinescu binterpt[1][1]=1.077350269189625764509148780502; 387c0cb691aSEmil Constantinescu binterpt[2][1]=-1.1547005383792515290182975610039; 388c0cb691aSEmil Constantinescu ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2,2,&binterpt[0][0]);CHKERRQ(ierr); 389b1c69cc3SEmil Constantinescu } 390b1c69cc3SEmil Constantinescu 391b1c69cc3SEmil Constantinescu { 392b1c69cc3SEmil Constantinescu const PetscReal 393b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 394b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 395b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 396b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 397b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 398b1c69cc3SEmil Constantinescu {0.0,1./4.,0,0}, 399b1c69cc3SEmil Constantinescu {-2.,-2./3.,2./3.,0}, 400b1c69cc3SEmil Constantinescu {1./2.,5./36.,-2./9,0}}, 401b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 402b1c69cc3SEmil Constantinescu b2[4] = {1./8.,3./4.,1./8.,0}; 403c0cb691aSEmil Constantinescu PetscReal binterpt[4][3]; 404c0cb691aSEmil Constantinescu binterpt[0][0]=6.25; 405c0cb691aSEmil Constantinescu binterpt[1][0]=-30.25; 406c0cb691aSEmil Constantinescu binterpt[2][0]=1.75; 407c0cb691aSEmil Constantinescu binterpt[3][0]=23.25; 408c0cb691aSEmil Constantinescu binterpt[0][1]=-9.75; 409c0cb691aSEmil Constantinescu binterpt[1][1]=58.75; 410c0cb691aSEmil Constantinescu binterpt[2][1]=-3.25; 411c0cb691aSEmil Constantinescu binterpt[3][1]=-45.75; 412c0cb691aSEmil Constantinescu binterpt[0][2]=3.6666666666666666666666666666667; 413c0cb691aSEmil Constantinescu binterpt[1][2]=-28.333333333333333333333333333333; 414c0cb691aSEmil Constantinescu binterpt[2][2]=1.6666666666666666666666666666667; 415c0cb691aSEmil Constantinescu binterpt[3][2]=23.; 416c0cb691aSEmil Constantinescu ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr); 417b1c69cc3SEmil Constantinescu } 418b1c69cc3SEmil Constantinescu 419b1c69cc3SEmil Constantinescu { 420b1c69cc3SEmil Constantinescu const PetscReal 421b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 422b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 423b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 424b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 425b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 426b1c69cc3SEmil Constantinescu {0.0,3./4.,0,0}, 427b1c69cc3SEmil Constantinescu {-2./3.,-23./9.,2./9.,0}, 428b1c69cc3SEmil Constantinescu {1./18.,65./108.,-2./27,0}}, 429b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 430b1c69cc3SEmil Constantinescu b2[4] = {3./16.,10./16.,3./16.,0}; 431c0cb691aSEmil Constantinescu 432c0cb691aSEmil Constantinescu PetscReal binterpt[4][3]; 433c0cb691aSEmil Constantinescu binterpt[0][0]=1.6911764705882352941176470588235; 434c0cb691aSEmil Constantinescu binterpt[1][0]=3.6813725490196078431372549019608; 435c0cb691aSEmil Constantinescu binterpt[2][0]=0.23039215686274509803921568627451; 436c0cb691aSEmil Constantinescu binterpt[3][0]=-4.6029411764705882352941176470588; 437c0cb691aSEmil Constantinescu binterpt[0][1]=-0.95588235294117647058823529411765; 438c0cb691aSEmil Constantinescu binterpt[1][1]=-6.2401960784313725490196078431373; 439c0cb691aSEmil Constantinescu binterpt[2][1]=-0.31862745098039215686274509803922; 440c0cb691aSEmil Constantinescu binterpt[3][1]=7.5147058823529411764705882352941; 441c0cb691aSEmil Constantinescu binterpt[0][2]=-0.56862745098039215686274509803922; 442c0cb691aSEmil Constantinescu binterpt[1][2]=2.7254901960784313725490196078431; 443c0cb691aSEmil Constantinescu binterpt[2][2]=0.25490196078431372549019607843137; 444c0cb691aSEmil Constantinescu binterpt[3][2]=-2.4117647058823529411764705882353; 445c0cb691aSEmil Constantinescu ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr); 446b1c69cc3SEmil Constantinescu } 447753f8adbSEmil Constantinescu 448753f8adbSEmil Constantinescu { 449753f8adbSEmil Constantinescu PetscReal A[4][4],Gamma[4][4],b[4],b2[4]; 4503ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 451753f8adbSEmil Constantinescu 452753f8adbSEmil Constantinescu Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816; 45305e8e825SJed Brown Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0; 454753f8adbSEmil Constantinescu Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476; 455753f8adbSEmil Constantinescu Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816; 45605e8e825SJed Brown Gamma[1][2]=0; Gamma[1][3]=0; 457753f8adbSEmil Constantinescu Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903; 458753f8adbSEmil Constantinescu Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131; 459753f8adbSEmil Constantinescu Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816; 46005e8e825SJed Brown Gamma[2][3]=0; 461753f8adbSEmil Constantinescu Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783; 462753f8adbSEmil Constantinescu Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984; 463753f8adbSEmil Constantinescu Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198; 464753f8adbSEmil Constantinescu Gamma[3][3]=0; 465753f8adbSEmil Constantinescu 46605e8e825SJed Brown A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0; 467753f8adbSEmil Constantinescu A[1][0]=0.8717330430169179988320388950590125027645343373957631; 46805e8e825SJed Brown A[1][1]=0; A[1][2]=0; A[1][3]=0; 469753f8adbSEmil Constantinescu A[2][0]=0.5275890119763004115618079766722914408876108660811028; 470753f8adbSEmil Constantinescu A[2][1]=0.07241098802369958843819203208518599088698057726988732; 47105e8e825SJed Brown A[2][2]=0; A[2][3]=0; 472753f8adbSEmil Constantinescu A[3][0]=0.3990960076760701320627260685975778145384666450351314; 473753f8adbSEmil Constantinescu A[3][1]=-0.4375576546135194437228463747348862825846903771419953; 474753f8adbSEmil Constantinescu A[3][2]=1.038461646937449311660120300601880176655352737312713; 47505e8e825SJed Brown A[3][3]=0; 476753f8adbSEmil Constantinescu 477753f8adbSEmil Constantinescu b[0]=0.1876410243467238251612921333138006734899663569186926; 478753f8adbSEmil Constantinescu b[1]=-0.5952974735769549480478230473706443582188442040780541; 479753f8adbSEmil Constantinescu b[2]=0.9717899277217721234705114616271378792182450260943198; 480753f8adbSEmil Constantinescu b[3]=0.4358665215084589994160194475295062513822671686978816; 481753f8adbSEmil Constantinescu 482753f8adbSEmil Constantinescu b2[0]=0.2147402862233891404862383521089097657790734483804460; 483753f8adbSEmil Constantinescu b2[1]=-0.4851622638849390928209050538171743017757490232519684; 484753f8adbSEmil Constantinescu b2[2]=0.8687250025203875511662123688667549217531982787600080; 485753f8adbSEmil Constantinescu b2[3]=0.4016969751411624011684543450940068201770721128357014; 486753f8adbSEmil Constantinescu 4873ca35412SEmil Constantinescu binterpt[0][0]=2.2565812720167954547104627844105; 4883ca35412SEmil Constantinescu binterpt[1][0]=1.349166413351089573796243820819; 4893ca35412SEmil Constantinescu binterpt[2][0]=-2.4695174540533503758652847586647; 4903ca35412SEmil Constantinescu binterpt[3][0]=-0.13623023131453465264142184656474; 4913ca35412SEmil Constantinescu binterpt[0][1]=-3.0826699111559187902922463354557; 4923ca35412SEmil Constantinescu binterpt[1][1]=-2.4689115685996042534544925650515; 4933ca35412SEmil Constantinescu binterpt[2][1]=5.7428279814696677152129332773553; 4943ca35412SEmil Constantinescu binterpt[3][1]=-0.19124650171414467146619437684812; 4953ca35412SEmil Constantinescu binterpt[0][2]=1.0137296634858471607430756831148; 4963ca35412SEmil Constantinescu binterpt[1][2]=0.52444768167155973161042570784064; 4973ca35412SEmil Constantinescu binterpt[2][2]=-2.3015205996945452158771370439586; 4983ca35412SEmil Constantinescu binterpt[3][2]=0.76334325453713832352363565300308; 499f4aed992SEmil Constantinescu 500f73f8d2cSSatish Balay ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2,3,&binterpt[0][0]);CHKERRQ(ierr); 501753f8adbSEmil Constantinescu } 502753f8adbSEmil Constantinescu 503e27a552bSJed Brown PetscFunctionReturn(0); 504e27a552bSJed Brown } 505e27a552bSJed Brown 506e27a552bSJed Brown #undef __FUNCT__ 507e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy" 508e27a552bSJed Brown /*@C 509e27a552bSJed Brown TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 510e27a552bSJed Brown 511e27a552bSJed Brown Not Collective 512e27a552bSJed Brown 513e27a552bSJed Brown Level: advanced 514e27a552bSJed Brown 515e27a552bSJed Brown .keywords: TSRosW, register, destroy 516e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 517e27a552bSJed Brown @*/ 518e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void) 519e27a552bSJed Brown { 520e27a552bSJed Brown PetscErrorCode ierr; 52161692a83SJed Brown RosWTableauLink link; 522e27a552bSJed Brown 523e27a552bSJed Brown PetscFunctionBegin; 52461692a83SJed Brown while ((link = RosWTableauList)) { 52561692a83SJed Brown RosWTableau t = &link->tab; 52661692a83SJed Brown RosWTableauList = link->next; 52761692a83SJed Brown ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr); 52843b21953SEmil Constantinescu ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr); 529fe7e6d57SJed Brown ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr); 530f4aed992SEmil Constantinescu ierr = PetscFree(t->binterpt);CHKERRQ(ierr); 531e27a552bSJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 532e27a552bSJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 533e27a552bSJed Brown } 534e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_FALSE; 535e27a552bSJed Brown PetscFunctionReturn(0); 536e27a552bSJed Brown } 537e27a552bSJed Brown 538e27a552bSJed Brown #undef __FUNCT__ 539e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage" 540e27a552bSJed Brown /*@C 541e27a552bSJed Brown TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called 542e27a552bSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW() 543e27a552bSJed Brown when using static libraries. 544e27a552bSJed Brown 545e27a552bSJed Brown Input Parameter: 546e27a552bSJed Brown path - The dynamic library path, or PETSC_NULL 547e27a552bSJed Brown 548e27a552bSJed Brown Level: developer 549e27a552bSJed Brown 550e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package 551e27a552bSJed Brown .seealso: PetscInitialize() 552e27a552bSJed Brown @*/ 553e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[]) 554e27a552bSJed Brown { 555e27a552bSJed Brown PetscErrorCode ierr; 556e27a552bSJed Brown 557e27a552bSJed Brown PetscFunctionBegin; 558e27a552bSJed Brown if (TSRosWPackageInitialized) PetscFunctionReturn(0); 559e27a552bSJed Brown TSRosWPackageInitialized = PETSC_TRUE; 560e27a552bSJed Brown ierr = TSRosWRegisterAll();CHKERRQ(ierr); 561e27a552bSJed Brown ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr); 562e27a552bSJed Brown PetscFunctionReturn(0); 563e27a552bSJed Brown } 564e27a552bSJed Brown 565e27a552bSJed Brown #undef __FUNCT__ 566e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage" 567e27a552bSJed Brown /*@C 568e27a552bSJed Brown TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is 569e27a552bSJed Brown called from PetscFinalize(). 570e27a552bSJed Brown 571e27a552bSJed Brown Level: developer 572e27a552bSJed Brown 573e27a552bSJed Brown .keywords: Petsc, destroy, package 574e27a552bSJed Brown .seealso: PetscFinalize() 575e27a552bSJed Brown @*/ 576e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void) 577e27a552bSJed Brown { 578e27a552bSJed Brown PetscErrorCode ierr; 579e27a552bSJed Brown 580e27a552bSJed Brown PetscFunctionBegin; 581e27a552bSJed Brown TSRosWPackageInitialized = PETSC_FALSE; 582e27a552bSJed Brown ierr = TSRosWRegisterDestroy();CHKERRQ(ierr); 583e27a552bSJed Brown PetscFunctionReturn(0); 584e27a552bSJed Brown } 585e27a552bSJed Brown 586e27a552bSJed Brown #undef __FUNCT__ 587e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister" 588e27a552bSJed Brown /*@C 58961692a83SJed Brown TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 590e27a552bSJed Brown 591e27a552bSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 592e27a552bSJed Brown 593e27a552bSJed Brown Input Parameters: 594e27a552bSJed Brown + name - identifier for method 595e27a552bSJed Brown . order - approximation order of method 596e27a552bSJed Brown . s - number of stages, this is the dimension of the matrices below 59761692a83SJed Brown . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 59861692a83SJed Brown . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 599fe7e6d57SJed Brown . b - Step completion table (dimension s) 600fe7e6d57SJed Brown - bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available) 601f4aed992SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt 602f4aed992SEmil Constantinescu . binterpt - Coefficients of the interpolation formula (dimension s*pinterp) 603e27a552bSJed Brown 604e27a552bSJed Brown Notes: 60561692a83SJed Brown Several Rosenbrock W methods are provided, this function is only needed to create new methods. 606e27a552bSJed Brown 607e27a552bSJed Brown Level: advanced 608e27a552bSJed Brown 609e27a552bSJed Brown .keywords: TS, register 610e27a552bSJed Brown 611e27a552bSJed Brown .seealso: TSRosW 612e27a552bSJed Brown @*/ 613e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s, 614f4aed992SEmil Constantinescu const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[], 615f4aed992SEmil Constantinescu PetscInt pinterp,const PetscReal binterpt[]) 616e27a552bSJed Brown { 617e27a552bSJed Brown PetscErrorCode ierr; 61861692a83SJed Brown RosWTableauLink link; 61961692a83SJed Brown RosWTableau t; 62061692a83SJed Brown PetscInt i,j,k; 62161692a83SJed Brown PetscScalar *GammaInv; 622e27a552bSJed Brown 623e27a552bSJed Brown PetscFunctionBegin; 624fe7e6d57SJed Brown PetscValidCharPointer(name,1); 625fe7e6d57SJed Brown PetscValidPointer(A,4); 626fe7e6d57SJed Brown PetscValidPointer(Gamma,5); 627fe7e6d57SJed Brown PetscValidPointer(b,6); 628fe7e6d57SJed Brown if (bembed) PetscValidPointer(bembed,7); 629fe7e6d57SJed Brown 630e27a552bSJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 631e27a552bSJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 632e27a552bSJed Brown t = &link->tab; 633e27a552bSJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 634e27a552bSJed Brown t->order = order; 635e27a552bSJed Brown t->s = s; 63661692a83SJed 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); 63743b21953SEmil 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); 638e27a552bSJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 63961692a83SJed Brown ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 64043b21953SEmil Constantinescu ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 64161692a83SJed Brown ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 642fe7e6d57SJed Brown if (bembed) { 643fe7e6d57SJed Brown ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr); 644fe7e6d57SJed Brown ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 645fe7e6d57SJed Brown } 64661692a83SJed Brown for (i=0; i<s; i++) { 64761692a83SJed Brown t->ASum[i] = 0; 64861692a83SJed Brown t->GammaSum[i] = 0; 64961692a83SJed Brown for (j=0; j<s; j++) { 65061692a83SJed Brown t->ASum[i] += A[i*s+j]; 651fe7e6d57SJed Brown t->GammaSum[i] += Gamma[i*s+j]; 65261692a83SJed Brown } 65361692a83SJed Brown } 65461692a83SJed Brown ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */ 65561692a83SJed Brown for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i]; 656fd96d5b0SEmil Constantinescu for (i=0; i<s; i++) { 657fd96d5b0SEmil Constantinescu if (Gamma[i*s+i] == 0.0) { 658fd96d5b0SEmil Constantinescu GammaInv[i*s+i] = 1.0; 659c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_TRUE; 660fd96d5b0SEmil Constantinescu } else { 661c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_FALSE; 662fd96d5b0SEmil Constantinescu } 663fd96d5b0SEmil Constantinescu } 664fd96d5b0SEmil Constantinescu 66561692a83SJed Brown switch (s) { 66661692a83SJed Brown case 1: GammaInv[0] = 1./GammaInv[0]; break; 66761692a83SJed Brown case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break; 66861692a83SJed Brown case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break; 66961692a83SJed Brown case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break; 67061692a83SJed Brown case 5: { 67161692a83SJed Brown PetscInt ipvt5[5]; 67261692a83SJed Brown MatScalar work5[5*5]; 67361692a83SJed Brown ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break; 67461692a83SJed Brown } 67561692a83SJed Brown case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break; 67661692a83SJed Brown case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break; 67761692a83SJed Brown default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s); 67861692a83SJed Brown } 67961692a83SJed Brown for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 68061692a83SJed Brown ierr = PetscFree(GammaInv);CHKERRQ(ierr); 68143b21953SEmil Constantinescu 68243b21953SEmil Constantinescu for (i=0; i<s; i++) { 68343b21953SEmil Constantinescu for (k=0; k<i+1; k++) { 68443b21953SEmil Constantinescu t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]); 68543b21953SEmil Constantinescu for (j=k+1; j<i+1; j++) { 68643b21953SEmil Constantinescu t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]); 68743b21953SEmil Constantinescu } 68843b21953SEmil Constantinescu } 68943b21953SEmil Constantinescu } 69043b21953SEmil Constantinescu 69161692a83SJed Brown for (i=0; i<s; i++) { 69261692a83SJed Brown for (j=0; j<s; j++) { 69361692a83SJed Brown t->At[i*s+j] = 0; 69461692a83SJed Brown for (k=0; k<s; k++) { 69561692a83SJed Brown t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j]; 69661692a83SJed Brown } 69761692a83SJed Brown } 69861692a83SJed Brown t->bt[i] = 0; 69961692a83SJed Brown for (j=0; j<s; j++) { 70061692a83SJed Brown t->bt[i] += t->b[j] * t->GammaInv[j*s+i]; 70161692a83SJed Brown } 702fe7e6d57SJed Brown if (bembed) { 703fe7e6d57SJed Brown t->bembedt[i] = 0; 704fe7e6d57SJed Brown for (j=0; j<s; j++) { 705fe7e6d57SJed Brown t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i]; 706fe7e6d57SJed Brown } 707fe7e6d57SJed Brown } 70861692a83SJed Brown } 7098d59e960SJed Brown t->ccfl = 1.0; /* Fix this */ 7108d59e960SJed Brown 711f4aed992SEmil Constantinescu t->pinterp = pinterp; 7123ca35412SEmil Constantinescu ierr = PetscMalloc(s*pinterp*sizeof(binterpt[0]),&t->binterpt);CHKERRQ(ierr); 7133ca35412SEmil Constantinescu ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 71461692a83SJed Brown link->next = RosWTableauList; 71561692a83SJed Brown RosWTableauList = link; 716e27a552bSJed Brown PetscFunctionReturn(0); 717e27a552bSJed Brown } 718e27a552bSJed Brown 719e27a552bSJed Brown #undef __FUNCT__ 7201c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW" 7211c3436cfSJed Brown /* 7221c3436cfSJed Brown The step completion formula is 7231c3436cfSJed Brown 7241c3436cfSJed Brown x1 = x0 + b^T Y 7251c3436cfSJed Brown 7261c3436cfSJed Brown where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 7271c3436cfSJed Brown updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 7281c3436cfSJed Brown 7291c3436cfSJed Brown x1e = x0 + be^T Y 7301c3436cfSJed Brown = x1 - b^T Y + be^T Y 7311c3436cfSJed Brown = x1 + (be - b)^T Y 7321c3436cfSJed Brown 7331c3436cfSJed Brown so we can evaluate the method of different order even after the step has been optimistically completed. 7341c3436cfSJed Brown */ 7351c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done) 7361c3436cfSJed Brown { 7371c3436cfSJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 7381c3436cfSJed Brown RosWTableau tab = ros->tableau; 7391c3436cfSJed Brown PetscScalar *w = ros->work; 7401c3436cfSJed Brown PetscInt i; 7411c3436cfSJed Brown PetscErrorCode ierr; 7421c3436cfSJed Brown 7431c3436cfSJed Brown PetscFunctionBegin; 7441c3436cfSJed Brown if (order == tab->order) { 745108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */ 7461c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 747de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bt[i]; 748de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 749108c343cSJed Brown } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 7501c3436cfSJed Brown if (done) *done = PETSC_TRUE; 7511c3436cfSJed Brown PetscFunctionReturn(0); 7521c3436cfSJed Brown } else if (order == tab->order-1) { 7531c3436cfSJed Brown if (!tab->bembedt) goto unavailable; 754108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */ 7551c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 756de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i]; 757de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 758108c343cSJed Brown } else { /* Use rollback-and-recomplete formula (bembedt - bt) */ 759108c343cSJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 760108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 761108c343cSJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 7621c3436cfSJed Brown } 7631c3436cfSJed Brown if (done) *done = PETSC_TRUE; 7641c3436cfSJed Brown PetscFunctionReturn(0); 7651c3436cfSJed Brown } 7661c3436cfSJed Brown unavailable: 7671c3436cfSJed Brown if (done) *done = PETSC_FALSE; 7681c3436cfSJed 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); 7691c3436cfSJed Brown PetscFunctionReturn(0); 7701c3436cfSJed Brown } 7711c3436cfSJed Brown 7721c3436cfSJed Brown #undef __FUNCT__ 773e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW" 774e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts) 775e27a552bSJed Brown { 77661692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 77761692a83SJed Brown RosWTableau tab = ros->tableau; 778e27a552bSJed Brown const PetscInt s = tab->s; 7791c3436cfSJed Brown const PetscReal *At = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv; 7800feba352SEmil Constantinescu const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr; 781c17803e7SJed Brown const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 78261692a83SJed Brown PetscScalar *w = ros->work; 7837d4bf2deSEmil Constantinescu Vec *Y = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage; 784e27a552bSJed Brown SNES snes; 7851c3436cfSJed Brown TSAdapt adapt; 7861c3436cfSJed Brown PetscInt i,j,its,lits,reject,next_scheme; 787cdbf8f93SLisandro Dalcin PetscReal next_time_step; 7881c3436cfSJed Brown PetscBool accept; 789e27a552bSJed Brown PetscErrorCode ierr; 7900feba352SEmil Constantinescu MatStructure str; 791e27a552bSJed Brown 792e27a552bSJed Brown PetscFunctionBegin; 793e27a552bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 794cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 7951c3436cfSJed Brown accept = PETSC_TRUE; 796108c343cSJed Brown ros->status = TS_STEP_INCOMPLETE; 797e27a552bSJed Brown 79897335746SJed Brown for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 7991c3436cfSJed Brown const PetscReal h = ts->time_step; 8003ca35412SEmil Constantinescu ierr = VecCopy(ts->vec_sol,ros->VecSolPrev);CHKERRQ(ierr);/*move this at the end*/ 801e27a552bSJed Brown for (i=0; i<s; i++) { 8021c3436cfSJed Brown ros->stage_time = ts->ptime + h*ASum[i]; 803c17803e7SJed Brown if (GammaZeroDiag[i]) { 804c17803e7SJed Brown ros->stage_explicit = PETSC_TRUE; 805fd96d5b0SEmil Constantinescu ros->shift = 1./h; 806c17803e7SJed Brown } else { 807c17803e7SJed Brown ros->stage_explicit = PETSC_FALSE; 80861692a83SJed Brown ros->shift = 1./(h*Gamma[i*s+i]); 809fd96d5b0SEmil Constantinescu } 81061692a83SJed Brown 81161692a83SJed Brown ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr); 812de19f811SJed Brown for (j=0; j<i; j++) w[j] = At[i*s+j]; 813de19f811SJed Brown ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 81461692a83SJed Brown 81561692a83SJed Brown for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j]; 81661692a83SJed Brown ierr = VecZeroEntries(Zdot);CHKERRQ(ierr); 81761692a83SJed Brown ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr); 81861692a83SJed Brown 819e27a552bSJed Brown /* Initial guess taken from last stage */ 82061692a83SJed Brown ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr); 82161692a83SJed Brown 8227d4bf2deSEmil Constantinescu if (!ros->stage_explicit) { 82361692a83SJed Brown if (!ros->recompute_jacobian && !i) { 82461692a83SJed Brown ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 82561692a83SJed Brown } 82661692a83SJed Brown ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr); 827e27a552bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 828e27a552bSJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 829e27a552bSJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 83097335746SJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 83197335746SJed Brown ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 83297335746SJed Brown if (!accept) goto reject_step; 8337d4bf2deSEmil Constantinescu } else { 8341ce71dffSSatish Balay Mat J,Jp; 8350feba352SEmil Constantinescu ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */ 8360feba352SEmil Constantinescu ierr = TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);CHKERRQ(ierr); 8370feba352SEmil Constantinescu ierr = VecScale(Y[i],-1.0); 8380feba352SEmil Constantinescu ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/ 8390feba352SEmil Constantinescu 8400feba352SEmil Constantinescu ierr = VecZeroEntries(Zstage);CHKERRQ(ierr); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */ 8410feba352SEmil Constantinescu for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j]; 8420feba352SEmil Constantinescu ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 8430feba352SEmil Constantinescu /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */ 8440feba352SEmil Constantinescu str = SAME_NONZERO_PATTERN; 845*ccbc64bcSJed Brown ierr = TSGetIJacobian(ts,&J,&Jp,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 8460feba352SEmil Constantinescu ierr = TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);CHKERRQ(ierr); 8470feba352SEmil Constantinescu ierr = MatMult(J,Zstage,Zdot); 8480feba352SEmil Constantinescu 8490feba352SEmil Constantinescu ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); 8500feba352SEmil Constantinescu ierr = VecScale(Y[i],h); 8517d4bf2deSEmil Constantinescu ts->linear_its += 1; 8527d4bf2deSEmil Constantinescu } 853e27a552bSJed Brown } 8541c3436cfSJed Brown ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 855108c343cSJed Brown ros->status = TS_STEP_PENDING; 856e27a552bSJed Brown 8571c3436cfSJed Brown /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 8581c3436cfSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 8591c3436cfSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 8608d59e960SJed Brown ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 8611c3436cfSJed Brown ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 8621c3436cfSJed Brown if (accept) { 8631c3436cfSJed Brown /* ignore next_scheme for now */ 864e27a552bSJed Brown ts->ptime += ts->time_step; 865cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 866e27a552bSJed Brown ts->steps++; 867108c343cSJed Brown ros->status = TS_STEP_COMPLETE; 8681c3436cfSJed Brown break; 8691c3436cfSJed Brown } else { /* Roll back the current step */ 8701c3436cfSJed Brown for (i=0; i<s; i++) w[i] = -tab->bt[i]; 8711c3436cfSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr); 8721c3436cfSJed Brown ts->time_step = next_time_step; 873108c343cSJed Brown ros->status = TS_STEP_INCOMPLETE; 8741c3436cfSJed Brown } 875476b6736SJed Brown reject_step: continue; 8761c3436cfSJed Brown } 877b2ce242eSJed Brown if (ros->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 878e27a552bSJed Brown PetscFunctionReturn(0); 879e27a552bSJed Brown } 880e27a552bSJed Brown 881e27a552bSJed Brown #undef __FUNCT__ 882e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW" 883e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X) 884e27a552bSJed Brown { 88561692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 886f4aed992SEmil Constantinescu PetscInt s = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j; 887f4aed992SEmil Constantinescu PetscReal h; 888f4aed992SEmil Constantinescu PetscReal tt,t; 889f4aed992SEmil Constantinescu PetscScalar *bt; 890f4aed992SEmil Constantinescu const PetscReal *Bt = ros->tableau->binterpt; 891f4aed992SEmil Constantinescu PetscErrorCode ierr; 892f4aed992SEmil Constantinescu const PetscReal *GammaInv = ros->tableau->GammaInv; 893f4aed992SEmil Constantinescu PetscScalar *w = ros->work; 894f4aed992SEmil Constantinescu Vec *Y = ros->Y; 895e27a552bSJed Brown 896e27a552bSJed Brown PetscFunctionBegin; 897f4aed992SEmil Constantinescu if (!Bt) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name); 898f4aed992SEmil Constantinescu 899f4aed992SEmil Constantinescu switch (ros->status) { 900f4aed992SEmil Constantinescu case TS_STEP_INCOMPLETE: 901f4aed992SEmil Constantinescu case TS_STEP_PENDING: 902f4aed992SEmil Constantinescu h = ts->time_step; 903f4aed992SEmil Constantinescu t = (itime - ts->ptime)/h; 904f4aed992SEmil Constantinescu break; 905f4aed992SEmil Constantinescu case TS_STEP_COMPLETE: 906f4aed992SEmil Constantinescu h = ts->time_step_prev; 907f4aed992SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 908f4aed992SEmil Constantinescu break; 909f4aed992SEmil Constantinescu default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus"); 910f4aed992SEmil Constantinescu } 9113ca35412SEmil Constantinescu ierr = PetscMalloc(s*sizeof(bt[0]),&bt);CHKERRQ(ierr); 912f4aed992SEmil Constantinescu for (i=0; i<s; i++) bt[i] = 0; 913f4aed992SEmil Constantinescu for (j=0,tt=t; j<pinterp; j++,tt*=t) { 914f4aed992SEmil Constantinescu for (i=0; i<s; i++) { 9153ca35412SEmil Constantinescu bt[i] += Bt[i*pinterp+j] * tt; 916f4aed992SEmil Constantinescu } 917f4aed992SEmil Constantinescu } 918f4aed992SEmil Constantinescu 919f4aed992SEmil Constantinescu /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */ 920f4aed992SEmil Constantinescu /*X<-0*/ 921f4aed992SEmil Constantinescu ierr = VecZeroEntries(X);CHKERRQ(ierr); 922f4aed992SEmil Constantinescu 923f4aed992SEmil Constantinescu /*X<- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */ 9243ca35412SEmil Constantinescu for (j=0; j<s; j++) w[j]=0; 9253ca35412SEmil Constantinescu for (j=0; j<s; j++) { 9263ca35412SEmil Constantinescu for (i=j; i<s; i++) { 9273ca35412SEmil Constantinescu w[j] += bt[i]*GammaInv[i*s+j]; 928f4aed992SEmil Constantinescu } 9293ca35412SEmil Constantinescu } 9303ca35412SEmil Constantinescu ierr = VecMAXPY(X,i,w,Y);CHKERRQ(ierr); 931f4aed992SEmil Constantinescu 932f4aed992SEmil Constantinescu /*X<-y(t) + X*/ 9333ca35412SEmil Constantinescu ierr = VecAXPY(X,1.0,ros->VecSolPrev);CHKERRQ(ierr); 934f4aed992SEmil Constantinescu 935f4aed992SEmil Constantinescu ierr = PetscFree(bt);CHKERRQ(ierr); 936f4aed992SEmil Constantinescu 937e27a552bSJed Brown PetscFunctionReturn(0); 938e27a552bSJed Brown } 939e27a552bSJed Brown 940e27a552bSJed Brown /*------------------------------------------------------------*/ 941e27a552bSJed Brown #undef __FUNCT__ 942e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW" 943e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts) 944e27a552bSJed Brown { 94561692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 946e27a552bSJed Brown PetscInt s; 947e27a552bSJed Brown PetscErrorCode ierr; 948e27a552bSJed Brown 949e27a552bSJed Brown PetscFunctionBegin; 95061692a83SJed Brown if (!ros->tableau) PetscFunctionReturn(0); 95161692a83SJed Brown s = ros->tableau->s; 95261692a83SJed Brown ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr); 95361692a83SJed Brown ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr); 95461692a83SJed Brown ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr); 95561692a83SJed Brown ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr); 95661692a83SJed Brown ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr); 9573ca35412SEmil Constantinescu ierr = VecDestroy(&ros->VecSolPrev);CHKERRQ(ierr); 95861692a83SJed Brown ierr = PetscFree(ros->work);CHKERRQ(ierr); 959e27a552bSJed Brown PetscFunctionReturn(0); 960e27a552bSJed Brown } 961e27a552bSJed Brown 962e27a552bSJed Brown #undef __FUNCT__ 963e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW" 964e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts) 965e27a552bSJed Brown { 966e27a552bSJed Brown PetscErrorCode ierr; 967e27a552bSJed Brown 968e27a552bSJed Brown PetscFunctionBegin; 969e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 970e27a552bSJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 971e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr); 972e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr); 97361692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr); 974e27a552bSJed Brown PetscFunctionReturn(0); 975e27a552bSJed Brown } 976e27a552bSJed Brown 977e27a552bSJed Brown /* 978e27a552bSJed Brown This defines the nonlinear equation that is to be solved with SNES 979e27a552bSJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 980e27a552bSJed Brown */ 981e27a552bSJed Brown #undef __FUNCT__ 982e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW" 983e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts) 984e27a552bSJed Brown { 98561692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 986e27a552bSJed Brown PetscErrorCode ierr; 987e27a552bSJed Brown 988e27a552bSJed Brown PetscFunctionBegin; 98961692a83SJed Brown ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */ 99061692a83SJed Brown ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr); /* Ystage = X + Zstage */ 99161692a83SJed Brown ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr); 992e27a552bSJed Brown PetscFunctionReturn(0); 993e27a552bSJed Brown } 994e27a552bSJed Brown 995e27a552bSJed Brown #undef __FUNCT__ 996e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW" 997e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 998e27a552bSJed Brown { 99961692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1000e27a552bSJed Brown PetscErrorCode ierr; 1001e27a552bSJed Brown 1002e27a552bSJed Brown PetscFunctionBegin; 100361692a83SJed Brown /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 100461692a83SJed Brown ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 1005e27a552bSJed Brown PetscFunctionReturn(0); 1006e27a552bSJed Brown } 1007e27a552bSJed Brown 1008e27a552bSJed Brown #undef __FUNCT__ 1009e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW" 1010e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts) 1011e27a552bSJed Brown { 101261692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 101361692a83SJed Brown RosWTableau tab = ros->tableau; 1014e27a552bSJed Brown PetscInt s = tab->s; 1015e27a552bSJed Brown PetscErrorCode ierr; 1016e27a552bSJed Brown 1017e27a552bSJed Brown PetscFunctionBegin; 101861692a83SJed Brown if (!ros->tableau) { 1019e27a552bSJed Brown ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr); 1020e27a552bSJed Brown } 102161692a83SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr); 102261692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr); 102361692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr); 102461692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr); 102561692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr); 10263ca35412SEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&ros->VecSolPrev);CHKERRQ(ierr); 102761692a83SJed Brown ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr); 1028e27a552bSJed Brown PetscFunctionReturn(0); 1029e27a552bSJed Brown } 1030e27a552bSJed Brown /*------------------------------------------------------------*/ 1031e27a552bSJed Brown 1032e27a552bSJed Brown #undef __FUNCT__ 1033e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW" 1034e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts) 1035e27a552bSJed Brown { 103661692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1037e27a552bSJed Brown PetscErrorCode ierr; 103861692a83SJed Brown char rostype[256]; 1039e27a552bSJed Brown 1040e27a552bSJed Brown PetscFunctionBegin; 1041e27a552bSJed Brown ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr); 1042e27a552bSJed Brown { 104361692a83SJed Brown RosWTableauLink link; 1044e27a552bSJed Brown PetscInt count,choice; 1045e27a552bSJed Brown PetscBool flg; 1046e27a552bSJed Brown const char **namelist; 104761692a83SJed Brown SNES snes; 104861692a83SJed Brown 104961692a83SJed Brown ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr); 105061692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) ; 1051e27a552bSJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 105261692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 105361692a83SJed Brown ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr); 105461692a83SJed Brown ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr); 1055e27a552bSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 105661692a83SJed Brown 105761692a83SJed Brown ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr); 105861692a83SJed Brown 105961692a83SJed Brown /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 106061692a83SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 106161692a83SJed Brown if (!((PetscObject)snes)->type_name) { 106261692a83SJed Brown ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 106361692a83SJed Brown } 106461692a83SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1065e27a552bSJed Brown } 1066e27a552bSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1067e27a552bSJed Brown PetscFunctionReturn(0); 1068e27a552bSJed Brown } 1069e27a552bSJed Brown 1070e27a552bSJed Brown #undef __FUNCT__ 1071e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray" 1072e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 1073e27a552bSJed Brown { 1074e27a552bSJed Brown PetscErrorCode ierr; 1075e408995aSJed Brown PetscInt i; 1076e408995aSJed Brown size_t left,count; 1077e27a552bSJed Brown char *p; 1078e27a552bSJed Brown 1079e27a552bSJed Brown PetscFunctionBegin; 1080e408995aSJed Brown for (i=0,p=buf,left=len; i<n; i++) { 1081e408995aSJed Brown ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 1082e27a552bSJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 1083e27a552bSJed Brown left -= count; 1084e27a552bSJed Brown p += count; 1085e27a552bSJed Brown *p++ = ' '; 1086e27a552bSJed Brown } 1087e27a552bSJed Brown p[i ? 0 : -1] = 0; 1088e27a552bSJed Brown PetscFunctionReturn(0); 1089e27a552bSJed Brown } 1090e27a552bSJed Brown 1091e27a552bSJed Brown #undef __FUNCT__ 1092e27a552bSJed Brown #define __FUNCT__ "TSView_RosW" 1093e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer) 1094e27a552bSJed Brown { 109561692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 109661692a83SJed Brown RosWTableau tab = ros->tableau; 1097e27a552bSJed Brown PetscBool iascii; 1098e27a552bSJed Brown PetscErrorCode ierr; 1099e27a552bSJed Brown 1100e27a552bSJed Brown PetscFunctionBegin; 1101e27a552bSJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1102e27a552bSJed Brown if (iascii) { 110361692a83SJed Brown const TSRosWType rostype; 1104e408995aSJed Brown PetscInt i; 1105e408995aSJed Brown PetscReal abscissa[512]; 1106e27a552bSJed Brown char buf[512]; 110761692a83SJed Brown ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr); 110861692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Rosenbrock-W %s\n",rostype);CHKERRQ(ierr); 1109e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr); 111061692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A = %s\n",buf);CHKERRQ(ierr); 1111e408995aSJed Brown for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 1112e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr); 1113e408995aSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr); 1114e27a552bSJed Brown } 1115e27a552bSJed Brown ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 1116e27a552bSJed Brown PetscFunctionReturn(0); 1117e27a552bSJed Brown } 1118e27a552bSJed Brown 1119e27a552bSJed Brown #undef __FUNCT__ 1120e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType" 1121e27a552bSJed Brown /*@C 112261692a83SJed Brown TSRosWSetType - Set the type of Rosenbrock-W scheme 1123e27a552bSJed Brown 1124e27a552bSJed Brown Logically collective 1125e27a552bSJed Brown 1126e27a552bSJed Brown Input Parameter: 1127e27a552bSJed Brown + ts - timestepping context 112861692a83SJed Brown - rostype - type of Rosenbrock-W scheme 1129e27a552bSJed Brown 1130020d8f30SJed Brown Level: beginner 1131e27a552bSJed Brown 1132020d8f30SJed Brown .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3 1133e27a552bSJed Brown @*/ 113461692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype) 1135e27a552bSJed Brown { 1136e27a552bSJed Brown PetscErrorCode ierr; 1137e27a552bSJed Brown 1138e27a552bSJed Brown PetscFunctionBegin; 1139e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 114061692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr); 1141e27a552bSJed Brown PetscFunctionReturn(0); 1142e27a552bSJed Brown } 1143e27a552bSJed Brown 1144e27a552bSJed Brown #undef __FUNCT__ 1145e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType" 1146e27a552bSJed Brown /*@C 114761692a83SJed Brown TSRosWGetType - Get the type of Rosenbrock-W scheme 1148e27a552bSJed Brown 1149e27a552bSJed Brown Logically collective 1150e27a552bSJed Brown 1151e27a552bSJed Brown Input Parameter: 1152e27a552bSJed Brown . ts - timestepping context 1153e27a552bSJed Brown 1154e27a552bSJed Brown Output Parameter: 115561692a83SJed Brown . rostype - type of Rosenbrock-W scheme 1156e27a552bSJed Brown 1157e27a552bSJed Brown Level: intermediate 1158e27a552bSJed Brown 1159e27a552bSJed Brown .seealso: TSRosWGetType() 1160e27a552bSJed Brown @*/ 116161692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype) 1162e27a552bSJed Brown { 1163e27a552bSJed Brown PetscErrorCode ierr; 1164e27a552bSJed Brown 1165e27a552bSJed Brown PetscFunctionBegin; 1166e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 116761692a83SJed Brown ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr); 1168e27a552bSJed Brown PetscFunctionReturn(0); 1169e27a552bSJed Brown } 1170e27a552bSJed Brown 1171e27a552bSJed Brown #undef __FUNCT__ 117261692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian" 1173e27a552bSJed Brown /*@C 117461692a83SJed Brown TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 1175e27a552bSJed Brown 1176e27a552bSJed Brown Logically collective 1177e27a552bSJed Brown 1178e27a552bSJed Brown Input Parameter: 1179e27a552bSJed Brown + ts - timestepping context 118061692a83SJed Brown - flg - PETSC_TRUE to recompute the Jacobian at each stage 1181e27a552bSJed Brown 1182e27a552bSJed Brown Level: intermediate 1183e27a552bSJed Brown 1184e27a552bSJed Brown .seealso: TSRosWGetType() 1185e27a552bSJed Brown @*/ 118661692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg) 1187e27a552bSJed Brown { 1188e27a552bSJed Brown PetscErrorCode ierr; 1189e27a552bSJed Brown 1190e27a552bSJed Brown PetscFunctionBegin; 1191e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 119261692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1193e27a552bSJed Brown PetscFunctionReturn(0); 1194e27a552bSJed Brown } 1195e27a552bSJed Brown 1196e27a552bSJed Brown EXTERN_C_BEGIN 1197e27a552bSJed Brown #undef __FUNCT__ 1198e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW" 119961692a83SJed Brown PetscErrorCode TSRosWGetType_RosW(TS ts,const TSRosWType *rostype) 1200e27a552bSJed Brown { 120161692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1202e27a552bSJed Brown PetscErrorCode ierr; 1203e27a552bSJed Brown 1204e27a552bSJed Brown PetscFunctionBegin; 120561692a83SJed Brown if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);} 120661692a83SJed Brown *rostype = ros->tableau->name; 1207e27a552bSJed Brown PetscFunctionReturn(0); 1208e27a552bSJed Brown } 1209e27a552bSJed Brown #undef __FUNCT__ 1210e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW" 121161692a83SJed Brown PetscErrorCode TSRosWSetType_RosW(TS ts,const TSRosWType rostype) 1212e27a552bSJed Brown { 121361692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1214e27a552bSJed Brown PetscErrorCode ierr; 1215e27a552bSJed Brown PetscBool match; 121661692a83SJed Brown RosWTableauLink link; 1217e27a552bSJed Brown 1218e27a552bSJed Brown PetscFunctionBegin; 121961692a83SJed Brown if (ros->tableau) { 122061692a83SJed Brown ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr); 1221e27a552bSJed Brown if (match) PetscFunctionReturn(0); 1222e27a552bSJed Brown } 122361692a83SJed Brown for (link = RosWTableauList; link; link=link->next) { 122461692a83SJed Brown ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr); 1225e27a552bSJed Brown if (match) { 1226e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 122761692a83SJed Brown ros->tableau = &link->tab; 1228e27a552bSJed Brown PetscFunctionReturn(0); 1229e27a552bSJed Brown } 1230e27a552bSJed Brown } 123161692a83SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype); 1232e27a552bSJed Brown PetscFunctionReturn(0); 1233e27a552bSJed Brown } 123461692a83SJed Brown 1235e27a552bSJed Brown #undef __FUNCT__ 123661692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW" 123761692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg) 1238e27a552bSJed Brown { 123961692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1240e27a552bSJed Brown 1241e27a552bSJed Brown PetscFunctionBegin; 124261692a83SJed Brown ros->recompute_jacobian = flg; 1243e27a552bSJed Brown PetscFunctionReturn(0); 1244e27a552bSJed Brown } 1245e27a552bSJed Brown EXTERN_C_END 1246e27a552bSJed Brown 1247e27a552bSJed Brown /* ------------------------------------------------------------ */ 1248e27a552bSJed Brown /*MC 1249020d8f30SJed Brown TSROSW - ODE solver using Rosenbrock-W schemes 1250e27a552bSJed Brown 1251e27a552bSJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1252e27a552bSJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1253e27a552bSJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1254e27a552bSJed Brown 1255e27a552bSJed Brown Notes: 125661692a83SJed Brown This method currently only works with autonomous ODE and DAE. 125761692a83SJed Brown 125861692a83SJed Brown Developer notes: 125961692a83SJed Brown Rosenbrock-W methods are typically specified for autonomous ODE 126061692a83SJed Brown 126161692a83SJed Brown $ xdot = f(x) 126261692a83SJed Brown 126361692a83SJed Brown by the stage equations 126461692a83SJed Brown 126561692a83SJed Brown $ k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 126661692a83SJed Brown 126761692a83SJed Brown and step completion formula 126861692a83SJed Brown 126961692a83SJed Brown $ x_1 = x_0 + sum_j b_j k_j 127061692a83SJed Brown 127161692a83SJed Brown with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x) 127261692a83SJed Brown and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 127361692a83SJed Brown we define new variables for the stage equations 127461692a83SJed Brown 127561692a83SJed Brown $ y_i = gamma_ij k_j 127661692a83SJed Brown 127761692a83SJed Brown The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 127861692a83SJed Brown 127961692a83SJed Brown $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i} 128061692a83SJed Brown 128161692a83SJed Brown to rewrite the method as 128261692a83SJed Brown 128361692a83SJed 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 128461692a83SJed Brown $ x_1 = x_0 + sum_j bt_j y_j 128561692a83SJed Brown 128661692a83SJed Brown where we have introduced the mass matrix M. Continue by defining 128761692a83SJed Brown 128861692a83SJed Brown $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 128961692a83SJed Brown 129061692a83SJed Brown or, more compactly in tensor notation 129161692a83SJed Brown 129261692a83SJed Brown $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 129361692a83SJed Brown 129461692a83SJed Brown Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 129561692a83SJed Brown stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 129661692a83SJed Brown equation 129761692a83SJed Brown 129861692a83SJed Brown $ g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 129961692a83SJed Brown 130061692a83SJed Brown with initial guess y_i = 0. 1301e27a552bSJed Brown 1302e27a552bSJed Brown Level: beginner 1303e27a552bSJed Brown 1304020d8f30SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister() 1305e27a552bSJed Brown 1306e27a552bSJed Brown M*/ 1307e27a552bSJed Brown EXTERN_C_BEGIN 1308e27a552bSJed Brown #undef __FUNCT__ 1309e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW" 1310e27a552bSJed Brown PetscErrorCode TSCreate_RosW(TS ts) 1311e27a552bSJed Brown { 131261692a83SJed Brown TS_RosW *ros; 1313e27a552bSJed Brown PetscErrorCode ierr; 1314e27a552bSJed Brown 1315e27a552bSJed Brown PetscFunctionBegin; 1316e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1317e27a552bSJed Brown ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1318e27a552bSJed Brown #endif 1319e27a552bSJed Brown 1320e27a552bSJed Brown ts->ops->reset = TSReset_RosW; 1321e27a552bSJed Brown ts->ops->destroy = TSDestroy_RosW; 1322e27a552bSJed Brown ts->ops->view = TSView_RosW; 1323e27a552bSJed Brown ts->ops->setup = TSSetUp_RosW; 1324e27a552bSJed Brown ts->ops->step = TSStep_RosW; 1325e27a552bSJed Brown ts->ops->interpolate = TSInterpolate_RosW; 13261c3436cfSJed Brown ts->ops->evaluatestep = TSEvaluateStep_RosW; 1327e27a552bSJed Brown ts->ops->setfromoptions = TSSetFromOptions_RosW; 1328e27a552bSJed Brown ts->ops->snesfunction = SNESTSFormFunction_RosW; 1329e27a552bSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1330e27a552bSJed Brown 133161692a83SJed Brown ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr); 133261692a83SJed Brown ts->data = (void*)ros; 1333e27a552bSJed Brown 1334e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr); 1335e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr); 133661692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr); 1337e27a552bSJed Brown PetscFunctionReturn(0); 1338e27a552bSJed Brown } 1339e27a552bSJed Brown EXTERN_C_END 1340