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 */ 2961692a83SJed Brown PetscReal *b; /* Step completion table */ 30fe7e6d57SJed Brown PetscReal *bembed; /* Step completion table for embedded method of order one less */ 3161692a83SJed Brown PetscReal *ASum; /* Row sum of A */ 3261692a83SJed Brown PetscReal *GammaSum; /* Row sum of Gamma, only needed for non-autonomous systems */ 3361692a83SJed Brown PetscReal *At; /* Propagation table in transformed variables */ 3461692a83SJed Brown PetscReal *bt; /* Step completion table in transformed variables */ 35fe7e6d57SJed Brown PetscReal *bembedt; /* Step completion table of order one less in transformed variables */ 3661692a83SJed Brown PetscReal *GammaInv; /* Inverse of Gamma, used for transformed variables */ 378d59e960SJed Brown PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 38e27a552bSJed Brown }; 3961692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink; 4061692a83SJed Brown struct _RosWTableauLink { 4161692a83SJed Brown struct _RosWTableau tab; 4261692a83SJed Brown RosWTableauLink next; 43e27a552bSJed Brown }; 4461692a83SJed Brown static RosWTableauLink RosWTableauList; 45e27a552bSJed Brown 46e27a552bSJed Brown typedef struct { 4761692a83SJed Brown RosWTableau tableau; 4861692a83SJed Brown Vec *Y; /* States computed during the step, used to complete the step */ 49e27a552bSJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 5061692a83SJed Brown Vec Ystage; /* Work vector for the state value at each stage */ 5161692a83SJed Brown Vec Zdot; /* Ydot = Zdot + shift*Y */ 5261692a83SJed Brown Vec Zstage; /* Y = Zstage + Y */ 531c3436cfSJed Brown PetscScalar *work; /* Scalar work space of length number of stages, used to prepare VecMAXPY() */ 54e27a552bSJed Brown PetscReal shift; 55e27a552bSJed Brown PetscReal stage_time; 56c17803e7SJed Brown PetscReal stage_explicit; /* Flag indicates that the current stage is explicit */ 5761692a83SJed Brown PetscBool recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */ 581c3436cfSJed Brown PetscBool step_taken; /* ts->vec_sol has been advanced to the end of the current time step */ 59e27a552bSJed Brown } TS_RosW; 60e27a552bSJed Brown 61fe7e6d57SJed Brown /*MC 62fe7e6d57SJed Brown TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme. 63fe7e6d57SJed Brown 64fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P. 65fe7e6d57SJed Brown 66fe7e6d57SJed Brown Level: intermediate 67fe7e6d57SJed Brown 68fe7e6d57SJed Brown .seealso: TSROSW 69fe7e6d57SJed Brown M*/ 70fe7e6d57SJed Brown 71fe7e6d57SJed Brown /*MC 72fe7e6d57SJed Brown TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme. 73fe7e6d57SJed Brown 74fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M. 75fe7e6d57SJed Brown 76fe7e6d57SJed Brown Level: intermediate 77fe7e6d57SJed Brown 78fe7e6d57SJed Brown .seealso: TSROSW 79fe7e6d57SJed Brown M*/ 80fe7e6d57SJed Brown 81fe7e6d57SJed Brown /*MC 82fe7e6d57SJed Brown TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1. 83fe7e6d57SJed Brown 84fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 85fe7e6d57SJed Brown 86fe7e6d57SJed 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. 87fe7e6d57SJed Brown 88fe7e6d57SJed Brown References: 89fe7e6d57SJed Brown Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 90fe7e6d57SJed Brown 91fe7e6d57SJed Brown Level: intermediate 92fe7e6d57SJed Brown 93fe7e6d57SJed Brown .seealso: TSROSW 94fe7e6d57SJed Brown M*/ 95fe7e6d57SJed Brown 96fe7e6d57SJed Brown /*MC 97fe7e6d57SJed Brown TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1. 98fe7e6d57SJed Brown 99fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 100fe7e6d57SJed Brown 101fe7e6d57SJed 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. 102fe7e6d57SJed Brown 103fe7e6d57SJed Brown References: 104fe7e6d57SJed Brown Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 105fe7e6d57SJed Brown 106fe7e6d57SJed Brown Level: intermediate 107fe7e6d57SJed Brown 108fe7e6d57SJed Brown .seealso: TSROSW 109fe7e6d57SJed Brown M*/ 110fe7e6d57SJed Brown 111ef3c5b88SJed Brown /*MC 112ef3c5b88SJed Brown TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme 113ef3c5b88SJed Brown 114ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 115ef3c5b88SJed Brown 116ef3c5b88SJed Brown Both the third order and embedded second order methods are stiffly accurate and L-stable. 117ef3c5b88SJed Brown 118ef3c5b88SJed Brown References: 119ef3c5b88SJed Brown Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 120ef3c5b88SJed Brown 121ef3c5b88SJed Brown Level: intermediate 122ef3c5b88SJed Brown 123ef3c5b88SJed Brown .seealso: TSROSW, TSROSWSANDU3 124ef3c5b88SJed Brown M*/ 125ef3c5b88SJed Brown 126ef3c5b88SJed Brown /*MC 127ef3c5b88SJed Brown TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme 128ef3c5b88SJed Brown 129ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 130ef3c5b88SJed Brown 131ef3c5b88SJed Brown The third order method is L-stable, but not stiffly accurate. 132ef3c5b88SJed Brown The second order embedded method is strongly A-stable with R(infty) = 0.5. 133ef3c5b88SJed Brown The internal stages are L-stable. 134ef3c5b88SJed Brown This method is called ROS3 in the paper. 135ef3c5b88SJed Brown 136ef3c5b88SJed Brown References: 137ef3c5b88SJed Brown Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 138ef3c5b88SJed Brown 139ef3c5b88SJed Brown Level: intermediate 140ef3c5b88SJed Brown 141ef3c5b88SJed Brown .seealso: TSROSW, TSROSWRODAS3 142ef3c5b88SJed Brown M*/ 143ef3c5b88SJed Brown 144961f28d0SJed Brown /*MC 145961f28d0SJed Brown TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages 146961f28d0SJed Brown 147961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 148961f28d0SJed Brown 149961f28d0SJed Brown A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3) 150961f28d0SJed Brown 151961f28d0SJed Brown References: 152961f28d0SJed Brown Emil Constantinescu 153961f28d0SJed Brown 154961f28d0SJed Brown Level: intermediate 155961f28d0SJed Brown 156961f28d0SJed Brown .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P3S2C, SSP 157961f28d0SJed Brown M*/ 158961f28d0SJed Brown 159961f28d0SJed Brown /*MC 160961f28d0SJed Brown TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 161961f28d0SJed Brown 162961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 163961f28d0SJed Brown 164961f28d0SJed Brown L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 165961f28d0SJed Brown 166961f28d0SJed Brown References: 167961f28d0SJed Brown Emil Constantinescu 168961f28d0SJed Brown 169961f28d0SJed Brown Level: intermediate 170961f28d0SJed Brown 171961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P3S2C, TSSSP 172961f28d0SJed Brown M*/ 173961f28d0SJed Brown 174961f28d0SJed Brown /*MC 175961f28d0SJed Brown TSROSWLLSSP3P3S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 176961f28d0SJed Brown 177961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 178961f28d0SJed Brown 179961f28d0SJed Brown L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 180961f28d0SJed Brown 181961f28d0SJed Brown References: 182961f28d0SJed Brown Emil Constantinescu 183961f28d0SJed Brown 184961f28d0SJed Brown Level: intermediate 185961f28d0SJed Brown 186961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP 187961f28d0SJed Brown M*/ 188961f28d0SJed Brown 189e27a552bSJed Brown #undef __FUNCT__ 190e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll" 191e27a552bSJed Brown /*@C 192e27a552bSJed Brown TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW 193e27a552bSJed Brown 194e27a552bSJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 195e27a552bSJed Brown 196e27a552bSJed Brown Level: advanced 197e27a552bSJed Brown 198e27a552bSJed Brown .keywords: TS, TSRosW, register, all 199e27a552bSJed Brown 200e27a552bSJed Brown .seealso: TSRosWRegisterDestroy() 201e27a552bSJed Brown @*/ 202e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void) 203e27a552bSJed Brown { 204e27a552bSJed Brown PetscErrorCode ierr; 205e27a552bSJed Brown 206e27a552bSJed Brown PetscFunctionBegin; 207e27a552bSJed Brown if (TSRosWRegisterAllCalled) PetscFunctionReturn(0); 208e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_TRUE; 209e27a552bSJed Brown 210e27a552bSJed Brown { 21161692a83SJed Brown const PetscReal g = 1. + 1./PetscSqrtReal(2.0); 212e27a552bSJed Brown const PetscReal 21361692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 21461692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2151c3436cfSJed Brown b[2] = {0.5,0.5}, 2161c3436cfSJed Brown b1[2] = {1.0,0.0}; 2171c3436cfSJed Brown ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr); 218e27a552bSJed Brown } 219e27a552bSJed Brown { 22061692a83SJed Brown const PetscReal g = 1. - 1./PetscSqrtReal(2.0); 221e27a552bSJed Brown const PetscReal 22261692a83SJed Brown A[2][2] = {{0,0}, {1.,0}}, 22361692a83SJed Brown Gamma[2][2] = {{g,0}, {-2.*g,g}}, 2241c3436cfSJed Brown b[2] = {0.5,0.5}, 2251c3436cfSJed Brown b1[2] = {1.0,0.0}; 2261c3436cfSJed Brown ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr); 227fe7e6d57SJed Brown } 228fe7e6d57SJed Brown { 229fe7e6d57SJed Brown const PetscReal g = 7.8867513459481287e-01; 230fe7e6d57SJed Brown const PetscReal 231fe7e6d57SJed Brown A[3][3] = {{0,0,0}, 232fe7e6d57SJed Brown {1.5773502691896257e+00,0,0}, 233fe7e6d57SJed Brown {0.5,0,0}}, 234fe7e6d57SJed Brown Gamma[3][3] = {{g,0,0}, 235fe7e6d57SJed Brown {-1.5773502691896257e+00,g,0}, 236fe7e6d57SJed Brown {-6.7075317547305480e-01,1.7075317547305482e-01,g}}, 237fe7e6d57SJed Brown b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01}, 238fe7e6d57SJed Brown b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01}; 239fe7e6d57SJed Brown ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 240fe7e6d57SJed Brown } 241fe7e6d57SJed Brown { 242fe7e6d57SJed Brown const PetscReal g = 4.3586652150845900e-01; 243fe7e6d57SJed Brown const PetscReal 244fe7e6d57SJed Brown A[4][4] = {{0,0,0,0}, 245fe7e6d57SJed Brown {8.7173304301691801e-01,0,0,0}, 246fe7e6d57SJed Brown {8.4457060015369423e-01,-1.1299064236484185e-01,0,0}, 247fe7e6d57SJed Brown {0,0,1.,0}}, 248fe7e6d57SJed Brown Gamma[4][4] = {{g,0,0,0}, 249fe7e6d57SJed Brown {-8.7173304301691801e-01,g,0,0}, 250fe7e6d57SJed Brown {-9.0338057013044082e-01,5.4180672388095326e-02,g,0}, 251fe7e6d57SJed Brown {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}}, 252fe7e6d57SJed Brown b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01}, 253fe7e6d57SJed Brown b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01}; 254fe7e6d57SJed Brown ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 255e27a552bSJed Brown } 256ef3c5b88SJed Brown { 257ef3c5b88SJed Brown const PetscReal g = 0.5; 258ef3c5b88SJed Brown const PetscReal 259ef3c5b88SJed Brown A[4][4] = {{0,0,0,0}, 260ef3c5b88SJed Brown {0,0,0,0}, 261ef3c5b88SJed Brown {1.,0,0,0}, 262ef3c5b88SJed Brown {0.75,-0.25,0.5,0}}, 263ef3c5b88SJed Brown Gamma[4][4] = {{g,0,0,0}, 264ef3c5b88SJed Brown {1.,g,0,0}, 265ef3c5b88SJed Brown {-0.25,-0.25,g,0}, 266ef3c5b88SJed Brown {1./12,1./12,-2./3,g}}, 267ef3c5b88SJed Brown b[4] = {5./6,-1./6,-1./6,0.5}, 268ef3c5b88SJed Brown b2[4] = {0.75,-0.25,0.5,0}; 269ef3c5b88SJed Brown ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 270ef3c5b88SJed Brown } 271ef3c5b88SJed Brown { 272ef3c5b88SJed Brown const PetscReal g = 0.43586652150845899941601945119356; 273ef3c5b88SJed Brown const PetscReal 274ef3c5b88SJed Brown A[3][3] = {{0,0,0}, 275ef3c5b88SJed Brown {g,0,0}, 276ef3c5b88SJed Brown {g,0,0}}, 277ef3c5b88SJed Brown Gamma[3][3] = {{g,0,0}, 278ef3c5b88SJed Brown {-0.19294655696029095575009695436041,g,0}, 279ef3c5b88SJed Brown {0,1.74927148125794685173529749738960,g}}, 280ef3c5b88SJed Brown b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829}, 281ef3c5b88SJed Brown b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619}; 282ef3c5b88SJed Brown ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 283ef3c5b88SJed Brown } 284b1c69cc3SEmil Constantinescu { 285b1c69cc3SEmil Constantinescu const PetscReal g = (3.0+sqrt(3.0))/6.0; 286b1c69cc3SEmil Constantinescu const PetscReal 287b1c69cc3SEmil Constantinescu A[3][3] = {{0,0,0}, 288b1c69cc3SEmil Constantinescu {1,0,0}, 289b1c69cc3SEmil Constantinescu {0.25,0.25,0}}, 290b1c69cc3SEmil Constantinescu Gamma[3][3] = {{0,0,0}, 291b1c69cc3SEmil Constantinescu {(-3.0-sqrt(3.0))/6.0,g,0}, 292b1c69cc3SEmil Constantinescu {(-3.0-sqrt(3.0))/24.0,(-3.0-sqrt(3.0))/8.0,g}}, 293b1c69cc3SEmil Constantinescu b[3] = {1./6.,1./6.,2./3.}, 294b1c69cc3SEmil Constantinescu b2[3] = {1./4.,1./4.,1./2.}; 295b1c69cc3SEmil Constantinescu ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 296b1c69cc3SEmil Constantinescu } 297b1c69cc3SEmil Constantinescu 298b1c69cc3SEmil Constantinescu { 299b1c69cc3SEmil Constantinescu const PetscReal 300b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 301b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 302b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 303b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 304b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 305b1c69cc3SEmil Constantinescu {0.0,1./4.,0,0}, 306b1c69cc3SEmil Constantinescu {-2.,-2./3.,2./3.,0}, 307b1c69cc3SEmil Constantinescu {1./2.,5./36.,-2./9,0}}, 308b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 309b1c69cc3SEmil Constantinescu b2[4] = {1./8.,3./4.,1./8.,0}; 310b1c69cc3SEmil Constantinescu ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 311b1c69cc3SEmil Constantinescu } 312b1c69cc3SEmil Constantinescu 313b1c69cc3SEmil Constantinescu { 314b1c69cc3SEmil Constantinescu const PetscReal 315b1c69cc3SEmil Constantinescu A[4][4] = {{0,0,0,0}, 316b1c69cc3SEmil Constantinescu {1./2.,0,0,0}, 317b1c69cc3SEmil Constantinescu {1./2.,1./2.,0,0}, 318b1c69cc3SEmil Constantinescu {1./6.,1./6.,1./6.,0}}, 319b1c69cc3SEmil Constantinescu Gamma[4][4] = {{1./2.,0,0,0}, 320b1c69cc3SEmil Constantinescu {0.0,3./4.,0,0}, 321b1c69cc3SEmil Constantinescu {-2./3.,-23./9.,2./9.,0}, 322b1c69cc3SEmil Constantinescu {1./18.,65./108.,-2./27,0}}, 323b1c69cc3SEmil Constantinescu b[4] = {1./6.,1./6.,1./6.,1./2.}, 324b1c69cc3SEmil Constantinescu b2[4] = {3./16.,10./16.,3./16.,0}; 325b1c69cc3SEmil Constantinescu ierr = TSRosWRegister(TSROSWLLSSP3P3S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 326b1c69cc3SEmil Constantinescu } 327753f8adbSEmil Constantinescu 328753f8adbSEmil Constantinescu { 329753f8adbSEmil Constantinescu PetscReal A[4][4],Gamma[4][4],b[4],b2[4]; 330753f8adbSEmil Constantinescu 331753f8adbSEmil Constantinescu Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816; 33205e8e825SJed Brown Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0; 333753f8adbSEmil Constantinescu Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476; 334753f8adbSEmil Constantinescu Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816; 33505e8e825SJed Brown Gamma[1][2]=0; Gamma[1][3]=0; 336753f8adbSEmil Constantinescu Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903; 337753f8adbSEmil Constantinescu Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131; 338753f8adbSEmil Constantinescu Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816; 33905e8e825SJed Brown Gamma[2][3]=0; 340753f8adbSEmil Constantinescu Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783; 341753f8adbSEmil Constantinescu Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984; 342753f8adbSEmil Constantinescu Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198; 343753f8adbSEmil Constantinescu Gamma[3][3]=0; 344753f8adbSEmil Constantinescu 34505e8e825SJed Brown A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0; 346753f8adbSEmil Constantinescu A[1][0]=0.8717330430169179988320388950590125027645343373957631; 34705e8e825SJed Brown A[1][1]=0; A[1][2]=0; A[1][3]=0; 348753f8adbSEmil Constantinescu A[2][0]=0.5275890119763004115618079766722914408876108660811028; 349753f8adbSEmil Constantinescu A[2][1]=0.07241098802369958843819203208518599088698057726988732; 35005e8e825SJed Brown A[2][2]=0; A[2][3]=0; 351753f8adbSEmil Constantinescu A[3][0]=0.3990960076760701320627260685975778145384666450351314; 352753f8adbSEmil Constantinescu A[3][1]=-0.4375576546135194437228463747348862825846903771419953; 353753f8adbSEmil Constantinescu A[3][2]=1.038461646937449311660120300601880176655352737312713; 35405e8e825SJed Brown A[3][3]=0; 355753f8adbSEmil Constantinescu 356753f8adbSEmil Constantinescu b[0]=0.1876410243467238251612921333138006734899663569186926; 357753f8adbSEmil Constantinescu b[1]=-0.5952974735769549480478230473706443582188442040780541; 358753f8adbSEmil Constantinescu b[2]=0.9717899277217721234705114616271378792182450260943198; 359753f8adbSEmil Constantinescu b[3]=0.4358665215084589994160194475295062513822671686978816; 360753f8adbSEmil Constantinescu 361753f8adbSEmil Constantinescu b2[0]=0.2147402862233891404862383521089097657790734483804460; 362753f8adbSEmil Constantinescu b2[1]=-0.4851622638849390928209050538171743017757490232519684; 363753f8adbSEmil Constantinescu b2[2]=0.8687250025203875511662123688667549217531982787600080; 364753f8adbSEmil Constantinescu b2[3]=0.4016969751411624011684543450940068201770721128357014; 365753f8adbSEmil Constantinescu 366753f8adbSEmil Constantinescu ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 367753f8adbSEmil Constantinescu } 368753f8adbSEmil Constantinescu 369e27a552bSJed Brown PetscFunctionReturn(0); 370e27a552bSJed Brown } 371e27a552bSJed Brown 372e27a552bSJed Brown #undef __FUNCT__ 373e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy" 374e27a552bSJed Brown /*@C 375e27a552bSJed Brown TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 376e27a552bSJed Brown 377e27a552bSJed Brown Not Collective 378e27a552bSJed Brown 379e27a552bSJed Brown Level: advanced 380e27a552bSJed Brown 381e27a552bSJed Brown .keywords: TSRosW, register, destroy 382e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 383e27a552bSJed Brown @*/ 384e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void) 385e27a552bSJed Brown { 386e27a552bSJed Brown PetscErrorCode ierr; 38761692a83SJed Brown RosWTableauLink link; 388e27a552bSJed Brown 389e27a552bSJed Brown PetscFunctionBegin; 39061692a83SJed Brown while ((link = RosWTableauList)) { 39161692a83SJed Brown RosWTableau t = &link->tab; 39261692a83SJed Brown RosWTableauList = link->next; 39361692a83SJed Brown ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr); 394c17803e7SJed Brown ierr = PetscFree4(t->At,t->bt,t->GammaInv,t->GammaZeroDiag);CHKERRQ(ierr); 395fe7e6d57SJed Brown ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr); 396e27a552bSJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 397e27a552bSJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 398e27a552bSJed Brown } 399e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_FALSE; 400e27a552bSJed Brown PetscFunctionReturn(0); 401e27a552bSJed Brown } 402e27a552bSJed Brown 403e27a552bSJed Brown #undef __FUNCT__ 404e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage" 405e27a552bSJed Brown /*@C 406e27a552bSJed Brown TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called 407e27a552bSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW() 408e27a552bSJed Brown when using static libraries. 409e27a552bSJed Brown 410e27a552bSJed Brown Input Parameter: 411e27a552bSJed Brown path - The dynamic library path, or PETSC_NULL 412e27a552bSJed Brown 413e27a552bSJed Brown Level: developer 414e27a552bSJed Brown 415e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package 416e27a552bSJed Brown .seealso: PetscInitialize() 417e27a552bSJed Brown @*/ 418e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[]) 419e27a552bSJed Brown { 420e27a552bSJed Brown PetscErrorCode ierr; 421e27a552bSJed Brown 422e27a552bSJed Brown PetscFunctionBegin; 423e27a552bSJed Brown if (TSRosWPackageInitialized) PetscFunctionReturn(0); 424e27a552bSJed Brown TSRosWPackageInitialized = PETSC_TRUE; 425e27a552bSJed Brown ierr = TSRosWRegisterAll();CHKERRQ(ierr); 426e27a552bSJed Brown ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr); 427e27a552bSJed Brown PetscFunctionReturn(0); 428e27a552bSJed Brown } 429e27a552bSJed Brown 430e27a552bSJed Brown #undef __FUNCT__ 431e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage" 432e27a552bSJed Brown /*@C 433e27a552bSJed Brown TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is 434e27a552bSJed Brown called from PetscFinalize(). 435e27a552bSJed Brown 436e27a552bSJed Brown Level: developer 437e27a552bSJed Brown 438e27a552bSJed Brown .keywords: Petsc, destroy, package 439e27a552bSJed Brown .seealso: PetscFinalize() 440e27a552bSJed Brown @*/ 441e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void) 442e27a552bSJed Brown { 443e27a552bSJed Brown PetscErrorCode ierr; 444e27a552bSJed Brown 445e27a552bSJed Brown PetscFunctionBegin; 446e27a552bSJed Brown TSRosWPackageInitialized = PETSC_FALSE; 447e27a552bSJed Brown ierr = TSRosWRegisterDestroy();CHKERRQ(ierr); 448e27a552bSJed Brown PetscFunctionReturn(0); 449e27a552bSJed Brown } 450e27a552bSJed Brown 451e27a552bSJed Brown #undef __FUNCT__ 452e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister" 453e27a552bSJed Brown /*@C 45461692a83SJed Brown TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 455e27a552bSJed Brown 456e27a552bSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 457e27a552bSJed Brown 458e27a552bSJed Brown Input Parameters: 459e27a552bSJed Brown + name - identifier for method 460e27a552bSJed Brown . order - approximation order of method 461e27a552bSJed Brown . s - number of stages, this is the dimension of the matrices below 46261692a83SJed Brown . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 46361692a83SJed Brown . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 464fe7e6d57SJed Brown . b - Step completion table (dimension s) 465fe7e6d57SJed Brown - bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available) 466e27a552bSJed Brown 467e27a552bSJed Brown Notes: 46861692a83SJed Brown Several Rosenbrock W methods are provided, this function is only needed to create new methods. 469e27a552bSJed Brown 470e27a552bSJed Brown Level: advanced 471e27a552bSJed Brown 472e27a552bSJed Brown .keywords: TS, register 473e27a552bSJed Brown 474e27a552bSJed Brown .seealso: TSRosW 475e27a552bSJed Brown @*/ 476e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s, 477fe7e6d57SJed Brown const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[]) 478e27a552bSJed Brown { 479e27a552bSJed Brown PetscErrorCode ierr; 48061692a83SJed Brown RosWTableauLink link; 48161692a83SJed Brown RosWTableau t; 48261692a83SJed Brown PetscInt i,j,k; 48361692a83SJed Brown PetscScalar *GammaInv; 484e27a552bSJed Brown 485e27a552bSJed Brown PetscFunctionBegin; 486fe7e6d57SJed Brown PetscValidCharPointer(name,1); 487fe7e6d57SJed Brown PetscValidPointer(A,4); 488fe7e6d57SJed Brown PetscValidPointer(Gamma,5); 489fe7e6d57SJed Brown PetscValidPointer(b,6); 490fe7e6d57SJed Brown if (bembed) PetscValidPointer(bembed,7); 491fe7e6d57SJed Brown 492e27a552bSJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 493e27a552bSJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 494e27a552bSJed Brown t = &link->tab; 495e27a552bSJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 496e27a552bSJed Brown t->order = order; 497e27a552bSJed Brown t->s = s; 49861692a83SJed 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); 499c17803e7SJed Brown ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag);CHKERRQ(ierr); 500e27a552bSJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 50161692a83SJed Brown ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 50261692a83SJed Brown ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 503fe7e6d57SJed Brown if (bembed) { 504fe7e6d57SJed Brown ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr); 505fe7e6d57SJed Brown ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 506fe7e6d57SJed Brown } 50761692a83SJed Brown for (i=0; i<s; i++) { 50861692a83SJed Brown t->ASum[i] = 0; 50961692a83SJed Brown t->GammaSum[i] = 0; 51061692a83SJed Brown for (j=0; j<s; j++) { 51161692a83SJed Brown t->ASum[i] += A[i*s+j]; 512fe7e6d57SJed Brown t->GammaSum[i] += Gamma[i*s+j]; 51361692a83SJed Brown } 51461692a83SJed Brown } 51561692a83SJed Brown ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */ 51661692a83SJed Brown for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i]; 517fd96d5b0SEmil Constantinescu for (i=0; i<s; i++) { 518fd96d5b0SEmil Constantinescu if (Gamma[i*s+i] == 0.0) { 519fd96d5b0SEmil Constantinescu GammaInv[i*s+i] = 1.0; 520c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_TRUE; 521fd96d5b0SEmil Constantinescu } else { 522c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_FALSE; 523fd96d5b0SEmil Constantinescu } 524fd96d5b0SEmil Constantinescu } 525fd96d5b0SEmil Constantinescu 52661692a83SJed Brown switch (s) { 52761692a83SJed Brown case 1: GammaInv[0] = 1./GammaInv[0]; break; 52861692a83SJed Brown case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break; 52961692a83SJed Brown case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break; 53061692a83SJed Brown case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break; 53161692a83SJed Brown case 5: { 53261692a83SJed Brown PetscInt ipvt5[5]; 53361692a83SJed Brown MatScalar work5[5*5]; 53461692a83SJed Brown ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break; 53561692a83SJed Brown } 53661692a83SJed Brown case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break; 53761692a83SJed Brown case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break; 53861692a83SJed Brown default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s); 53961692a83SJed Brown } 54061692a83SJed Brown for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 54161692a83SJed Brown ierr = PetscFree(GammaInv);CHKERRQ(ierr); 54261692a83SJed Brown for (i=0; i<s; i++) { 54361692a83SJed Brown for (j=0; j<s; j++) { 54461692a83SJed Brown t->At[i*s+j] = 0; 54561692a83SJed Brown for (k=0; k<s; k++) { 54661692a83SJed Brown t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j]; 54761692a83SJed Brown } 54861692a83SJed Brown } 54961692a83SJed Brown t->bt[i] = 0; 55061692a83SJed Brown for (j=0; j<s; j++) { 55161692a83SJed Brown t->bt[i] += t->b[j] * t->GammaInv[j*s+i]; 55261692a83SJed Brown } 553fe7e6d57SJed Brown if (bembed) { 554fe7e6d57SJed Brown t->bembedt[i] = 0; 555fe7e6d57SJed Brown for (j=0; j<s; j++) { 556fe7e6d57SJed Brown t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i]; 557fe7e6d57SJed Brown } 558fe7e6d57SJed Brown } 55961692a83SJed Brown } 5608d59e960SJed Brown t->ccfl = 1.0; /* Fix this */ 5618d59e960SJed Brown 56261692a83SJed Brown link->next = RosWTableauList; 56361692a83SJed Brown RosWTableauList = link; 564e27a552bSJed Brown PetscFunctionReturn(0); 565e27a552bSJed Brown } 566e27a552bSJed Brown 567e27a552bSJed Brown #undef __FUNCT__ 5681c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW" 5691c3436cfSJed Brown /* 5701c3436cfSJed Brown The step completion formula is 5711c3436cfSJed Brown 5721c3436cfSJed Brown x1 = x0 + b^T Y 5731c3436cfSJed Brown 5741c3436cfSJed Brown where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 5751c3436cfSJed Brown updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 5761c3436cfSJed Brown 5771c3436cfSJed Brown x1e = x0 + be^T Y 5781c3436cfSJed Brown = x1 - b^T Y + be^T Y 5791c3436cfSJed Brown = x1 + (be - b)^T Y 5801c3436cfSJed Brown 5811c3436cfSJed Brown so we can evaluate the method of different order even after the step has been optimistically completed. 5821c3436cfSJed Brown */ 5831c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done) 5841c3436cfSJed Brown { 5851c3436cfSJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 5861c3436cfSJed Brown RosWTableau tab = ros->tableau; 5871c3436cfSJed Brown PetscScalar *w = ros->work; 5881c3436cfSJed Brown PetscInt i; 5891c3436cfSJed Brown PetscErrorCode ierr; 5901c3436cfSJed Brown 5911c3436cfSJed Brown PetscFunctionBegin; 5921c3436cfSJed Brown if (order == tab->order) { 5931c3436cfSJed Brown if (ros->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 5941c3436cfSJed Brown else { 5951c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 596*de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bt[i]; 597*de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 5981c3436cfSJed Brown } 5991c3436cfSJed Brown if (done) *done = PETSC_TRUE; 6001c3436cfSJed Brown PetscFunctionReturn(0); 6011c3436cfSJed Brown } else if (order == tab->order-1) { 6021c3436cfSJed Brown if (!tab->bembedt) goto unavailable; 6031c3436cfSJed Brown if (ros->step_taken) { 6041c3436cfSJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 6051c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 6061c3436cfSJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 6071c3436cfSJed Brown } else { 6081c3436cfSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 609*de19f811SJed Brown for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i]; 610*de19f811SJed Brown ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 6111c3436cfSJed Brown } 6121c3436cfSJed Brown if (done) *done = PETSC_TRUE; 6131c3436cfSJed Brown PetscFunctionReturn(0); 6141c3436cfSJed Brown } 6151c3436cfSJed Brown unavailable: 6161c3436cfSJed Brown if (done) *done = PETSC_FALSE; 6171c3436cfSJed 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); 6181c3436cfSJed Brown PetscFunctionReturn(0); 6191c3436cfSJed Brown } 6201c3436cfSJed Brown 6211c3436cfSJed Brown #undef __FUNCT__ 622e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW" 623e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts) 624e27a552bSJed Brown { 62561692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 62661692a83SJed Brown RosWTableau tab = ros->tableau; 627e27a552bSJed Brown const PetscInt s = tab->s; 6281c3436cfSJed Brown const PetscReal *At = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv; 629c17803e7SJed Brown const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 63061692a83SJed Brown PetscScalar *w = ros->work; 63161692a83SJed Brown Vec *Y = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage; 632e27a552bSJed Brown SNES snes; 6331c3436cfSJed Brown TSAdapt adapt; 6341c3436cfSJed Brown PetscInt i,j,its,lits,reject,next_scheme; 635cdbf8f93SLisandro Dalcin PetscReal next_time_step; 6361c3436cfSJed Brown PetscBool accept; 637e27a552bSJed Brown PetscErrorCode ierr; 638e27a552bSJed Brown 639e27a552bSJed Brown PetscFunctionBegin; 640e27a552bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 641cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 6421c3436cfSJed Brown accept = PETSC_TRUE; 6431c3436cfSJed Brown ros->step_taken = PETSC_FALSE; 644e27a552bSJed Brown 6451c3436cfSJed Brown for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 6461c3436cfSJed Brown const PetscReal h = ts->time_step; 647e27a552bSJed Brown for (i=0; i<s; i++) { 6481c3436cfSJed Brown ros->stage_time = ts->ptime + h*ASum[i]; 649c17803e7SJed Brown if (GammaZeroDiag[i]) { 650c17803e7SJed Brown ros->stage_explicit = PETSC_TRUE; 651fd96d5b0SEmil Constantinescu ros->shift = 1./h; 652c17803e7SJed Brown } else { 653c17803e7SJed Brown ros->stage_explicit = PETSC_FALSE; 65461692a83SJed Brown ros->shift = 1./(h*Gamma[i*s+i]); 655fd96d5b0SEmil Constantinescu } 65661692a83SJed Brown 65761692a83SJed Brown ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr); 658*de19f811SJed Brown for (j=0; j<i; j++) w[j] = At[i*s+j]; 659*de19f811SJed Brown ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 66061692a83SJed Brown 66161692a83SJed Brown for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j]; 66261692a83SJed Brown ierr = VecZeroEntries(Zdot);CHKERRQ(ierr); 66361692a83SJed Brown ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr); 66461692a83SJed Brown 665e27a552bSJed Brown /* Initial guess taken from last stage */ 66661692a83SJed Brown ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr); 66761692a83SJed Brown 66861692a83SJed Brown if (!ros->recompute_jacobian && !i) { 66961692a83SJed Brown ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 67061692a83SJed Brown } 67161692a83SJed Brown 67261692a83SJed Brown ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr); 673e27a552bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 674e27a552bSJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 675e27a552bSJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 676e27a552bSJed Brown } 6771c3436cfSJed Brown ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 6781c3436cfSJed Brown ros->step_taken = PETSC_TRUE; 679e27a552bSJed Brown 6801c3436cfSJed Brown /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 6811c3436cfSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 6821c3436cfSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 6838d59e960SJed Brown ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 6841c3436cfSJed Brown ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 6851c3436cfSJed Brown if (accept) { 6861c3436cfSJed Brown /* ignore next_scheme for now */ 687e27a552bSJed Brown ts->ptime += ts->time_step; 688cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 689e27a552bSJed Brown ts->steps++; 6901c3436cfSJed Brown break; 6911c3436cfSJed Brown } else { /* Roll back the current step */ 6921c3436cfSJed Brown for (i=0; i<s; i++) w[i] = -tab->bt[i]; 6931c3436cfSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr); 6941c3436cfSJed Brown ts->time_step = next_time_step; 6951c3436cfSJed Brown ros->step_taken = PETSC_FALSE; 6961c3436cfSJed Brown } 6971c3436cfSJed Brown } 6981c3436cfSJed Brown 699e27a552bSJed Brown PetscFunctionReturn(0); 700e27a552bSJed Brown } 701e27a552bSJed Brown 702e27a552bSJed Brown #undef __FUNCT__ 703e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW" 704e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X) 705e27a552bSJed Brown { 70661692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 707e27a552bSJed Brown 708e27a552bSJed Brown PetscFunctionBegin; 70961692a83SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name); 710e27a552bSJed Brown PetscFunctionReturn(0); 711e27a552bSJed Brown } 712e27a552bSJed Brown 713e27a552bSJed Brown /*------------------------------------------------------------*/ 714e27a552bSJed Brown #undef __FUNCT__ 715e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW" 716e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts) 717e27a552bSJed Brown { 71861692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 719e27a552bSJed Brown PetscInt s; 720e27a552bSJed Brown PetscErrorCode ierr; 721e27a552bSJed Brown 722e27a552bSJed Brown PetscFunctionBegin; 72361692a83SJed Brown if (!ros->tableau) PetscFunctionReturn(0); 72461692a83SJed Brown s = ros->tableau->s; 72561692a83SJed Brown ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr); 72661692a83SJed Brown ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr); 72761692a83SJed Brown ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr); 72861692a83SJed Brown ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr); 72961692a83SJed Brown ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr); 73061692a83SJed Brown ierr = PetscFree(ros->work);CHKERRQ(ierr); 731e27a552bSJed Brown PetscFunctionReturn(0); 732e27a552bSJed Brown } 733e27a552bSJed Brown 734e27a552bSJed Brown #undef __FUNCT__ 735e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW" 736e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts) 737e27a552bSJed Brown { 738e27a552bSJed Brown PetscErrorCode ierr; 739e27a552bSJed Brown 740e27a552bSJed Brown PetscFunctionBegin; 741e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 742e27a552bSJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 743e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr); 744e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr); 74561692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr); 746e27a552bSJed Brown PetscFunctionReturn(0); 747e27a552bSJed Brown } 748e27a552bSJed Brown 749e27a552bSJed Brown /* 750e27a552bSJed Brown This defines the nonlinear equation that is to be solved with SNES 751e27a552bSJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 752e27a552bSJed Brown */ 753e27a552bSJed Brown #undef __FUNCT__ 754e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW" 755e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts) 756e27a552bSJed Brown { 75761692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 758e27a552bSJed Brown PetscErrorCode ierr; 759e27a552bSJed Brown 760e27a552bSJed Brown PetscFunctionBegin; 761c17803e7SJed Brown if (ros->stage_explicit) { 762c17803e7SJed Brown ierr = VecAXPBY(ros->Ydot,ros->shift,0.0,X);CHKERRQ(ierr); /* Ydot = shift*X*/ 763c17803e7SJed Brown } else { 76461692a83SJed Brown ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */ 765c17803e7SJed Brown } 76661692a83SJed Brown ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr); /* Ystage = X + Zstage */ 76761692a83SJed Brown ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr); 768e27a552bSJed Brown PetscFunctionReturn(0); 769e27a552bSJed Brown } 770e27a552bSJed Brown 771e27a552bSJed Brown #undef __FUNCT__ 772e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW" 773e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 774e27a552bSJed Brown { 77561692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 776e27a552bSJed Brown PetscErrorCode ierr; 777e27a552bSJed Brown 778e27a552bSJed Brown PetscFunctionBegin; 77961692a83SJed Brown /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 78061692a83SJed Brown ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 781e27a552bSJed Brown PetscFunctionReturn(0); 782e27a552bSJed Brown } 783e27a552bSJed Brown 784e27a552bSJed Brown #undef __FUNCT__ 785e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW" 786e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts) 787e27a552bSJed Brown { 78861692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 78961692a83SJed Brown RosWTableau tab = ros->tableau; 790e27a552bSJed Brown PetscInt s = tab->s; 791e27a552bSJed Brown PetscErrorCode ierr; 792e27a552bSJed Brown 793e27a552bSJed Brown PetscFunctionBegin; 79461692a83SJed Brown if (!ros->tableau) { 795e27a552bSJed Brown ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr); 796e27a552bSJed Brown } 79761692a83SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr); 79861692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr); 79961692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr); 80061692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr); 80161692a83SJed Brown ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr); 80261692a83SJed Brown ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr); 803e27a552bSJed Brown PetscFunctionReturn(0); 804e27a552bSJed Brown } 805e27a552bSJed Brown /*------------------------------------------------------------*/ 806e27a552bSJed Brown 807e27a552bSJed Brown #undef __FUNCT__ 808e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW" 809e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts) 810e27a552bSJed Brown { 81161692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 812e27a552bSJed Brown PetscErrorCode ierr; 81361692a83SJed Brown char rostype[256]; 814e27a552bSJed Brown 815e27a552bSJed Brown PetscFunctionBegin; 816e27a552bSJed Brown ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr); 817e27a552bSJed Brown { 81861692a83SJed Brown RosWTableauLink link; 819e27a552bSJed Brown PetscInt count,choice; 820e27a552bSJed Brown PetscBool flg; 821e27a552bSJed Brown const char **namelist; 82261692a83SJed Brown SNES snes; 82361692a83SJed Brown 82461692a83SJed Brown ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr); 82561692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) ; 826e27a552bSJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 82761692a83SJed Brown for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 82861692a83SJed Brown ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr); 82961692a83SJed Brown ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr); 830e27a552bSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 83161692a83SJed Brown 83261692a83SJed Brown ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr); 83361692a83SJed Brown 83461692a83SJed Brown /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 83561692a83SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 83661692a83SJed Brown if (!((PetscObject)snes)->type_name) { 83761692a83SJed Brown ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 83861692a83SJed Brown } 83961692a83SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 840e27a552bSJed Brown } 841e27a552bSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 842e27a552bSJed Brown PetscFunctionReturn(0); 843e27a552bSJed Brown } 844e27a552bSJed Brown 845e27a552bSJed Brown #undef __FUNCT__ 846e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray" 847e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 848e27a552bSJed Brown { 849e27a552bSJed Brown PetscErrorCode ierr; 850e408995aSJed Brown PetscInt i; 851e408995aSJed Brown size_t left,count; 852e27a552bSJed Brown char *p; 853e27a552bSJed Brown 854e27a552bSJed Brown PetscFunctionBegin; 855e408995aSJed Brown for (i=0,p=buf,left=len; i<n; i++) { 856e408995aSJed Brown ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 857e27a552bSJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 858e27a552bSJed Brown left -= count; 859e27a552bSJed Brown p += count; 860e27a552bSJed Brown *p++ = ' '; 861e27a552bSJed Brown } 862e27a552bSJed Brown p[i ? 0 : -1] = 0; 863e27a552bSJed Brown PetscFunctionReturn(0); 864e27a552bSJed Brown } 865e27a552bSJed Brown 866e27a552bSJed Brown #undef __FUNCT__ 867e27a552bSJed Brown #define __FUNCT__ "TSView_RosW" 868e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer) 869e27a552bSJed Brown { 87061692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 87161692a83SJed Brown RosWTableau tab = ros->tableau; 872e27a552bSJed Brown PetscBool iascii; 873e27a552bSJed Brown PetscErrorCode ierr; 874e27a552bSJed Brown 875e27a552bSJed Brown PetscFunctionBegin; 876e27a552bSJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 877e27a552bSJed Brown if (iascii) { 87861692a83SJed Brown const TSRosWType rostype; 879e408995aSJed Brown PetscInt i; 880e408995aSJed Brown PetscReal abscissa[512]; 881e27a552bSJed Brown char buf[512]; 88261692a83SJed Brown ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr); 88361692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Rosenbrock-W %s\n",rostype);CHKERRQ(ierr); 884e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr); 88561692a83SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A = %s\n",buf);CHKERRQ(ierr); 886e408995aSJed Brown for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 887e408995aSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr); 888e408995aSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr); 889e27a552bSJed Brown } 890e27a552bSJed Brown ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 891e27a552bSJed Brown PetscFunctionReturn(0); 892e27a552bSJed Brown } 893e27a552bSJed Brown 894e27a552bSJed Brown #undef __FUNCT__ 895e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType" 896e27a552bSJed Brown /*@C 89761692a83SJed Brown TSRosWSetType - Set the type of Rosenbrock-W scheme 898e27a552bSJed Brown 899e27a552bSJed Brown Logically collective 900e27a552bSJed Brown 901e27a552bSJed Brown Input Parameter: 902e27a552bSJed Brown + ts - timestepping context 90361692a83SJed Brown - rostype - type of Rosenbrock-W scheme 904e27a552bSJed Brown 905e27a552bSJed Brown Level: intermediate 906e27a552bSJed Brown 907e27a552bSJed Brown .seealso: TSRosWGetType() 908e27a552bSJed Brown @*/ 90961692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype) 910e27a552bSJed Brown { 911e27a552bSJed Brown PetscErrorCode ierr; 912e27a552bSJed Brown 913e27a552bSJed Brown PetscFunctionBegin; 914e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 91561692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr); 916e27a552bSJed Brown PetscFunctionReturn(0); 917e27a552bSJed Brown } 918e27a552bSJed Brown 919e27a552bSJed Brown #undef __FUNCT__ 920e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType" 921e27a552bSJed Brown /*@C 92261692a83SJed Brown TSRosWGetType - Get the type of Rosenbrock-W scheme 923e27a552bSJed Brown 924e27a552bSJed Brown Logically collective 925e27a552bSJed Brown 926e27a552bSJed Brown Input Parameter: 927e27a552bSJed Brown . ts - timestepping context 928e27a552bSJed Brown 929e27a552bSJed Brown Output Parameter: 93061692a83SJed Brown . rostype - type of Rosenbrock-W scheme 931e27a552bSJed Brown 932e27a552bSJed Brown Level: intermediate 933e27a552bSJed Brown 934e27a552bSJed Brown .seealso: TSRosWGetType() 935e27a552bSJed Brown @*/ 93661692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype) 937e27a552bSJed Brown { 938e27a552bSJed Brown PetscErrorCode ierr; 939e27a552bSJed Brown 940e27a552bSJed Brown PetscFunctionBegin; 941e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 94261692a83SJed Brown ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr); 943e27a552bSJed Brown PetscFunctionReturn(0); 944e27a552bSJed Brown } 945e27a552bSJed Brown 946e27a552bSJed Brown #undef __FUNCT__ 94761692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian" 948e27a552bSJed Brown /*@C 94961692a83SJed Brown TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 950e27a552bSJed Brown 951e27a552bSJed Brown Logically collective 952e27a552bSJed Brown 953e27a552bSJed Brown Input Parameter: 954e27a552bSJed Brown + ts - timestepping context 95561692a83SJed Brown - flg - PETSC_TRUE to recompute the Jacobian at each stage 956e27a552bSJed Brown 957e27a552bSJed Brown Level: intermediate 958e27a552bSJed Brown 959e27a552bSJed Brown .seealso: TSRosWGetType() 960e27a552bSJed Brown @*/ 96161692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg) 962e27a552bSJed Brown { 963e27a552bSJed Brown PetscErrorCode ierr; 964e27a552bSJed Brown 965e27a552bSJed Brown PetscFunctionBegin; 966e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 96761692a83SJed Brown ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 968e27a552bSJed Brown PetscFunctionReturn(0); 969e27a552bSJed Brown } 970e27a552bSJed Brown 971e27a552bSJed Brown EXTERN_C_BEGIN 972e27a552bSJed Brown #undef __FUNCT__ 973e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW" 97461692a83SJed Brown PetscErrorCode TSRosWGetType_RosW(TS ts,const TSRosWType *rostype) 975e27a552bSJed Brown { 97661692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 977e27a552bSJed Brown PetscErrorCode ierr; 978e27a552bSJed Brown 979e27a552bSJed Brown PetscFunctionBegin; 98061692a83SJed Brown if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);} 98161692a83SJed Brown *rostype = ros->tableau->name; 982e27a552bSJed Brown PetscFunctionReturn(0); 983e27a552bSJed Brown } 984e27a552bSJed Brown #undef __FUNCT__ 985e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW" 98661692a83SJed Brown PetscErrorCode TSRosWSetType_RosW(TS ts,const TSRosWType rostype) 987e27a552bSJed Brown { 98861692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 989e27a552bSJed Brown PetscErrorCode ierr; 990e27a552bSJed Brown PetscBool match; 99161692a83SJed Brown RosWTableauLink link; 992e27a552bSJed Brown 993e27a552bSJed Brown PetscFunctionBegin; 99461692a83SJed Brown if (ros->tableau) { 99561692a83SJed Brown ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr); 996e27a552bSJed Brown if (match) PetscFunctionReturn(0); 997e27a552bSJed Brown } 99861692a83SJed Brown for (link = RosWTableauList; link; link=link->next) { 99961692a83SJed Brown ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr); 1000e27a552bSJed Brown if (match) { 1001e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 100261692a83SJed Brown ros->tableau = &link->tab; 1003e27a552bSJed Brown PetscFunctionReturn(0); 1004e27a552bSJed Brown } 1005e27a552bSJed Brown } 100661692a83SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype); 1007e27a552bSJed Brown PetscFunctionReturn(0); 1008e27a552bSJed Brown } 100961692a83SJed Brown 1010e27a552bSJed Brown #undef __FUNCT__ 101161692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW" 101261692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg) 1013e27a552bSJed Brown { 101461692a83SJed Brown TS_RosW *ros = (TS_RosW*)ts->data; 1015e27a552bSJed Brown 1016e27a552bSJed Brown PetscFunctionBegin; 101761692a83SJed Brown ros->recompute_jacobian = flg; 1018e27a552bSJed Brown PetscFunctionReturn(0); 1019e27a552bSJed Brown } 1020e27a552bSJed Brown EXTERN_C_END 1021e27a552bSJed Brown 1022e27a552bSJed Brown /* ------------------------------------------------------------ */ 1023e27a552bSJed Brown /*MC 1024e27a552bSJed Brown TSRosW - ODE solver using Rosenbrock-W schemes 1025e27a552bSJed Brown 1026e27a552bSJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1027e27a552bSJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1028e27a552bSJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1029e27a552bSJed Brown 1030e27a552bSJed Brown Notes: 103161692a83SJed Brown This method currently only works with autonomous ODE and DAE. 103261692a83SJed Brown 103361692a83SJed Brown Developer notes: 103461692a83SJed Brown Rosenbrock-W methods are typically specified for autonomous ODE 103561692a83SJed Brown 103661692a83SJed Brown $ xdot = f(x) 103761692a83SJed Brown 103861692a83SJed Brown by the stage equations 103961692a83SJed Brown 104061692a83SJed Brown $ k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 104161692a83SJed Brown 104261692a83SJed Brown and step completion formula 104361692a83SJed Brown 104461692a83SJed Brown $ x_1 = x_0 + sum_j b_j k_j 104561692a83SJed Brown 104661692a83SJed Brown with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x) 104761692a83SJed Brown and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 104861692a83SJed Brown we define new variables for the stage equations 104961692a83SJed Brown 105061692a83SJed Brown $ y_i = gamma_ij k_j 105161692a83SJed Brown 105261692a83SJed Brown The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 105361692a83SJed Brown 105461692a83SJed Brown $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i} 105561692a83SJed Brown 105661692a83SJed Brown to rewrite the method as 105761692a83SJed Brown 105861692a83SJed 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 105961692a83SJed Brown $ x_1 = x_0 + sum_j bt_j y_j 106061692a83SJed Brown 106161692a83SJed Brown where we have introduced the mass matrix M. Continue by defining 106261692a83SJed Brown 106361692a83SJed Brown $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 106461692a83SJed Brown 106561692a83SJed Brown or, more compactly in tensor notation 106661692a83SJed Brown 106761692a83SJed Brown $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 106861692a83SJed Brown 106961692a83SJed Brown Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 107061692a83SJed Brown stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 107161692a83SJed Brown equation 107261692a83SJed Brown 107361692a83SJed Brown $ g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 107461692a83SJed Brown 107561692a83SJed Brown with initial guess y_i = 0. 1076e27a552bSJed Brown 1077e27a552bSJed Brown Level: beginner 1078e27a552bSJed Brown 1079e27a552bSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSRosWRegister() 1080e27a552bSJed Brown 1081e27a552bSJed Brown M*/ 1082e27a552bSJed Brown EXTERN_C_BEGIN 1083e27a552bSJed Brown #undef __FUNCT__ 1084e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW" 1085e27a552bSJed Brown PetscErrorCode TSCreate_RosW(TS ts) 1086e27a552bSJed Brown { 108761692a83SJed Brown TS_RosW *ros; 1088e27a552bSJed Brown PetscErrorCode ierr; 1089e27a552bSJed Brown 1090e27a552bSJed Brown PetscFunctionBegin; 1091e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1092e27a552bSJed Brown ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1093e27a552bSJed Brown #endif 1094e27a552bSJed Brown 1095e27a552bSJed Brown ts->ops->reset = TSReset_RosW; 1096e27a552bSJed Brown ts->ops->destroy = TSDestroy_RosW; 1097e27a552bSJed Brown ts->ops->view = TSView_RosW; 1098e27a552bSJed Brown ts->ops->setup = TSSetUp_RosW; 1099e27a552bSJed Brown ts->ops->step = TSStep_RosW; 1100e27a552bSJed Brown ts->ops->interpolate = TSInterpolate_RosW; 11011c3436cfSJed Brown ts->ops->evaluatestep = TSEvaluateStep_RosW; 1102e27a552bSJed Brown ts->ops->setfromoptions = TSSetFromOptions_RosW; 1103e27a552bSJed Brown ts->ops->snesfunction = SNESTSFormFunction_RosW; 1104e27a552bSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1105e27a552bSJed Brown 110661692a83SJed Brown ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr); 110761692a83SJed Brown ts->data = (void*)ros; 1108e27a552bSJed Brown 1109e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr); 1110e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr); 111161692a83SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr); 1112e27a552bSJed Brown PetscFunctionReturn(0); 1113e27a552bSJed Brown } 1114e27a552bSJed Brown EXTERN_C_END 1115