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 1761692a83SJed Brown static const TSRosWType TSRosWDefault = TSROSW2P; 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 */ 2661692a83SJed Brown PetscReal *A; /* Propagation table, strictly lower triangular */ 2761692a83SJed Brown PetscReal *Gamma; /* Stage table, lower triangular with nonzero diagonal */ 28c17803e7SJed Brown PetscBool *GammaZeroDiag; /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */ 2943b21953SEmil Constantinescu PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/ 3061692a83SJed Brown PetscReal *b; /* Step completion table */ 31fe7e6d57SJed Brown PetscReal *bembed; /* Step completion table for embedded method of order one less */ 3261692a83SJed Brown PetscReal *ASum; /* Row sum of A */ 3361692a83SJed Brown PetscReal *GammaSum; /* Row sum of Gamma, only needed for non-autonomous systems */ 3461692a83SJed Brown PetscReal *At; /* Propagation table in transformed variables */ 3561692a83SJed Brown PetscReal *bt; /* Step completion table in transformed variables */ 36fe7e6d57SJed Brown PetscReal *bembedt; /* Step completion table of order one less in transformed variables */ 3761692a83SJed Brown PetscReal *GammaInv; /* Inverse of Gamma, used for transformed variables */ 388d59e960SJed Brown PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 39e27a552bSJed Brown }; 4061692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink; 4161692a83SJed Brown struct _RosWTableauLink { 4261692a83SJed Brown struct _RosWTableau tab; 4361692a83SJed Brown RosWTableauLink next; 44e27a552bSJed Brown }; 4561692a83SJed Brown static RosWTableauLink RosWTableauList; 46e27a552bSJed Brown 47e27a552bSJed Brown typedef struct { 4861692a83SJed Brown RosWTableau tableau; 4961692a83SJed Brown Vec *Y; /* States computed during the step, used to complete the step */ 50e27a552bSJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 5161692a83SJed Brown Vec Ystage; /* Work vector for the state value at each stage */ 5261692a83SJed Brown Vec Zdot; /* Ydot = Zdot + shift*Y */ 5361692a83SJed Brown Vec Zstage; /* Y = Zstage + Y */ 541c3436cfSJed Brown PetscScalar *work; /* Scalar work space of length number of stages, used to prepare VecMAXPY() */ 55e27a552bSJed Brown PetscReal shift; 56e27a552bSJed Brown PetscReal stage_time; 57c17803e7SJed Brown PetscReal stage_explicit; /* Flag indicates that the current stage is explicit */ 5861692a83SJed Brown PetscBool recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */ 59108c343cSJed Brown TSStepStatus status; 60e27a552bSJed Brown } TS_RosW; 61e27a552bSJed Brown 62fe7e6d57SJed Brown /*MC 63fe7e6d57SJed Brown TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme. 64fe7e6d57SJed Brown 65fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P. 66fe7e6d57SJed Brown 67fe7e6d57SJed Brown Level: intermediate 68fe7e6d57SJed Brown 69fe7e6d57SJed Brown .seealso: TSROSW 70fe7e6d57SJed Brown M*/ 71fe7e6d57SJed Brown 72fe7e6d57SJed Brown /*MC 73fe7e6d57SJed Brown TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme. 74fe7e6d57SJed Brown 75fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M. 76fe7e6d57SJed Brown 77fe7e6d57SJed Brown Level: intermediate 78fe7e6d57SJed Brown 79fe7e6d57SJed Brown .seealso: TSROSW 80fe7e6d57SJed Brown M*/ 81fe7e6d57SJed Brown 82fe7e6d57SJed Brown /*MC 83fe7e6d57SJed Brown TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1. 84fe7e6d57SJed Brown 85fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 86fe7e6d57SJed Brown 87fe7e6d57SJed 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. 88fe7e6d57SJed Brown 89fe7e6d57SJed Brown References: 90fe7e6d57SJed Brown Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 91fe7e6d57SJed Brown 92fe7e6d57SJed Brown Level: intermediate 93fe7e6d57SJed Brown 94fe7e6d57SJed Brown .seealso: TSROSW 95fe7e6d57SJed Brown M*/ 96fe7e6d57SJed Brown 97fe7e6d57SJed Brown /*MC 98fe7e6d57SJed Brown TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1. 99fe7e6d57SJed Brown 100fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 101fe7e6d57SJed Brown 102fe7e6d57SJed 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. 103fe7e6d57SJed Brown 104fe7e6d57SJed Brown References: 105fe7e6d57SJed Brown Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 106fe7e6d57SJed Brown 107fe7e6d57SJed Brown Level: intermediate 108fe7e6d57SJed Brown 109fe7e6d57SJed Brown .seealso: TSROSW 110fe7e6d57SJed Brown M*/ 111fe7e6d57SJed Brown 112ef3c5b88SJed Brown /*MC 113ef3c5b88SJed Brown TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme 114ef3c5b88SJed Brown 115ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 116ef3c5b88SJed Brown 117ef3c5b88SJed Brown Both the third order and embedded second order methods are stiffly accurate and L-stable. 118ef3c5b88SJed Brown 119ef3c5b88SJed Brown References: 120ef3c5b88SJed Brown Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 121ef3c5b88SJed Brown 122ef3c5b88SJed Brown Level: intermediate 123ef3c5b88SJed Brown 124ef3c5b88SJed Brown .seealso: TSROSW, TSROSWSANDU3 125ef3c5b88SJed Brown M*/ 126ef3c5b88SJed Brown 127ef3c5b88SJed Brown /*MC 128ef3c5b88SJed Brown TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme 129ef3c5b88SJed Brown 130ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 131ef3c5b88SJed Brown 132ef3c5b88SJed Brown The third order method is L-stable, but not stiffly accurate. 133ef3c5b88SJed Brown The second order embedded method is strongly A-stable with R(infty) = 0.5. 134ef3c5b88SJed Brown The internal stages are L-stable. 135ef3c5b88SJed Brown This method is called ROS3 in the paper. 136ef3c5b88SJed Brown 137ef3c5b88SJed Brown References: 138ef3c5b88SJed Brown Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 139ef3c5b88SJed Brown 140ef3c5b88SJed Brown Level: intermediate 141ef3c5b88SJed Brown 142ef3c5b88SJed Brown .seealso: TSROSW, TSROSWRODAS3 143ef3c5b88SJed Brown M*/ 144ef3c5b88SJed Brown 145961f28d0SJed Brown /*MC 146961f28d0SJed Brown TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages 147961f28d0SJed Brown 148961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 149961f28d0SJed Brown 150961f28d0SJed Brown A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3) 151961f28d0SJed Brown 152961f28d0SJed Brown References: 153961f28d0SJed Brown Emil Constantinescu 154961f28d0SJed Brown 155961f28d0SJed Brown Level: intermediate 156961f28d0SJed Brown 15743b21953SEmil Constantinescu .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP 158961f28d0SJed Brown M*/ 159961f28d0SJed Brown 160961f28d0SJed Brown /*MC 161961f28d0SJed Brown TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 162961f28d0SJed Brown 163961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 164961f28d0SJed Brown 165961f28d0SJed Brown L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 166961f28d0SJed Brown 167961f28d0SJed Brown References: 168961f28d0SJed Brown Emil Constantinescu 169961f28d0SJed Brown 170961f28d0SJed Brown Level: intermediate 171961f28d0SJed Brown 17243b21953SEmil Constantinescu .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP 173961f28d0SJed Brown M*/ 174961f28d0SJed Brown 175961f28d0SJed Brown /*MC 17643b21953SEmil Constantinescu TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 177961f28d0SJed Brown 178961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 179961f28d0SJed Brown 180961f28d0SJed Brown L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 181961f28d0SJed Brown 182961f28d0SJed Brown References: 183961f28d0SJed Brown Emil Constantinescu 184961f28d0SJed Brown 185961f28d0SJed Brown Level: intermediate 186961f28d0SJed Brown 187961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP 188961f28d0SJed Brown M*/ 189961f28d0SJed Brown 190e27a552bSJed Brown #undef __FUNCT__ 191e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll" 192e27a552bSJed Brown /*@C 193e27a552bSJed Brown TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW 194e27a552bSJed Brown 195e27a552bSJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 196e27a552bSJed Brown 197e27a552bSJed Brown Level: advanced 198e27a552bSJed Brown 199e27a552bSJed Brown .keywords: TS, TSRosW, register, all 200e27a552bSJed Brown 201e27a552bSJed Brown .seealso: TSRosWRegisterDestroy() 202e27a552bSJed Brown @*/ 203e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void) 204e27a552bSJed Brown { 205e27a552bSJed Brown PetscErrorCode ierr; 206e27a552bSJed Brown 207e27a552bSJed Brown PetscFunctionBegin; 208e27a552bSJed Brown if (TSRosWRegisterAllCalled) PetscFunctionReturn(0); 209e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_TRUE; 210e27a552bSJed Brown 211e27a552bSJed Brown { 21261692a83SJed Brown const PetscReal g = 1. + 1./PetscSqrtReal(2.0); 213e27a552bSJed Brown const PetscReal 21461692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 21561692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2161c3436cfSJed Brown b[2] = {0.5,0.5}, 2171c3436cfSJed Brown b1[2] = {1.0,0.0}; 2181c3436cfSJed Brown ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr); 219e27a552bSJed Brown } 220e27a552bSJed Brown { 22161692a83SJed Brown const PetscReal g = 1. - 1./PetscSqrtReal(2.0); 222e27a552bSJed Brown const PetscReal 22361692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 22461692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2251c3436cfSJed Brown b[2] = {0.5,0.5}, 2261c3436cfSJed Brown b1[2] = {1.0,0.0}; 2271c3436cfSJed Brown ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr); 228fe7e6d57SJed Brown } 229fe7e6d57SJed Brown { 230fe7e6d57SJed Brown const PetscReal g = 7.8867513459481287e-01; 231fe7e6d57SJed Brown const PetscReal 232fe7e6d57SJed Brown A[3][3] = {{0,0,0}, 233fe7e6d57SJed Brown {1.5773502691896257e+00,0,0}, 234fe7e6d57SJed Brown {0.5,0,0}}, 235fe7e6d57SJed Brown Gamma[3][3] = {{g,0,0}, 236fe7e6d57SJed Brown {-1.5773502691896257e+00,g,0}, 237fe7e6d57SJed Brown {-6.7075317547305480e-01,1.7075317547305482e-01,g}}, 238fe7e6d57SJed Brown b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01}, 239fe7e6d57SJed Brown b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01}; 240fe7e6d57SJed Brown ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 241fe7e6d57SJed Brown } 242fe7e6d57SJed Brown { 243fe7e6d57SJed Brown const PetscReal g = 4.3586652150845900e-01; 244fe7e6d57SJed Brown const PetscReal 245fe7e6d57SJed Brown A[4][4] = {{0,0,0,0}, 246fe7e6d57SJed Brown {8.7173304301691801e-01,0,0,0}, 247fe7e6d57SJed Brown {8.4457060015369423e-01,-1.1299064236484185e-01,0,0}, 248fe7e6d57SJed Brown {0,0,1.,0}}, 249fe7e6d57SJed Brown Gamma[4][4] = {{g,0,0,0}, 250fe7e6d57SJed Brown {-8.7173304301691801e-01,g,0,0}, 251fe7e6d57SJed Brown {-9.0338057013044082e-01,5.4180672388095326e-02,g,0}, 252fe7e6d57SJed Brown {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}}, 253fe7e6d57SJed Brown b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01}, 254fe7e6d57SJed Brown b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01}; 255fe7e6d57SJed Brown ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 256e27a552bSJed Brown } 257ef3c5b88SJed Brown { 258ef3c5b88SJed Brown const PetscReal g = 0.5; 259ef3c5b88SJed Brown const PetscReal 260ef3c5b88SJed Brown A[4][4] = {{0,0,0,0}, 261ef3c5b88SJed Brown {0,0,0,0}, 262ef3c5b88SJed Brown {1.,0,0,0}, 263ef3c5b88SJed Brown {0.75,-0.25,0.5,0}}, 264ef3c5b88SJed Brown Gamma[4][4] = {{g,0,0,0}, 265ef3c5b88SJed Brown {1.,g,0,0}, 266ef3c5b88SJed Brown {-0.25,-0.25,g,0}, 267ef3c5b88SJed Brown {1./12,1./12,-2./3,g}}, 268ef3c5b88SJed Brown b[4] = {5./6,-1./6,-1./6,0.5}, 269ef3c5b88SJed Brown b2[4] = {0.75,-0.25,0.5,0}; 270ef3c5b88SJed Brown ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 271ef3c5b88SJed Brown } 272ef3c5b88SJed Brown { 273ef3c5b88SJed Brown const PetscReal g = 0.43586652150845899941601945119356; 274ef3c5b88SJed Brown const PetscReal 275ef3c5b88SJed Brown A[3][3] = {{0,0,0}, 276ef3c5b88SJed Brown {g,0,0}, 277ef3c5b88SJed Brown {g,0,0}}, 278ef3c5b88SJed Brown Gamma[3][3] = {{g,0,0}, 279ef3c5b88SJed Brown {-0.19294655696029095575009695436041,g,0}, 280ef3c5b88SJed Brown {0,1.74927148125794685173529749738960,g}}, 281ef3c5b88SJed Brown b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829}, 282ef3c5b88SJed Brown b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619}; 283ef3c5b88SJed Brown ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 284ef3c5b88SJed Brown } 285b1c69cc3SEmil Constantinescu { 286*aaf9cf16SJed Brown const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0; 287b1c69cc3SEmil Constantinescu const PetscReal 288b1c69cc3SEmil Constantinescu A[3][3] = {{0,0,0}, 289b1c69cc3SEmil Constantinescu {1,0,0}, 290b1c69cc3SEmil Constantinescu {0.25,0.25,0}}, 291b1c69cc3SEmil Constantinescu Gamma[3][3] = {{0,0,0}, 292*aaf9cf16SJed Brown {(-3.0-s3)/6.0,g,0}, 293*aaf9cf16SJed Brown {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}}, 294b1c69cc3SEmil Constantinescu b[3] = {1./6.,1./6.,2./3.}, 295b1c69cc3SEmil Constantinescu b2[3] = {1./4.,1./4.,1./2.}; 296b1c69cc3SEmil Constantinescu ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 297b1c69cc3SEmil Constantinescu } 298b1c69cc3SEmil Constantinescu 299b1c69cc3SEmil Constantinescu { 300b1c69cc3SEmil Constantinescu const PetscReal 301b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 302b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 303b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 304b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 305b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 306b1c69cc3SEmil Constantinescu {0.0,1./4.,0,0}, 307b1c69cc3SEmil Constantinescu {-2.,-2./3.,2./3.,0}, 308b1c69cc3SEmil Constantinescu {1./2.,5./36.,-2./9,0}}, 309b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 310b1c69cc3SEmil Constantinescu b2[4] = {1./8.,3./4.,1./8.,0}; 311b1c69cc3SEmil Constantinescu ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 312b1c69cc3SEmil Constantinescu } 313b1c69cc3SEmil Constantinescu 314b1c69cc3SEmil Constantinescu { 315b1c69cc3SEmil Constantinescu const PetscReal 316b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 317b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 318b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 319b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 320b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 321b1c69cc3SEmil Constantinescu {0.0,3./4.,0,0}, 322b1c69cc3SEmil Constantinescu {-2./3.,-23./9.,2./9.,0}, 323b1c69cc3SEmil Constantinescu {1./18.,65./108.,-2./27,0}}, 324b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 325b1c69cc3SEmil Constantinescu b2[4] = {3./16.,10./16.,3./16.,0}; 32643b21953SEmil Constantinescu ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 327b1c69cc3SEmil Constantinescu } 328753f8adbSEmil Constantinescu 329753f8adbSEmil Constantinescu { 330753f8adbSEmil Constantinescu PetscReal A[4][4],Gamma[4][4],b[4],b2[4]; 331753f8adbSEmil Constantinescu 332753f8adbSEmil Constantinescu Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816; 33305e8e825SJed Brown Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0; 334753f8adbSEmil Constantinescu Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476; 335753f8adbSEmil Constantinescu Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816; 33605e8e825SJed Brown Gamma[1][2]=0; Gamma[1][3]=0; 337753f8adbSEmil Constantinescu Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903; 338753f8adbSEmil Constantinescu Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131; 339753f8adbSEmil Constantinescu Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816; 34005e8e825SJed Brown Gamma[2][3]=0; 341753f8adbSEmil Constantinescu Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783; 342753f8adbSEmil Constantinescu Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984; 343753f8adbSEmil Constantinescu Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198; 344753f8adbSEmil Constantinescu Gamma[3][3]=0; 345753f8adbSEmil Constantinescu 34605e8e825SJed Brown A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0; 347753f8adbSEmil Constantinescu A[1][0]=0.8717330430169179988320388950590125027645343373957631; 34805e8e825SJed Brown A[1][1]=0; A[1][2]=0; A[1][3]=0; 349753f8adbSEmil Constantinescu A[2][0]=0.5275890119763004115618079766722914408876108660811028; 350753f8adbSEmil Constantinescu A[2][1]=0.07241098802369958843819203208518599088698057726988732; 35105e8e825SJed Brown A[2][2]=0; A[2][3]=0; 352753f8adbSEmil Constantinescu A[3][0]=0.3990960076760701320627260685975778145384666450351314; 353753f8adbSEmil Constantinescu A[3][1]=-0.4375576546135194437228463747348862825846903771419953; 354753f8adbSEmil Constantinescu A[3][2]=1.038461646937449311660120300601880176655352737312713; 35505e8e825SJed Brown A[3][3]=0; 356753f8adbSEmil Constantinescu 357753f8adbSEmil Constantinescu b[0]=0.1876410243467238251612921333138006734899663569186926; 358753f8adbSEmil Constantinescu b[1]=-0.5952974735769549480478230473706443582188442040780541; 359753f8adbSEmil Constantinescu b[2]=0.9717899277217721234705114616271378792182450260943198; 360753f8adbSEmil Constantinescu b[3]=0.4358665215084589994160194475295062513822671686978816; 361753f8adbSEmil Constantinescu 362753f8adbSEmil Constantinescu b2[0]=0.2147402862233891404862383521089097657790734483804460; 363753f8adbSEmil Constantinescu b2[1]=-0.4851622638849390928209050538171743017757490232519684; 364753f8adbSEmil Constantinescu b2[2]=0.8687250025203875511662123688667549217531982787600080; 365753f8adbSEmil Constantinescu b2[3]=0.4016969751411624011684543450940068201770721128357014; 366753f8adbSEmil Constantinescu 367753f8adbSEmil Constantinescu ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 368753f8adbSEmil Constantinescu } 369753f8adbSEmil Constantinescu 370e27a552bSJed Brown PetscFunctionReturn(0); 371e27a552bSJed Brown } 372e27a552bSJed Brown 373e27a552bSJed Brown #undef __FUNCT__ 374e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy" 375e27a552bSJed Brown /*@C 376e27a552bSJed Brown TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 377e27a552bSJed Brown 378e27a552bSJed Brown Not Collective 379e27a552bSJed Brown 380e27a552bSJed Brown Level: advanced 381e27a552bSJed Brown 382e27a552bSJed Brown .keywords: TSRosW, register, destroy 383e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 384e27a552bSJed Brown @*/ 385e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void) 386e27a552bSJed Brown { 387e27a552bSJed Brown PetscErrorCode ierr; 38861692a83SJed Brown RosWTableauLink link; 389e27a552bSJed Brown 390e27a552bSJed Brown PetscFunctionBegin; 39161692a83SJed Brown while ((link = RosWTableauList)) { 39261692a83SJed Brown RosWTableau t = &link->tab; 39361692a83SJed Brown RosWTableauList = link->next; 39461692a83SJed Brown ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr); 39543b21953SEmil Constantinescu ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr); 396fe7e6d57SJed Brown ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr); 397e27a552bSJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 398e27a552bSJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 399e27a552bSJed Brown } 400e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_FALSE; 401e27a552bSJed Brown PetscFunctionReturn(0); 402e27a552bSJed Brown } 403e27a552bSJed Brown 404e27a552bSJed Brown #undef __FUNCT__ 405e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage" 406e27a552bSJed Brown /*@C 407e27a552bSJed Brown TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called 408e27a552bSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW() 409e27a552bSJed Brown when using static libraries. 410e27a552bSJed Brown 411e27a552bSJed Brown Input Parameter: 412e27a552bSJed Brown path - The dynamic library path, or PETSC_NULL 413e27a552bSJed Brown 414e27a552bSJed Brown Level: developer 415e27a552bSJed Brown 416e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package 417e27a552bSJed Brown .seealso: PetscInitialize() 418e27a552bSJed Brown @*/ 419e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[]) 420e27a552bSJed Brown { 421e27a552bSJed Brown PetscErrorCode ierr; 422e27a552bSJed Brown 423e27a552bSJed Brown PetscFunctionBegin; 424e27a552bSJed Brown if (TSRosWPackageInitialized) PetscFunctionReturn(0); 425e27a552bSJed Brown TSRosWPackageInitialized = PETSC_TRUE; 426e27a552bSJed Brown ierr = TSRosWRegisterAll();CHKERRQ(ierr); 427e27a552bSJed Brown ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr); 428e27a552bSJed Brown PetscFunctionReturn(0); 429e27a552bSJed Brown } 430e27a552bSJed Brown 431e27a552bSJed Brown #undef __FUNCT__ 432e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage" 433e27a552bSJed Brown /*@C 434e27a552bSJed Brown TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is 435e27a552bSJed Brown called from PetscFinalize(). 436e27a552bSJed Brown 437e27a552bSJed Brown Level: developer 438e27a552bSJed Brown 439e27a552bSJed Brown .keywords: Petsc, destroy, package 440e27a552bSJed Brown .seealso: PetscFinalize() 441e27a552bSJed Brown @*/ 442e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void) 443e27a552bSJed Brown { 444e27a552bSJed Brown PetscErrorCode ierr; 445e27a552bSJed Brown 446e27a552bSJed Brown PetscFunctionBegin; 447e27a552bSJed Brown TSRosWPackageInitialized = PETSC_FALSE; 448e27a552bSJed Brown ierr = TSRosWRegisterDestroy();CHKERRQ(ierr); 449e27a552bSJed Brown PetscFunctionReturn(0); 450e27a552bSJed Brown } 451e27a552bSJed Brown 452e27a552bSJed Brown #undef __FUNCT__ 453e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister" 454e27a552bSJed Brown /*@C 45561692a83SJed Brown TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 456e27a552bSJed Brown 457e27a552bSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 458e27a552bSJed Brown 459e27a552bSJed Brown Input Parameters: 460e27a552bSJed Brown + name - identifier for method 461e27a552bSJed Brown . order - approximation order of method 462e27a552bSJed Brown . s - number of stages, this is the dimension of the matrices below 46361692a83SJed Brown . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 46461692a83SJed Brown . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 465fe7e6d57SJed Brown . b - Step completion table (dimension s) 466fe7e6d57SJed Brown - bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available) 467e27a552bSJed Brown 468e27a552bSJed Brown Notes: 46961692a83SJed Brown Several Rosenbrock W methods are provided, this function is only needed to create new methods. 470e27a552bSJed Brown 471e27a552bSJed Brown Level: advanced 472e27a552bSJed Brown 473e27a552bSJed Brown .keywords: TS, register 474e27a552bSJed Brown 475e27a552bSJed Brown .seealso: TSRosW 476e27a552bSJed Brown @*/ 477e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s, 478fe7e6d57SJed Brown const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[]) 479e27a552bSJed Brown { 480e27a552bSJed Brown PetscErrorCode ierr; 48161692a83SJed Brown RosWTableauLink link; 48261692a83SJed Brown RosWTableau t; 48361692a83SJed Brown PetscInt i,j,k; 48461692a83SJed Brown PetscScalar *GammaInv; 485e27a552bSJed Brown 486e27a552bSJed Brown PetscFunctionBegin; 487fe7e6d57SJed Brown PetscValidCharPointer(name,1); 488fe7e6d57SJed Brown PetscValidPointer(A,4); 489fe7e6d57SJed Brown PetscValidPointer(Gamma,5); 490fe7e6d57SJed Brown PetscValidPointer(b,6); 491fe7e6d57SJed Brown if (bembed) PetscValidPointer(bembed,7); 492fe7e6d57SJed Brown 493e27a552bSJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 494e27a552bSJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 495e27a552bSJed Brown t = &link->tab; 496e27a552bSJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 497e27a552bSJed Brown t->order = order; 498e27a552bSJed Brown t->s = s; 49961692a83SJed 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); 50043b21953SEmil 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); 501e27a552bSJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 50261692a83SJed Brown ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 50343b21953SEmil Constantinescu ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 50461692a83SJed Brown ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 505fe7e6d57SJed Brown if (bembed) { 506fe7e6d57SJed Brown ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr); 507fe7e6d57SJed Brown ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 508fe7e6d57SJed Brown } 50961692a83SJed Brown for (i=0; i<s; i++) { 51061692a83SJed Brown t->ASum[i] = 0; 51161692a83SJed Brown t->GammaSum[i] = 0; 51261692a83SJed Brown for (j=0; j<s; j++) { 51361692a83SJed Brown t->ASum[i] += A[i*s+j]; 514fe7e6d57SJed Brown t->GammaSum[i] += Gamma[i*s+j]; 51561692a83SJed Brown } 51661692a83SJed Brown } 51761692a83SJed Brown ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */ 51861692a83SJed Brown for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i]; 519fd96d5b0SEmil Constantinescu for (i=0; i<s; i++) { 520fd96d5b0SEmil Constantinescu if (Gamma[i*s+i] == 0.0) { 521fd96d5b0SEmil Constantinescu GammaInv[i*s+i] = 1.0; 522c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_TRUE; 523fd96d5b0SEmil Constantinescu } else { 524c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_FALSE; 525fd96d5b0SEmil Constantinescu } 526fd96d5b0SEmil Constantinescu } 527fd96d5b0SEmil Constantinescu 52861692a83SJed Brown switch (s) { 52961692a83SJed Brown case 1: GammaInv[0] = 1./GammaInv[0]; break; 53061692a83SJed Brown case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break; 53161692a83SJed Brown case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break; 53261692a83SJed Brown case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break; 53361692a83SJed Brown case 5: { 53461692a83SJed Brown PetscInt ipvt5[5]; 53561692a83SJed Brown MatScalar work5[5*5]; 53661692a83SJed Brown ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break; 53761692a83SJed Brown } 53861692a83SJed Brown case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break; 53961692a83SJed Brown case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break; 54061692a83SJed Brown default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s); 54161692a83SJed Brown } 54261692a83SJed Brown for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 54361692a83SJed Brown ierr = PetscFree(GammaInv);CHKERRQ(ierr); 54443b21953SEmil Constantinescu 54543b21953SEmil Constantinescu for (i=0; i<s; i++) { 54643b21953SEmil Constantinescu for (k=0; k<i+1; k++) { 54743b21953SEmil Constantinescu t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]); 54843b21953SEmil Constantinescu for (j=k+1; j<i+1; j++) { 54943b21953SEmil Constantinescu t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]); 55043b21953SEmil Constantinescu } 55143b21953SEmil Constantinescu } 55243b21953SEmil Constantinescu } 55343b21953SEmil Constantinescu 55461692a83SJed Brown for (i=0; i<s; i++) { 55561692a83SJed Brown for (j=0; j<s; j++) { 55661692a83SJed Brown t->At[i*s+j] = 0; 55761692a83SJed Brown for (k=0; k<s; k++) { 55861692a83SJed Brown t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j]; 55961692a83SJed Brown } 56061692a83SJed Brown } 56161692a83SJed Brown t->bt[i] = 0; 56261692a83SJed Brown for (j=0; j<s; j++) { 56361692a83SJed Brown t->bt[i] += t->b[j] * t->GammaInv[j*s+i]; 56461692a83SJed Brown } 565fe7e6d57SJed Brown if (bembed) { 566fe7e6d57SJed Brown t->bembedt[i] = 0; 567fe7e6d57SJed Brown for (j=0; j<s; j++) { 568fe7e6d57SJed Brown t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i]; 569fe7e6d57SJed Brown } 570fe7e6d57SJed Brown } 57161692a83SJed Brown } 5728d59e960SJed Brown t->ccfl = 1.0; /* Fix this */ 5738d59e960SJed Brown 57461692a83SJed Brown link->next = RosWTableauList; 57561692a83SJed Brown RosWTableauList = link; 576e27a552bSJed Brown PetscFunctionReturn(0); 577e27a552bSJed Brown } 578e27a552bSJed Brown 579e27a552bSJed Brown #undef __FUNCT__ 5801c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW" 5811c3436cfSJed Brown /* 5821c3436cfSJed Brown The step completion formula is 5831c3436cfSJed Brown 5841c3436cfSJed Brown x1 = x0 + b^T Y 5851c3436cfSJed Brown 5861c3436cfSJed Brown where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 5871c3436cfSJed Brown updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 5881c3436cfSJed Brown 5891c3436cfSJed Brown x1e = x0 + be^T Y 5901c3436cfSJed Brown = x1 - b^T Y + be^T Y 5911c3436cfSJed Brown = x1 + (be - b)^T Y 5921c3436cfSJed Brown 5931c3436cfSJed Brown so we can evaluate the method of different order even after the step has been optimistically completed. 5941c3436cfSJed Brown */ 5951c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done) 5961c3436cfSJed Brown { 5971c3436cfSJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 5981c3436cfSJed Brown RosWTableau tab = ros->tableau; 5991c3436cfSJed Brown PetscScalar *w = ros->work; 6001c3436cfSJed Brown PetscInt i; 6011c3436cfSJed Brown PetscErrorCode ierr; 6021c3436cfSJed Brown 6031c3436cfSJed Brown PetscFunctionBegin; 6041c3436cfSJed Brown if (order == tab->order) { 605108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */ 6061c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 607de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bt[i]; 608de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 609108c343cSJed Brown } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 6101c3436cfSJed Brown if (done) *done = PETSC_TRUE; 6111c3436cfSJed Brown PetscFunctionReturn(0); 6121c3436cfSJed Brown } else if (order == tab->order-1) { 6131c3436cfSJed Brown if (!tab->bembedt) goto unavailable; 614108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */ 6151c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 616de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i]; 617de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 618108c343cSJed Brown } else { /* Use rollback-and-recomplete formula (bembedt - bt) */ 619108c343cSJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 620108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 621108c343cSJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 6221c3436cfSJed Brown } 6231c3436cfSJed Brown if (done) *done = PETSC_TRUE; 6241c3436cfSJed Brown PetscFunctionReturn(0); 6251c3436cfSJed Brown } 6261c3436cfSJed Brown unavailable: 6271c3436cfSJed Brown if (done) *done = PETSC_FALSE; 6281c3436cfSJed 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); 6291c3436cfSJed Brown PetscFunctionReturn(0); 6301c3436cfSJed Brown } 6311c3436cfSJed Brown 6321c3436cfSJed Brown #undef __FUNCT__ 633e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW" 634e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts) 635e27a552bSJed Brown { 63661692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 63761692a83SJed Brown RosWTableau tab = ros->tableau; 638e27a552bSJed Brown const PetscInt s = tab->s; 6391c3436cfSJed Brown const PetscReal *At = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv; 640c17803e7SJed Brown const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 64161692a83SJed Brown PetscScalar *w = ros->work; 6427d4bf2deSEmil Constantinescu Vec *Y = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage; 643e27a552bSJed Brown SNES snes; 6441c3436cfSJed Brown TSAdapt adapt; 6451c3436cfSJed Brown PetscInt i,j,its,lits,reject,next_scheme; 646cdbf8f93SLisandro Dalcin PetscReal next_time_step; 6471c3436cfSJed Brown PetscBool accept; 648e27a552bSJed Brown PetscErrorCode ierr; 649e27a552bSJed Brown 650e27a552bSJed Brown PetscFunctionBegin; 651e27a552bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 652cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 6531c3436cfSJed Brown accept = PETSC_TRUE; 654108c343cSJed Brown ros->status = TS_STEP_INCOMPLETE; 655e27a552bSJed Brown 6561c3436cfSJed Brown for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 6571c3436cfSJed Brown const PetscReal h = ts->time_step; 658e27a552bSJed Brown for (i=0; i<s; i++) { 6591c3436cfSJed Brown ros->stage_time = ts->ptime + h*ASum[i]; 660c17803e7SJed Brown if (GammaZeroDiag[i]) { 661c17803e7SJed Brown ros->stage_explicit = PETSC_TRUE; 662fd96d5b0SEmil Constantinescu ros->shift = 1./h; 663c17803e7SJed Brown } else { 664c17803e7SJed Brown ros->stage_explicit = PETSC_FALSE; 66561692a83SJed Brown ros->shift = 1./(h*Gamma[i*s+i]); 666fd96d5b0SEmil Constantinescu } 66761692a83SJed Brown 66861692a83SJed Brown ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr); 669de19f811SJed Brown for (j=0; j<i; j++) w[j] = At[i*s+j]; 670de19f811SJed Brown ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 67161692a83SJed Brown 67261692a83SJed Brown for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j]; 67361692a83SJed Brown ierr = VecZeroEntries(Zdot);CHKERRQ(ierr); 67461692a83SJed Brown ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr); 67561692a83SJed Brown 676e27a552bSJed Brown /* Initial guess taken from last stage */ 67761692a83SJed Brown ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr); 67861692a83SJed Brown 6797d4bf2deSEmil Constantinescu if (!ros->stage_explicit) { 68061692a83SJed Brown if (!ros->recompute_jacobian && !i) { 68161692a83SJed Brown ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 68261692a83SJed Brown } 68361692a83SJed Brown ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr); 684e27a552bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 685e27a552bSJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 686e27a552bSJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 6877d4bf2deSEmil Constantinescu } else { 6887d4bf2deSEmil Constantinescu ierr = VecWAXPY(Ydot,1,ts->vec_sol,Zdot);CHKERRQ(ierr); /* Ydot = x0 + Zdot */ 6897d4bf2deSEmil Constantinescu ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,Ydot,Zdot,PETSC_FALSE);CHKERRQ(ierr); 6907d4bf2deSEmil Constantinescu ierr = VecWAXPY(ros->Ystage,1.0,Zdot,ros->Zstage);CHKERRQ(ierr); /* Ystage = F + Zstage */ 6917d4bf2deSEmil Constantinescu ts->linear_its += 1; 6927d4bf2deSEmil Constantinescu } 693e27a552bSJed Brown } 6941c3436cfSJed Brown ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 695108c343cSJed Brown ros->status = TS_STEP_PENDING; 696e27a552bSJed Brown 6971c3436cfSJed Brown /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 6981c3436cfSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 6991c3436cfSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 7008d59e960SJed Brown ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 7011c3436cfSJed Brown ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 7021c3436cfSJed Brown if (accept) { 7031c3436cfSJed Brown /* ignore next_scheme for now */ 704e27a552bSJed Brown ts->ptime += ts->time_step; 705cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 706e27a552bSJed Brown ts->steps++; 707108c343cSJed Brown ros->status = TS_STEP_COMPLETE; 7081c3436cfSJed Brown break; 7091c3436cfSJed Brown } else { /* Roll back the current step */ 7101c3436cfSJed Brown for (i=0; i<s; i++) w[i] = -tab->bt[i]; 7111c3436cfSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr); 7121c3436cfSJed Brown ts->time_step = next_time_step; 713108c343cSJed Brown ros->status = TS_STEP_INCOMPLETE; 7141c3436cfSJed Brown } 7151c3436cfSJed Brown } 716e27a552bSJed Brown PetscFunctionReturn(0); 717e27a552bSJed Brown } 718e27a552bSJed Brown 719e27a552bSJed Brown #undef __FUNCT__ 720e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW" 721e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X) 722e27a552bSJed Brown { 72361692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 724e27a552bSJed Brown 725e27a552bSJed Brown PetscFunctionBegin; 72661692a83SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name); 727e27a552bSJed Brown PetscFunctionReturn(0); 728e27a552bSJed Brown } 729e27a552bSJed Brown 730e27a552bSJed Brown /*------------------------------------------------------------*/ 731e27a552bSJed Brown #undef __FUNCT__ 732e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW" 733e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts) 734e27a552bSJed Brown { 73561692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 736e27a552bSJed Brown PetscInt s; 737e27a552bSJed Brown PetscErrorCode ierr; 738e27a552bSJed Brown 739e27a552bSJed Brown PetscFunctionBegin; 74061692a83SJed Brown if (!ros->tableau) PetscFunctionReturn(0); 74161692a83SJed Brown s = ros->tableau->s; 74261692a83SJed Brown ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr); 74361692a83SJed Brown ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr); 74461692a83SJed Brown ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr); 74561692a83SJed Brown ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr); 74661692a83SJed Brown ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr); 74761692a83SJed Brown ierr = PetscFree(ros->work);CHKERRQ(ierr); 748e27a552bSJed Brown PetscFunctionReturn(0); 749e27a552bSJed Brown } 750e27a552bSJed Brown 751e27a552bSJed Brown #undef __FUNCT__ 752e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW" 753e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts) 754e27a552bSJed Brown { 755e27a552bSJed Brown PetscErrorCode ierr; 756e27a552bSJed Brown 757e27a552bSJed Brown PetscFunctionBegin; 758e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 759e27a552bSJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 760e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr); 761e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr); 76261692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr); 763e27a552bSJed Brown PetscFunctionReturn(0); 764e27a552bSJed Brown } 765e27a552bSJed Brown 766e27a552bSJed Brown /* 767e27a552bSJed Brown This defines the nonlinear equation that is to be solved with SNES 768e27a552bSJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 769e27a552bSJed Brown */ 770e27a552bSJed Brown #undef __FUNCT__ 771e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW" 772e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts) 773e27a552bSJed Brown { 77461692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 775e27a552bSJed Brown PetscErrorCode ierr; 776e27a552bSJed Brown 777e27a552bSJed Brown PetscFunctionBegin; 77861692a83SJed Brown ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */ 77961692a83SJed Brown ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr); /* Ystage = X + Zstage */ 78061692a83SJed Brown ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr); 781e27a552bSJed Brown PetscFunctionReturn(0); 782e27a552bSJed Brown } 783e27a552bSJed Brown 784e27a552bSJed Brown #undef __FUNCT__ 785e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW" 786e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 787e27a552bSJed Brown { 78861692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 789e27a552bSJed Brown PetscErrorCode ierr; 790e27a552bSJed Brown 791e27a552bSJed Brown PetscFunctionBegin; 79261692a83SJed Brown /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 79361692a83SJed Brown ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 794e27a552bSJed Brown PetscFunctionReturn(0); 795e27a552bSJed Brown } 796e27a552bSJed Brown 797e27a552bSJed Brown #undef __FUNCT__ 798e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW" 799e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts) 800e27a552bSJed Brown { 80161692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 80261692a83SJed Brown RosWTableau tab = ros->tableau; 803e27a552bSJed Brown PetscInt s = tab->s; 804e27a552bSJed Brown PetscErrorCode ierr; 805e27a552bSJed Brown 806e27a552bSJed Brown PetscFunctionBegin; 80761692a83SJed Brown if (!ros->tableau) { 808e27a552bSJed Brown ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr); 809e27a552bSJed Brown } 81061692a83SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr); 81161692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr); 81261692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr); 81361692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr); 81461692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr); 81561692a83SJed Brown ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr); 816e27a552bSJed Brown PetscFunctionReturn(0); 817e27a552bSJed Brown } 818e27a552bSJed Brown /*------------------------------------------------------------*/ 819e27a552bSJed Brown 820e27a552bSJed Brown #undef __FUNCT__ 821e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW" 822e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts) 823e27a552bSJed Brown { 82461692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 825e27a552bSJed Brown PetscErrorCode ierr; 82661692a83SJed Brown char rostype[256]; 827e27a552bSJed Brown 828e27a552bSJed Brown PetscFunctionBegin; 829e27a552bSJed Brown ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr); 830e27a552bSJed Brown { 83161692a83SJed Brown RosWTableauLink link; 832e27a552bSJed Brown PetscInt count,choice; 833e27a552bSJed Brown PetscBool flg; 834e27a552bSJed Brown const char **namelist; 83561692a83SJed Brown SNES snes; 83661692a83SJed Brown 83761692a83SJed Brown ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr); 83861692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) ; 839e27a552bSJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 84061692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 84161692a83SJed Brown ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr); 84261692a83SJed Brown ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr); 843e27a552bSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 84461692a83SJed Brown 84561692a83SJed Brown ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr); 84661692a83SJed Brown 84761692a83SJed Brown /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 84861692a83SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 84961692a83SJed Brown if (!((PetscObject)snes)->type_name) { 85061692a83SJed Brown ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 85161692a83SJed Brown } 85261692a83SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 853e27a552bSJed Brown } 854e27a552bSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 855e27a552bSJed Brown PetscFunctionReturn(0); 856e27a552bSJed Brown } 857e27a552bSJed Brown 858e27a552bSJed Brown #undef __FUNCT__ 859e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray" 860e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 861e27a552bSJed Brown { 862e27a552bSJed Brown PetscErrorCode ierr; 863e408995aSJed Brown PetscInt i; 864e408995aSJed Brown size_t left,count; 865e27a552bSJed Brown char *p; 866e27a552bSJed Brown 867e27a552bSJed Brown PetscFunctionBegin; 868e408995aSJed Brown for (i=0,p=buf,left=len; i<n; i++) { 869e408995aSJed Brown ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 870e27a552bSJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 871e27a552bSJed Brown left -= count; 872e27a552bSJed Brown p += count; 873e27a552bSJed Brown *p++ = ' '; 874e27a552bSJed Brown } 875e27a552bSJed Brown p[i ? 0 : -1] = 0; 876e27a552bSJed Brown PetscFunctionReturn(0); 877e27a552bSJed Brown } 878e27a552bSJed Brown 879e27a552bSJed Brown #undef __FUNCT__ 880e27a552bSJed Brown #define __FUNCT__ "TSView_RosW" 881e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer) 882e27a552bSJed Brown { 88361692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 88461692a83SJed Brown RosWTableau tab = ros->tableau; 885e27a552bSJed Brown PetscBool iascii; 886e27a552bSJed Brown PetscErrorCode ierr; 887e27a552bSJed Brown 888e27a552bSJed Brown PetscFunctionBegin; 889e27a552bSJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 890e27a552bSJed Brown if (iascii) { 89161692a83SJed Brown const TSRosWType rostype; 892e408995aSJed Brown PetscInt i; 893e408995aSJed Brown PetscReal abscissa[512]; 894e27a552bSJed Brown char buf[512]; 89561692a83SJed Brown ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr); 89661692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Rosenbrock-W %s\n",rostype);CHKERRQ(ierr); 897e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr); 89861692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A = %s\n",buf);CHKERRQ(ierr); 899e408995aSJed Brown for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 900e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr); 901e408995aSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr); 902e27a552bSJed Brown } 903e27a552bSJed Brown ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 904e27a552bSJed Brown PetscFunctionReturn(0); 905e27a552bSJed Brown } 906e27a552bSJed Brown 907e27a552bSJed Brown #undef __FUNCT__ 908e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType" 909e27a552bSJed Brown /*@C 91061692a83SJed Brown TSRosWSetType - Set the type of Rosenbrock-W scheme 911e27a552bSJed Brown 912e27a552bSJed Brown Logically collective 913e27a552bSJed Brown 914e27a552bSJed Brown Input Parameter: 915e27a552bSJed Brown + ts - timestepping context 91661692a83SJed Brown - rostype - type of Rosenbrock-W scheme 917e27a552bSJed Brown 918e27a552bSJed Brown Level: intermediate 919e27a552bSJed Brown 920e27a552bSJed Brown .seealso: TSRosWGetType() 921e27a552bSJed Brown @*/ 92261692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype) 923e27a552bSJed Brown { 924e27a552bSJed Brown PetscErrorCode ierr; 925e27a552bSJed Brown 926e27a552bSJed Brown PetscFunctionBegin; 927e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 92861692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr); 929e27a552bSJed Brown PetscFunctionReturn(0); 930e27a552bSJed Brown } 931e27a552bSJed Brown 932e27a552bSJed Brown #undef __FUNCT__ 933e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType" 934e27a552bSJed Brown /*@C 93561692a83SJed Brown TSRosWGetType - Get the type of Rosenbrock-W scheme 936e27a552bSJed Brown 937e27a552bSJed Brown Logically collective 938e27a552bSJed Brown 939e27a552bSJed Brown Input Parameter: 940e27a552bSJed Brown . ts - timestepping context 941e27a552bSJed Brown 942e27a552bSJed Brown Output Parameter: 94361692a83SJed Brown . rostype - type of Rosenbrock-W scheme 944e27a552bSJed Brown 945e27a552bSJed Brown Level: intermediate 946e27a552bSJed Brown 947e27a552bSJed Brown .seealso: TSRosWGetType() 948e27a552bSJed Brown @*/ 94961692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype) 950e27a552bSJed Brown { 951e27a552bSJed Brown PetscErrorCode ierr; 952e27a552bSJed Brown 953e27a552bSJed Brown PetscFunctionBegin; 954e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 95561692a83SJed Brown ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr); 956e27a552bSJed Brown PetscFunctionReturn(0); 957e27a552bSJed Brown } 958e27a552bSJed Brown 959e27a552bSJed Brown #undef __FUNCT__ 96061692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian" 961e27a552bSJed Brown /*@C 96261692a83SJed Brown TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 963e27a552bSJed Brown 964e27a552bSJed Brown Logically collective 965e27a552bSJed Brown 966e27a552bSJed Brown Input Parameter: 967e27a552bSJed Brown + ts - timestepping context 96861692a83SJed Brown - flg - PETSC_TRUE to recompute the Jacobian at each stage 969e27a552bSJed Brown 970e27a552bSJed Brown Level: intermediate 971e27a552bSJed Brown 972e27a552bSJed Brown .seealso: TSRosWGetType() 973e27a552bSJed Brown @*/ 97461692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg) 975e27a552bSJed Brown { 976e27a552bSJed Brown PetscErrorCode ierr; 977e27a552bSJed Brown 978e27a552bSJed Brown PetscFunctionBegin; 979e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 98061692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 981e27a552bSJed Brown PetscFunctionReturn(0); 982e27a552bSJed Brown } 983e27a552bSJed Brown 984e27a552bSJed Brown EXTERN_C_BEGIN 985e27a552bSJed Brown #undef __FUNCT__ 986e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW" 98761692a83SJed Brown PetscErrorCode TSRosWGetType_RosW(TS ts,const TSRosWType *rostype) 988e27a552bSJed Brown { 98961692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 990e27a552bSJed Brown PetscErrorCode ierr; 991e27a552bSJed Brown 992e27a552bSJed Brown PetscFunctionBegin; 99361692a83SJed Brown if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);} 99461692a83SJed Brown *rostype = ros->tableau->name; 995e27a552bSJed Brown PetscFunctionReturn(0); 996e27a552bSJed Brown } 997e27a552bSJed Brown #undef __FUNCT__ 998e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW" 99961692a83SJed Brown PetscErrorCode TSRosWSetType_RosW(TS ts,const TSRosWType rostype) 1000e27a552bSJed Brown { 100161692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1002e27a552bSJed Brown PetscErrorCode ierr; 1003e27a552bSJed Brown PetscBool match; 100461692a83SJed Brown RosWTableauLink link; 1005e27a552bSJed Brown 1006e27a552bSJed Brown PetscFunctionBegin; 100761692a83SJed Brown if (ros->tableau) { 100861692a83SJed Brown ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr); 1009e27a552bSJed Brown if (match) PetscFunctionReturn(0); 1010e27a552bSJed Brown } 101161692a83SJed Brown for (link = RosWTableauList; link; link=link->next) { 101261692a83SJed Brown ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr); 1013e27a552bSJed Brown if (match) { 1014e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 101561692a83SJed Brown ros->tableau = &link->tab; 1016e27a552bSJed Brown PetscFunctionReturn(0); 1017e27a552bSJed Brown } 1018e27a552bSJed Brown } 101961692a83SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype); 1020e27a552bSJed Brown PetscFunctionReturn(0); 1021e27a552bSJed Brown } 102261692a83SJed Brown 1023e27a552bSJed Brown #undef __FUNCT__ 102461692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW" 102561692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg) 1026e27a552bSJed Brown { 102761692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1028e27a552bSJed Brown 1029e27a552bSJed Brown PetscFunctionBegin; 103061692a83SJed Brown ros->recompute_jacobian = flg; 1031e27a552bSJed Brown PetscFunctionReturn(0); 1032e27a552bSJed Brown } 1033e27a552bSJed Brown EXTERN_C_END 1034e27a552bSJed Brown 1035e27a552bSJed Brown /* ------------------------------------------------------------ */ 1036e27a552bSJed Brown /*MC 1037e27a552bSJed Brown TSRosW - ODE solver using Rosenbrock-W schemes 1038e27a552bSJed Brown 1039e27a552bSJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1040e27a552bSJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1041e27a552bSJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1042e27a552bSJed Brown 1043e27a552bSJed Brown Notes: 104461692a83SJed Brown This method currently only works with autonomous ODE and DAE. 104561692a83SJed Brown 104661692a83SJed Brown Developer notes: 104761692a83SJed Brown Rosenbrock-W methods are typically specified for autonomous ODE 104861692a83SJed Brown 104961692a83SJed Brown $ xdot = f(x) 105061692a83SJed Brown 105161692a83SJed Brown by the stage equations 105261692a83SJed Brown 105361692a83SJed Brown $ k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 105461692a83SJed Brown 105561692a83SJed Brown and step completion formula 105661692a83SJed Brown 105761692a83SJed Brown $ x_1 = x_0 + sum_j b_j k_j 105861692a83SJed Brown 105961692a83SJed Brown with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x) 106061692a83SJed Brown and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 106161692a83SJed Brown we define new variables for the stage equations 106261692a83SJed Brown 106361692a83SJed Brown $ y_i = gamma_ij k_j 106461692a83SJed Brown 106561692a83SJed Brown The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 106661692a83SJed Brown 106761692a83SJed Brown $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i} 106861692a83SJed Brown 106961692a83SJed Brown to rewrite the method as 107061692a83SJed Brown 107161692a83SJed 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 107261692a83SJed Brown $ x_1 = x_0 + sum_j bt_j y_j 107361692a83SJed Brown 107461692a83SJed Brown where we have introduced the mass matrix M. Continue by defining 107561692a83SJed Brown 107661692a83SJed Brown $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 107761692a83SJed Brown 107861692a83SJed Brown or, more compactly in tensor notation 107961692a83SJed Brown 108061692a83SJed Brown $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 108161692a83SJed Brown 108261692a83SJed Brown Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 108361692a83SJed Brown stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 108461692a83SJed Brown equation 108561692a83SJed Brown 108661692a83SJed Brown $ g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 108761692a83SJed Brown 108861692a83SJed Brown with initial guess y_i = 0. 1089e27a552bSJed Brown 1090e27a552bSJed Brown Level: beginner 1091e27a552bSJed Brown 1092e27a552bSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSRosWRegister() 1093e27a552bSJed Brown 1094e27a552bSJed Brown M*/ 1095e27a552bSJed Brown EXTERN_C_BEGIN 1096e27a552bSJed Brown #undef __FUNCT__ 1097e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW" 1098e27a552bSJed Brown PetscErrorCode TSCreate_RosW(TS ts) 1099e27a552bSJed Brown { 110061692a83SJed Brown TS_RosW *ros; 1101e27a552bSJed Brown PetscErrorCode ierr; 1102e27a552bSJed Brown 1103e27a552bSJed Brown PetscFunctionBegin; 1104e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1105e27a552bSJed Brown ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1106e27a552bSJed Brown #endif 1107e27a552bSJed Brown 1108e27a552bSJed Brown ts->ops->reset = TSReset_RosW; 1109e27a552bSJed Brown ts->ops->destroy = TSDestroy_RosW; 1110e27a552bSJed Brown ts->ops->view = TSView_RosW; 1111e27a552bSJed Brown ts->ops->setup = TSSetUp_RosW; 1112e27a552bSJed Brown ts->ops->step = TSStep_RosW; 1113e27a552bSJed Brown ts->ops->interpolate = TSInterpolate_RosW; 11141c3436cfSJed Brown ts->ops->evaluatestep = TSEvaluateStep_RosW; 1115e27a552bSJed Brown ts->ops->setfromoptions = TSSetFromOptions_RosW; 1116e27a552bSJed Brown ts->ops->snesfunction = SNESTSFormFunction_RosW; 1117e27a552bSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1118e27a552bSJed Brown 111961692a83SJed Brown ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr); 112061692a83SJed Brown ts->data = (void*)ros; 1121e27a552bSJed Brown 1122e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr); 1123e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr); 112461692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr); 1125e27a552bSJed Brown PetscFunctionReturn(0); 1126e27a552bSJed Brown } 1127e27a552bSJed Brown EXTERN_C_END 1128