Feagin's RK14(12) - a 14th-order method with an embedded 12th-order method.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rkf1412_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t |
initial time |
||
real(kind=wp), | intent(in), | dimension(me%n) | :: | x |
initial state |
|
real(kind=wp), | intent(in) | :: | h |
time step |
||
real(kind=wp), | intent(out), | dimension(me%n) | :: | xf |
state at time |
|
real(kind=wp), | intent(out), | dimension(me%n) | :: | terr |
truncation error estimate |
subroutine rkf1412(me,t,x,h,xf,terr) implicit none class(rkf1412_class),intent(inout) :: me real(wp),intent(in) :: t !! initial time real(wp),dimension(me%n),intent(in) :: x !! initial state real(wp),intent(in) :: h !! time step real(wp),dimension(me%n),intent(out) :: xf !! state at time `t+h` real(wp),dimension(me%n),intent(out) :: terr !! truncation error estimate real(wp),parameter :: a0 = 0.000000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: a1 = 0.111111111111111111111111111111111111111111111111111111111111_wp real(wp),parameter :: a2 = 0.555555555555555555555555555555555555555555555555555555555556_wp real(wp),parameter :: a3 = 0.833333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: a4 = 0.333333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: a5 = 1.00000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: a6 = 0.669986979272772921764683785505998513938845229638460353285142_wp real(wp),parameter :: a7 = 0.297068384213818357389584716808219413223332094698915687379168_wp real(wp),parameter :: a8 = 0.727272727272727272727272727272727272727272727272727272727273_wp real(wp),parameter :: a9 = 0.140152799042188765276187487966946717629806463082532936287323_wp real(wp),parameter :: a10 = 0.700701039770150737151099854830749337941407049265546408969222_wp real(wp),parameter :: a11 = 0.363636363636363636363636363636363636363636363636363636363636_wp real(wp),parameter :: a12 = 0.263157894736842105263157894736842105263157894736842105263158_wp real(wp),parameter :: a13 = 0.0392172246650270859125196642501208648863714315266128052078483_wp real(wp),parameter :: a14 = 0.812917502928376762983393159278036506189612372617238550774312_wp real(wp),parameter :: a15 = 0.166666666666666666666666666666666666666666666666666666666667_wp real(wp),parameter :: a16 = 0.900000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: a17 = 0.0641299257451966923312771193896682809481096651615083225402924_wp real(wp),parameter :: a18 = 0.204149909283428848927744634301023405027149505241333751628870_wp real(wp),parameter :: a19 = 0.395350391048760565615671369827324372352227297456659450554577_wp real(wp),parameter :: a20 = 0.604649608951239434384328630172675627647772702543340549445423_wp real(wp),parameter :: a21 = 0.795850090716571151072255365698976594972850494758666248371130_wp real(wp),parameter :: a22 = 0.935870074254803307668722880610331719051890334838491677459708_wp real(wp),parameter :: a23 = 0.166666666666666666666666666666666666666666666666666666666667_wp real(wp),parameter :: a24 = 0.812917502928376762983393159278036506189612372617238550774312_wp real(wp),parameter :: a25 = 0.0392172246650270859125196642501208648863714315266128052078483_wp real(wp),parameter :: a26 = 0.363636363636363636363636363636363636363636363636363636363636_wp real(wp),parameter :: a27 = 0.700701039770150737151099854830749337941407049265546408969222_wp real(wp),parameter :: a28 = 0.140152799042188765276187487966946717629806463082532936287323_wp real(wp),parameter :: a29 = 0.297068384213818357389584716808219413223332094698915687379168_wp real(wp),parameter :: a30 = 0.669986979272772921764683785505998513938845229638460353285142_wp real(wp),parameter :: a31 = 0.333333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: a32 = 0.555555555555555555555555555555555555555555555555555555555556_wp real(wp),parameter :: a33 = 0.111111111111111111111111111111111111111111111111111111111111_wp real(wp),parameter :: a34 = 1.00000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c0 = 0.0178571428571428571428571428571428571428571428571428571428571_wp real(wp),parameter :: c1 = 0.00585937500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c2 = 0.0117187500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c3 = 0.0_wp real(wp),parameter :: c4 = 0.0175781250000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c5 = 0.0_wp real(wp),parameter :: c6 = 0.0234375000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c7 = 0.0292968750000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c8 = 0.0_wp real(wp),parameter :: c9 = 0.0351562500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c10 = 0.0410156250000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c11 = 0.0468750000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c12 = 0.0_wp real(wp),parameter :: c13 = 0.0527343750000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c14 = 0.0585937500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c15 = 0.0644531250000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c16 = 0.0_wp real(wp),parameter :: c17 = 0.105352113571753019691496032887878162227673083080523884041670_wp real(wp),parameter :: c18 = 0.170561346241752182382120338553874085887555487802790804737501_wp real(wp),parameter :: c19 = 0.206229397329351940783526485701104894741914286259542454077972_wp real(wp),parameter :: c20 = 0.206229397329351940783526485701104894741914286259542454077972_wp real(wp),parameter :: c21 = 0.170561346241752182382120338553874085887555487802790804737501_wp real(wp),parameter :: c22 = 0.105352113571753019691496032887878162227673083080523884041670_wp real(wp),parameter :: c23 = -0.0644531250000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c24 = -0.0585937500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c25 = -0.0527343750000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c26 = -0.0468750000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c27 = -0.0410156250000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c28 = -0.0351562500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c29 = -0.0292968750000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c30 = -0.0234375000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c31 = -0.0175781250000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c32 = -0.0117187500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c33 = -0.00585937500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: c34 = 0.0178571428571428571428571428571428571428571428571428571428571_wp real(wp),parameter :: b10 = 0.111111111111111111111111111111111111111111111111111111111111_wp real(wp),parameter :: b20 = -0.833333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: b21 = 1.38888888888888888888888888888888888888888888888888888888889_wp real(wp),parameter :: b30 = 0.208333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: b31 = 0.0_wp real(wp),parameter :: b32 = 0.625000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b40 = 0.193333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: b41 = 0.0_wp real(wp),parameter :: b42 = 0.220000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b43 = -0.0800000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b50 = 0.100000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b51 = 0.0_wp real(wp),parameter :: b52 = 0.0_wp real(wp),parameter :: b53 = 0.400000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b54 = 0.500000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b60 = 0.103484561636679776672993546511910344499744798201971316606663_wp real(wp),parameter :: b61 = 0.0_wp real(wp),parameter :: b62 = 0.0_wp real(wp),parameter :: b63 = 0.122068887306407222589644082868962077139592714834162134741275_wp real(wp),parameter :: b64 = 0.482574490331246622475134780125688112865919023850168049679402_wp real(wp),parameter :: b65 = -0.0381409600015606999730886240005620205664113072478411477421970_wp real(wp),parameter :: b70 = 0.124380526654094412881516420868799316268491466359671423163289_wp real(wp),parameter :: b71 = 0.0_wp real(wp),parameter :: b72 = 0.0_wp real(wp),parameter :: b73 = 0.0_wp real(wp),parameter :: b74 = 0.226120282197584301422238662979202901196752320742633143965145_wp real(wp),parameter :: b75 = 0.0137885887618080880607695837016477814530969417491493385363543_wp real(wp),parameter :: b76 = -0.0672210133996684449749399507414305856950086341525382182856200_wp real(wp),parameter :: b80 = 0.0936919065659673815530885456083005933866349695217750085655603_wp real(wp),parameter :: b81 = 0.0_wp real(wp),parameter :: b82 = 0.0_wp real(wp),parameter :: b83 = 0.0_wp real(wp),parameter :: b84 = 0.0_wp real(wp),parameter :: b85 = -0.00613406843450510987229498995641664735620914507128858871007099_wp real(wp),parameter :: b86 = 0.216019825625503063708860097659866573490979433278117320188668_wp real(wp),parameter :: b87 = 0.423695063515761937337619073960976753205867469544123532683116_wp real(wp),parameter :: b90 = 0.0838479812409052664616968791372814085980533139224911131069335_wp real(wp),parameter :: b91 = 0.0_wp real(wp),parameter :: b92 = 0.0_wp real(wp),parameter :: b93 = 0.0_wp real(wp),parameter :: b94 = 0.0_wp real(wp),parameter :: b95 = -0.0117949367100973814319755056031295775367961960590736150777613_wp real(wp),parameter :: b96 = -0.247299020568812652339473838743194598325992840353340132697498_wp real(wp),parameter :: b97 = 0.0978080858367729012259313014081291665503740655476733940756599_wp real(wp),parameter :: b98 = 0.217590689243420631360008651767860318344168120024782176879989_wp real(wp),parameter :: b100 = 0.0615255359769428227954562389614314714333423969064821107453940_wp real(wp),parameter :: b101 = 0.0_wp real(wp),parameter :: b102 = 0.0_wp real(wp),parameter :: b103 = 0.0_wp real(wp),parameter :: b104 = 0.0_wp real(wp),parameter :: b105 = 0.00592232780324503308042990005798046524738389560444257136834990_wp real(wp),parameter :: b106 = 0.470326159963841112217224303205894113455362530746108825010848_wp real(wp),parameter :: b107 = 0.299688863848679000853981837096192399136831121671781279184194_wp real(wp),parameter :: b108 = -0.247656877593994914689992276329810825853958069263947095548189_wp real(wp),parameter :: b109 = 0.110895029771437682893999851839061714522445173600678718208625_wp real(wp),parameter :: b110 = 0.0419700073362782579861792864787277787213483656543104611245994_wp real(wp),parameter :: b111 = 0.0_wp real(wp),parameter :: b112 = 0.0_wp real(wp),parameter :: b113 = 0.0_wp real(wp),parameter :: b114 = 0.0_wp real(wp),parameter :: b115 = -0.00317987696266205093901912847692712407988609169703103952205634_wp real(wp),parameter :: b116 = 0.806397714906192077260821711520379506393543111567419750119748_wp real(wp),parameter :: b117 = 0.0975983126412388979093522850684288851314672048003054550357187_wp real(wp),parameter :: b118 = 0.778575578158398909027512446452927238999763460594181964958853_wp real(wp),parameter :: b119 = 0.204890423831599428189499202098105603312029235081420653574829_wp real(wp),parameter :: b1110 = -1.56261579627468188307070943950527825211462892236424360892806_wp real(wp),parameter :: b120 = 0.0437726782233730163574465242495339811688214967071614123256973_wp real(wp),parameter :: b121 = 0.0_wp real(wp),parameter :: b122 = 0.0_wp real(wp),parameter :: b123 = 0.0_wp real(wp),parameter :: b124 = 0.0_wp real(wp),parameter :: b125 = 0.0_wp real(wp),parameter :: b126 = 0.0_wp real(wp),parameter :: b127 = 0.0_wp real(wp),parameter :: b128 = 0.00624365027520195208794358628580933625281631216903095917201250_wp real(wp),parameter :: b129 = 0.200043097109577314994435165469647856829066232218264969608768_wp real(wp),parameter :: b1210 = -0.00805328367804983036823857162048902911923392887337029314844206_wp real(wp),parameter :: b1211 = 0.0211517528067396521915711903523399601316877825157550573051221_wp real(wp),parameter :: b130 = 0.0283499250363514563095023591920717312247137654896477097768495_wp real(wp),parameter :: b131 = 0.0_wp real(wp),parameter :: b132 = 0.0_wp real(wp),parameter :: b133 = 0.0_wp real(wp),parameter :: b134 = 0.0_wp real(wp),parameter :: b135 = 0.0_wp real(wp),parameter :: b136 = 0.0_wp real(wp),parameter :: b137 = 0.0_wp real(wp),parameter :: b138 = 0.00249163204855817407538949148805995149459884653585417680098222_wp real(wp),parameter :: b139 = 0.0230138787854593149638399846373742768772087122638142234223658_wp real(wp),parameter :: b1310 = -0.00322155956692977098724476092467120878189463604760620461043308_wp real(wp),parameter :: b1311 = 0.00988442549447664668946335414487885256040819982786014648129297_wp real(wp),parameter :: b1312 = -0.0213010771328887351384307642875927384886634565429572466632092_wp real(wp),parameter :: b140 = 0.343511894290243001049432234735147943083353174980701426268122_wp real(wp),parameter :: b141 = 0.0_wp real(wp),parameter :: b142 = 0.0_wp real(wp),parameter :: b143 = 0.0_wp real(wp),parameter :: b144 = 0.0_wp real(wp),parameter :: b145 = 0.0_wp real(wp),parameter :: b146 = 0.0_wp real(wp),parameter :: b147 = 0.0_wp real(wp),parameter :: b148 = 0.210451912023627385609097011999010655788807405225626700040882_wp real(wp),parameter :: b149 = 1.03427452057230411936482926828825709938667999698324740166559_wp real(wp),parameter :: b1410 = 0.00600303645864422487051240448206640574939078092406156945568306_wp real(wp),parameter :: b1411 = 0.855938125099619537578012106002407728915062652616416005816477_wp real(wp),parameter :: b1412 = -0.977235005036766810872264852372525633013107656892839677696022_wp real(wp),parameter :: b1413 = -0.660026980479294694616225013856327693720573981219974874776419_wp real(wp),parameter :: b150 = -0.0143574001672168069538206399935076366657755954378399880691949_wp real(wp),parameter :: b151 = 0.0_wp real(wp),parameter :: b152 = 0.0_wp real(wp),parameter :: b153 = 0.0_wp real(wp),parameter :: b154 = 0.0_wp real(wp),parameter :: b155 = 0.0_wp real(wp),parameter :: b156 = 0.0_wp real(wp),parameter :: b157 = 0.0_wp real(wp),parameter :: b158 = -0.0366253270049039970293685796848974791733119081733552207318285_wp real(wp),parameter :: b159 = 0.0350254975636213681976849406979846524346789082471103574920148_wp real(wp),parameter :: b1510 = 0.0360946016362113508931786658758335239823689929864237671348749_wp real(wp),parameter :: b1511 = -0.0265219967553681106351595946834601923649627012457464284442911_wp real(wp),parameter :: b1512 = 0.0445699011305698119638911537508839908104336323082226770910408_wp real(wp),parameter :: b1513 = 0.124343093331358243286225595741786448038973408895106741855721_wp real(wp),parameter :: b1514 = 0.00413829693239480694403512496204335960426192908674476033832967_wp real(wp),parameter :: b160 = 0.356032404425120290975609116398089176264106222379748802654822_wp real(wp),parameter :: b161 = 0.0_wp real(wp),parameter :: b162 = 0.0_wp real(wp),parameter :: b163 = 0.0_wp real(wp),parameter :: b164 = 0.0_wp real(wp),parameter :: b165 = 0.0_wp real(wp),parameter :: b166 = 0.0_wp real(wp),parameter :: b167 = 0.0_wp real(wp),parameter :: b168 = -0.450192758947562595966821779075956175110645100214763601190349_wp real(wp),parameter :: b169 = 0.430527907083710898626656292808782917793030154094709462877146_wp real(wp),parameter :: b1610 = 0.511973029011022237668556960394071692077125787030651386389972_wp real(wp),parameter :: b1611 = 0.908303638886404260390159124638110213997496214819904630546596_wp real(wp),parameter :: b1612 = -1.23921093371933931757372469151534028854413889248605726186520_wp real(wp),parameter :: b1613 = -0.649048661671761465141672348879062553905402831967191097656668_wp real(wp),parameter :: b1614 = 0.251708904586819292210480529948970541404887852931447491219418_wp real(wp),parameter :: b1615 = 0.779906470345586398810756795282334476023540593411550187024263_wp real(wp),parameter :: b170 = 0.0130935687406513066406881206418834980127470438213192487844956_wp real(wp),parameter :: b171 = 0.0_wp real(wp),parameter :: b172 = 0.0_wp real(wp),parameter :: b173 = 0.0_wp real(wp),parameter :: b174 = 0.0_wp real(wp),parameter :: b175 = 0.0_wp real(wp),parameter :: b176 = 0.0_wp real(wp),parameter :: b177 = 0.0_wp real(wp),parameter :: b178 = 0.0_wp real(wp),parameter :: b179 = 0.0_wp real(wp),parameter :: b1710 = 0.0_wp real(wp),parameter :: b1711 = 0.0_wp real(wp),parameter :: b1712 = -0.0000932053067985113945908461962767108237858631509684667142124826_wp real(wp),parameter :: b1713 = 0.0505374334262299359640090443138590726770942344716122381702746_wp real(wp),parameter :: b1714 = 8.04470341944487979109579109610197797641311868930865361048975e-7_wp real(wp),parameter :: b1715 = 0.000591726029494171190528755742777717259844340971924321528178248_wp real(wp),parameter :: b1716 = -4.01614722154557337064691684906375587732264247950093804676867e-7_wp real(wp),parameter :: b180 = 0.0207926484466053012541944544000765652167255206144373407979758_wp real(wp),parameter :: b181 = 0.0_wp real(wp),parameter :: b182 = 0.0_wp real(wp),parameter :: b183 = 0.0_wp real(wp),parameter :: b184 = 0.0_wp real(wp),parameter :: b185 = 0.0_wp real(wp),parameter :: b186 = 0.0_wp real(wp),parameter :: b187 = 0.0_wp real(wp),parameter :: b188 = 0.0_wp real(wp),parameter :: b189 = 0.0_wp real(wp),parameter :: b1810 = 0.0_wp real(wp),parameter :: b1811 = 0.0_wp real(wp),parameter :: b1812 = 0.000582695918800085915101902697837284108951406103029871570103075_wp real(wp),parameter :: b1813 = -0.00801700732358815939083342186525852746640558465919633524655451_wp real(wp),parameter :: b1814 = 4.03847643847136940375170821743560570484117290330895506618968e-6_wp real(wp),parameter :: b1815 = 0.0854609998055506144225056114567535602510114622033622491802597_wp real(wp),parameter :: b1816 = -2.04486480935804242706707569691004307904442837552677456232848e-6_wp real(wp),parameter :: b1817 = 0.105328578824431893399799402979093997354240904235172843146582_wp real(wp),parameter :: b190 = 1.40153449795736021415446247355771306718486452917597731683689_wp real(wp),parameter :: b191 = 0.0_wp real(wp),parameter :: b192 = 0.0_wp real(wp),parameter :: b193 = 0.0_wp real(wp),parameter :: b194 = 0.0_wp real(wp),parameter :: b195 = 0.0_wp real(wp),parameter :: b196 = 0.0_wp real(wp),parameter :: b197 = 0.0_wp real(wp),parameter :: b198 = 0.0_wp real(wp),parameter :: b199 = 0.0_wp real(wp),parameter :: b1910 = 0.0_wp real(wp),parameter :: b1911 = 0.0_wp real(wp),parameter :: b1912 = -0.230252000984221261616272410367415621261130298274455611733277_wp real(wp),parameter :: b1913 = -7.21106840466912905659582237106874247165856493509961561958267_wp real(wp),parameter :: b1914 = 0.00372901560694836335236995327852132340217759566678662385552634_wp real(wp),parameter :: b1915 = -4.71415495727125020678778179392224757011323373221820091641216_wp real(wp),parameter :: b1916 = -0.00176367657545349242053841995032797673574903886695600132759652_wp real(wp),parameter :: b1917 = 7.64130548038698765563029310880237651185173367813936997648198_wp real(wp),parameter :: b1918 = 3.50602043659751834989896082949744710968212949893375368243588_wp real(wp),parameter :: b200 = 11.9514650694120686799372385830716401674473610826553517297976_wp real(wp),parameter :: b201 = 0.0_wp real(wp),parameter :: b202 = 0.0_wp real(wp),parameter :: b203 = 0.0_wp real(wp),parameter :: b204 = 0.0_wp real(wp),parameter :: b205 = 0.0_wp real(wp),parameter :: b206 = 0.0_wp real(wp),parameter :: b207 = 0.0_wp real(wp),parameter :: b208 = 0.0_wp real(wp),parameter :: b209 = 0.0_wp real(wp),parameter :: b2010 = 0.0_wp real(wp),parameter :: b2011 = 0.0_wp real(wp),parameter :: b2012 = 7.79480932108175968783516700231764388220284279598980948538579_wp real(wp),parameter :: b2013 = -56.4501393867325792523560991120904281440468100061340556540132_wp real(wp),parameter :: b2014 = 0.0912376306930644901344530449290276645709607450403673704844997_wp real(wp),parameter :: b2015 = -12.7336279925434886201945524309199275038162717529918963305155_wp real(wp),parameter :: b2016 = -0.0396895921904719712313542810939736674712383070433147873009352_wp real(wp),parameter :: b2017 = 54.4392141883570886996225765155307791861438378423305337073797_wp real(wp),parameter :: b2018 = -3.64411637921569236846406990361350645806721478409266709351203_wp real(wp),parameter :: b2019 = -0.804503249910509910899030787958579499315694913210787878260459_wp real(wp),parameter :: b210 = -148.809426507100488427838868268647625561930612082148597076690_wp real(wp),parameter :: b211 = 0.0_wp real(wp),parameter :: b212 = 0.0_wp real(wp),parameter :: b213 = 0.0_wp real(wp),parameter :: b214 = 0.0_wp real(wp),parameter :: b215 = 0.0_wp real(wp),parameter :: b216 = 0.0_wp real(wp),parameter :: b217 = 0.0_wp real(wp),parameter :: b218 = 0.0_wp real(wp),parameter :: b219 = 0.0_wp real(wp),parameter :: b2110 = 0.0_wp real(wp),parameter :: b2111 = 0.0_wp real(wp),parameter :: b2112 = -91.7295278291256484357935662402321623495228729036354276506427_wp real(wp),parameter :: b2113 = 707.656144971598359834575719286335716154821128966649565194286_wp real(wp),parameter :: b2114 = -1.10563611857482440905296961311590930801338308942637769555540_wp real(wp),parameter :: b2115 = 176.134591883811372587859898076055660406999516762301689616841_wp real(wp),parameter :: b2116 = 0.491384824214880662268898345164454557416884631402764792538746_wp real(wp),parameter :: b2117 = -684.278000449814944358237535610895081956077167893600278300805_wp real(wp),parameter :: b2118 = 27.9910604998398258984224332124380407446002518400668657974589_wp real(wp),parameter :: b2119 = 13.1939710030282333443670964371153238435064159623744975073252_wp real(wp),parameter :: b2120 = 1.25128781283980445450114974148056006317268830077396406361417_wp real(wp),parameter :: b220 = -9.67307946948196763644126118433219395839951408571877262880482_wp real(wp),parameter :: b221 = 0.0_wp real(wp),parameter :: b222 = 0.0_wp real(wp),parameter :: b223 = 0.0_wp real(wp),parameter :: b224 = 0.0_wp real(wp),parameter :: b225 = 0.0_wp real(wp),parameter :: b226 = 0.0_wp real(wp),parameter :: b227 = 0.0_wp real(wp),parameter :: b228 = 0.0_wp real(wp),parameter :: b229 = 0.0_wp real(wp),parameter :: b2210 = 0.0_wp real(wp),parameter :: b2211 = 0.0_wp real(wp),parameter :: b2212 = -4.46990150858505531443846227701960360497830681408751431146712_wp real(wp),parameter :: b2213 = 45.5127128690952681968241950400052751178905907817398483534845_wp real(wp),parameter :: b2214 = -0.0713085086183826912791492024438246129930559805352394367050813_wp real(wp),parameter :: b2215 = 11.2273614068412741582590624479939384207826800776794485051540_wp real(wp),parameter :: b2216 = 0.126244376717622724516237912909138809361786889819105426371393_wp real(wp),parameter :: b2217 = -43.5439339549483313605810624907242107623814304467621407753424_wp real(wp),parameter :: b2218 = 0.787174307543058978398792994996550902064546091443233850464377_wp real(wp),parameter :: b2219 = 0.532264696744684215669300708603886690785395776821503851830821_wp real(wp),parameter :: b2220 = 0.422422733996325326010225127471388772575086538809603346825334_wp real(wp),parameter :: b2221 = 0.0859131249503067107308438031499859443441115056294154956487671_wp real(wp),parameter :: b230 = -10.0664032447054702403396606900426891472202824757968765569183_wp real(wp),parameter :: b231 = 0.0_wp real(wp),parameter :: b232 = 0.0_wp real(wp),parameter :: b233 = 0.0_wp real(wp),parameter :: b234 = 0.0_wp real(wp),parameter :: b235 = 0.0_wp real(wp),parameter :: b236 = 0.0_wp real(wp),parameter :: b237 = 0.0_wp real(wp),parameter :: b238 = -0.0366253270049039970293685796848974791733119081733552207318285_wp real(wp),parameter :: b239 = 0.0350254975636213681976849406979846524346789082471103574920148_wp real(wp),parameter :: b2310 = 0.0360946016362113508931786658758335239823689929864237671348749_wp real(wp),parameter :: b2311 = -0.0265219967553681106351595946834601923649627012457464284442911_wp real(wp),parameter :: b2312 = -6.27088972181464143590553149478871603839356122957396018530209_wp real(wp),parameter :: b2313 = 48.2079237442562989090702103008195063923492593141636117832993_wp real(wp),parameter :: b2314 = -0.0694471689136165640882395180583732834557754169149088630301342_wp real(wp),parameter :: b2315 = 12.6810690204850295698341370913609807066108483811412127009785_wp real(wp),parameter :: b2316 = 0.0119671168968323754838161435501011294100927813964199613229864_wp real(wp),parameter :: b2317 = -46.7249764992482408003358268242662695593201321659795608950429_wp real(wp),parameter :: b2318 = 1.33029613326626711314710039298216591399033511191227101321435_wp real(wp),parameter :: b2319 = 1.00766787503398298353438903619926657771162717793661719708370_wp real(wp),parameter :: b2320 = 0.0209512051933665091664122388475480702892770753864487241177616_wp real(wp),parameter :: b2321 = 0.0210134706331264177317735424331396407424412188443757490871603_wp real(wp),parameter :: b2322 = 0.00952196014417121794175101542454575907376360233658356240547761_wp real(wp),parameter :: b240 = -409.478081677743708772589097409370357624424341606752069725341_wp real(wp),parameter :: b241 = 0.0_wp real(wp),parameter :: b242 = 0.0_wp real(wp),parameter :: b243 = 0.0_wp real(wp),parameter :: b244 = 0.0_wp real(wp),parameter :: b245 = 0.0_wp real(wp),parameter :: b246 = 0.0_wp real(wp),parameter :: b247 = 0.0_wp real(wp),parameter :: b248 = 0.210451912023627385609097011999010655788807405225626700040882_wp real(wp),parameter :: b249 = 1.03427452057230411936482926828825709938667999698324740166559_wp real(wp),parameter :: b2410 = 0.00600303645864422487051240448206640574939078092406156945568306_wp real(wp),parameter :: b2411 = 0.855938125099619537578012106002407728915062652616416005816477_wp real(wp),parameter :: b2412 = -250.516998547447860492777657729316130386584050420782075966990_wp real(wp),parameter :: b2413 = 1946.42466652388427766053750328264758595829850895761428240231_wp real(wp),parameter :: b2414 = -3.04503882102310365506105809086860882786950544097602101685174_wp real(wp),parameter :: b2415 = 490.626379528281713521208265299168083841598542274061671576230_wp real(wp),parameter :: b2416 = 1.56647589531270907115484067013597445739595615245966775329993_wp real(wp),parameter :: b2417 = -1881.97428994011173362217267377035870619215906638453056643641_wp real(wp),parameter :: b2418 = 75.2592224724847175278837713643303149821620618914245864351135_wp real(wp),parameter :: b2419 = 34.5734356980331067622434344736554689696728644793551014989002_wp real(wp),parameter :: b2420 = 3.21147679440968961435417361847073755169022966748891627882572_wp real(wp),parameter :: b2421 = -0.460408041738414391307201404237058848867245095265382820823055_wp real(wp),parameter :: b2422 = -0.0870718339841810522431884137957986245724252047388936572215438_wp real(wp),parameter :: b2423 = -7.39351814158303067567016952195521063999185773249132944724553_wp real(wp),parameter :: b250 = 3.43347475853550878921093496257596781120623891072008459930197_wp real(wp),parameter :: b251 = 0.0_wp real(wp),parameter :: b252 = 0.0_wp real(wp),parameter :: b253 = 0.0_wp real(wp),parameter :: b254 = 0.0_wp real(wp),parameter :: b255 = 0.0_wp real(wp),parameter :: b256 = 0.0_wp real(wp),parameter :: b257 = 0.0_wp real(wp),parameter :: b258 = 0.00249163204855817407538949148805995149459884653585417680098222_wp real(wp),parameter :: b259 = 0.0230138787854593149638399846373742768772087122638142234223658_wp real(wp),parameter :: b2510 = -0.00322155956692977098724476092467120878189463604760620461043308_wp real(wp),parameter :: b2511 = 0.00988442549447664668946335414487885256040819982786014648129297_wp real(wp),parameter :: b2512 = 2.16252799377922507788307841904757354045759225335732707916530_wp real(wp),parameter :: b2513 = -16.2699864546457421328065640660139489006987552040228852402716_wp real(wp),parameter :: b2514 = -0.128534502120524552843583417470935010538029037542654506231743_wp real(wp),parameter :: b2515 = -8.98915042666504253089307820833379330486511746063552853023189_wp real(wp),parameter :: b2516 = -0.00348595363232025333387080201851013650192401767250513765000963_wp real(wp),parameter :: b2517 = 15.7936194113339807536235187388695574135853387025139738341334_wp real(wp),parameter :: b2518 = -0.574403330914095065628165482017335820148383663195675408024658_wp real(wp),parameter :: b2519 = -0.345602039021393296692722496608124982535237228827655306030152_wp real(wp),parameter :: b2520 = -0.00662241490206585091731619991383757781133067992707418687587487_wp real(wp),parameter :: b2521 = -0.00777788129242204164032546458607364309759347209626759111946150_wp real(wp),parameter :: b2522 = -0.00356084192402274913338827232697437364675240818791706587952939_wp real(wp),parameter :: b2523 = 4.79282506449930799649797749629840189457296934139359048988332_wp real(wp),parameter :: b2524 = 0.153725464873068577844576387402512082757034273069877432944621_wp real(wp),parameter :: b260 = 32.3038520871985442326994734440031535091364975047784630088983_wp real(wp),parameter :: b261 = 0.0_wp real(wp),parameter :: b262 = 0.0_wp real(wp),parameter :: b263 = 0.0_wp real(wp),parameter :: b264 = 0.0_wp real(wp),parameter :: b265 = -0.00317987696266205093901912847692712407988609169703103952205634_wp real(wp),parameter :: b266 = 0.806397714906192077260821711520379506393543111567419750119748_wp real(wp),parameter :: b267 = 0.0975983126412388979093522850684288851314672048003054550357187_wp real(wp),parameter :: b268 = 0.778575578158398909027512446452927238999763460594181964958853_wp real(wp),parameter :: b269 = 0.204890423831599428189499202098105603312029235081420653574829_wp real(wp),parameter :: b2610 = -1.56261579627468188307070943950527825211462892236424360892806_wp real(wp),parameter :: b2611 = 0.0_wp real(wp),parameter :: b2612 = 16.3429891882310570648504243973927174708753353504154550405647_wp real(wp),parameter :: b2613 = -154.544555293543621230730189631471036399316683669609116705323_wp real(wp),parameter :: b2614 = 1.56971088703334872692034283417621761466263593582497085955201_wp real(wp),parameter :: b2615 = 3.27685545087248131321429817269900731165522404974733504794135_wp real(wp),parameter :: b2616 = -0.0503489245193653176348040727199783626534081095691632396802451_wp real(wp),parameter :: b2617 = 153.321151858041665070593767885914694011224363102594556731397_wp real(wp),parameter :: b2618 = 7.17568186327720495846766484814784143567826308034865369443637_wp real(wp),parameter :: b2619 = -2.94036748675300481945917659896930989215320594380777597403592_wp real(wp),parameter :: b2620 = -0.0665845946076803144470749676022628870281920493197256887985612_wp real(wp),parameter :: b2621 = -0.0462346054990843661229248668562217261176966514016859284197145_wp real(wp),parameter :: b2622 = -0.0204198733585679401539388228617269778848579774821581777675337_wp real(wp),parameter :: b2623 = -53.3523106438735850515953441165998107974045090495791591218714_wp real(wp),parameter :: b2624 = -1.35548714715078654978732186705996404017554501614191325114947_wp real(wp),parameter :: b2625 = -1.57196275801232751882901735171459249177687219114442583461866_wp real(wp),parameter :: b270 = -16.6451467486341512872031294403931758764560371130818978459405_wp real(wp),parameter :: b271 = 0.0_wp real(wp),parameter :: b272 = 0.0_wp real(wp),parameter :: b273 = 0.0_wp real(wp),parameter :: b274 = 0.0_wp real(wp),parameter :: b275 = 0.00592232780324503308042990005798046524738389560444257136834990_wp real(wp),parameter :: b276 = 0.470326159963841112217224303205894113455362530746108825010848_wp real(wp),parameter :: b277 = 0.299688863848679000853981837096192399136831121671781279184194_wp real(wp),parameter :: b278 = -0.247656877593994914689992276329810825853958069263947095548189_wp real(wp),parameter :: b279 = 0.110895029771437682893999851839061714522445173600678718208625_wp real(wp),parameter :: b2710 = 0.0_wp real(wp),parameter :: b2711 = -0.491719043846229147070666628704194097678081907210673044988866_wp real(wp),parameter :: b2712 = -11.4743154427289496968389492564352536350842454130853175250727_wp real(wp),parameter :: b2713 = 80.2593166576230272541702485886484400152793366623589989106256_wp real(wp),parameter :: b2714 = -0.384132303980042847625312526759029103746926841342088219165648_wp real(wp),parameter :: b2715 = 7.28147667468107583471326950926136115767612581862877764249646_wp real(wp),parameter :: b2716 = -0.132699384612248379510571708176035274836827341616751884314074_wp real(wp),parameter :: b2717 = -81.0799832525730726674679289752255240006070716633632990308935_wp real(wp),parameter :: b2718 = -1.25037492835620639521768185656179119962253747492403205797494_wp real(wp),parameter :: b2719 = 2.59263594969543681023776379504377324994226447359296887778718_wp real(wp),parameter :: b2720 = -0.301440298346404539830163997260526875264431537275641495291993_wp real(wp),parameter :: b2721 = 0.221384460789832337451706451572773791695246839057318414301020_wp real(wp),parameter :: b2722 = 0.0827577274771892931955989870974693152996276435429809890551210_wp real(wp),parameter :: b2723 = 18.9960662040611520464672450037243263998175161412237156872211_wp real(wp),parameter :: b2724 = 0.269231946409639685623468015128334167460051910348912845121977_wp real(wp),parameter :: b2725 = 1.62674827447066537462989364929628933988125029284183680279020_wp real(wp),parameter :: b2726 = 0.491719043846229147070666628704194097678081907210673044988866_wp real(wp),parameter :: b280 = 0.0838479812409052664616968791372814085980533139224911131069335_wp real(wp),parameter :: b281 = 0.0_wp real(wp),parameter :: b282 = 0.0_wp real(wp),parameter :: b283 = 0.0_wp real(wp),parameter :: b284 = 0.0_wp real(wp),parameter :: b285 = -0.0117949367100973814319755056031295775367961960590736150777613_wp real(wp),parameter :: b286 = -0.247299020568812652339473838743194598325992840353340132697498_wp real(wp),parameter :: b287 = 0.0978080858367729012259313014081291665503740655476733940756599_wp real(wp),parameter :: b288 = 0.217590689243420631360008651767860318344168120024782176879989_wp real(wp),parameter :: b289 = 0.0_wp real(wp),parameter :: b2810 = 0.137585606763325224865659632196787746647447222975084865975440_wp real(wp),parameter :: b2811 = 0.0439870229715046685058790092341545026046103890294261359042581_wp real(wp),parameter :: b2812 = 0.0_wp real(wp),parameter :: b2813 = -0.513700813768193341957004456618630303738757363641964030086972_wp real(wp),parameter :: b2814 = 0.826355691151315508644211308399153458701423158616168576922372_wp real(wp),parameter :: b2815 = 25.7018139719811832625873882972519939511136556341960074626615_wp real(wp),parameter :: b2816 = 0.0_wp real(wp),parameter :: b2817 = 0.0_wp real(wp),parameter :: b2818 = 0.0_wp real(wp),parameter :: b2819 = 0.0_wp real(wp),parameter :: b2820 = 0.0_wp real(wp),parameter :: b2821 = 0.0_wp real(wp),parameter :: b2822 = 0.0_wp real(wp),parameter :: b2823 = -25.7018139719811832625873882972519939511136556341960074626615_wp real(wp),parameter :: b2824 = -0.826355691151315508644211308399153458701423158616168576922372_wp real(wp),parameter :: b2825 = 0.513700813768193341957004456618630303738757363641964030086972_wp real(wp),parameter :: b2826 = -0.0439870229715046685058790092341545026046103890294261359042581_wp real(wp),parameter :: b2827 = -0.137585606763325224865659632196787746647447222975084865975440_wp real(wp),parameter :: b290 = 0.124380526654094412881516420868799316268491466359671423163289_wp real(wp),parameter :: b291 = 0.0_wp real(wp),parameter :: b292 = 0.0_wp real(wp),parameter :: b293 = 0.0_wp real(wp),parameter :: b294 = 0.226120282197584301422238662979202901196752320742633143965145_wp real(wp),parameter :: b295 = 0.0137885887618080880607695837016477814530969417491493385363543_wp real(wp),parameter :: b296 = -0.0672210133996684449749399507414305856950086341525382182856200_wp real(wp),parameter :: b297 = 0.0_wp real(wp),parameter :: b298 = 0.0_wp real(wp),parameter :: b299 = -0.856238975085428354755349769879501772112121597411563802855067_wp real(wp),parameter :: b2910 = -1.96337522866858908928262850028093813988180440518267404553576_wp real(wp),parameter :: b2911 = -0.232332822724119401237246257308921847250108199230419994978218_wp real(wp),parameter :: b2912 = 0.0_wp real(wp),parameter :: b2913 = 4.30660719086453349461668936876562947772432562053478092626764_wp real(wp),parameter :: b2914 = -2.92722963249465482659787911202390446687687394950633612630592_wp real(wp),parameter :: b2915 = -82.3131666397858944454492334105458707735761966428138676971041_wp real(wp),parameter :: b2916 = 0.0_wp real(wp),parameter :: b2917 = 0.0_wp real(wp),parameter :: b2918 = 0.0_wp real(wp),parameter :: b2919 = 0.0_wp real(wp),parameter :: b2920 = 0.0_wp real(wp),parameter :: b2921 = 0.0_wp real(wp),parameter :: b2922 = 0.0_wp real(wp),parameter :: b2923 = 82.3131666397858944454492334105458707735761966428138676971041_wp real(wp),parameter :: b2924 = 2.92722963249465482659787911202390446687687394950633612630592_wp real(wp),parameter :: b2925 = -4.30660719086453349461668936876562947772432562053478092626764_wp real(wp),parameter :: b2926 = 0.232332822724119401237246257308921847250108199230419994978218_wp real(wp),parameter :: b2927 = 1.96337522866858908928262850028093813988180440518267404553576_wp real(wp),parameter :: b2928 = 0.856238975085428354755349769879501772112121597411563802855067_wp real(wp),parameter :: b300 = 0.103484561636679776672993546511910344499744798201971316606663_wp real(wp),parameter :: b301 = 0.0_wp real(wp),parameter :: b302 = 0.0_wp real(wp),parameter :: b303 = 0.122068887306407222589644082868962077139592714834162134741275_wp real(wp),parameter :: b304 = 0.482574490331246622475134780125688112865919023850168049679402_wp real(wp),parameter :: b305 = -0.0381409600015606999730886240005620205664113072478411477421970_wp real(wp),parameter :: b306 = 0.0_wp real(wp),parameter :: b307 = -0.550499525310802324138388507020508177411414311000037561712836_wp real(wp),parameter :: b308 = 0.0_wp real(wp),parameter :: b309 = -0.711915811585189227887648262043794387578291882406745570495765_wp real(wp),parameter :: b3010 = -0.584129605671551340432988730158480872095335329645227595707052_wp real(wp),parameter :: b3011 = 0.0_wp real(wp),parameter :: b3012 = 0.0_wp real(wp),parameter :: b3013 = 2.11046308125864932128717300046622750300375054278936987850718_wp real(wp),parameter :: b3014 = -0.0837494736739572135525742023001037992695260175335123517729291_wp real(wp),parameter :: b3015 = 5.10021499072320914075295969043344113107545060862804249161191_wp real(wp),parameter :: b3016 = 0.0_wp real(wp),parameter :: b3017 = 0.0_wp real(wp),parameter :: b3018 = 0.0_wp real(wp),parameter :: b3019 = 0.0_wp real(wp),parameter :: b3020 = 0.0_wp real(wp),parameter :: b3021 = 0.0_wp real(wp),parameter :: b3022 = 0.0_wp real(wp),parameter :: b3023 = -5.10021499072320914075295969043344113107545060862804249161191_wp real(wp),parameter :: b3024 = 0.0837494736739572135525742023001037992695260175335123517729291_wp real(wp),parameter :: b3025 = -2.11046308125864932128717300046622750300375054278936987850718_wp real(wp),parameter :: b3026 = 0.0_wp real(wp),parameter :: b3027 = 0.584129605671551340432988730158480872095335329645227595707052_wp real(wp),parameter :: b3028 = 0.711915811585189227887648262043794387578291882406745570495765_wp real(wp),parameter :: b3029 = 0.550499525310802324138388507020508177411414311000037561712836_wp real(wp),parameter :: b310 = 0.193333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: b311 = 0.0_wp real(wp),parameter :: b312 = 0.220000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b313 = -0.0800000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b314 = 0.0_wp real(wp),parameter :: b315 = 0.0_wp real(wp),parameter :: b316 = 0.109993425580724703919462404865068340845119058295846426463652_wp real(wp),parameter :: b317 = -0.254297048076270161384068506997153122141835626976703920846242_wp real(wp),parameter :: b318 = 0.0_wp real(wp),parameter :: b319 = 0.865570777116694254343770343821098281832847401233011859346737_wp real(wp),parameter :: b3110 = 3.32416449114093083106799552786572018336860092936986407160200_wp real(wp),parameter :: b3111 = 0.0_wp real(wp),parameter :: b3112 = 0.0_wp real(wp),parameter :: b3113 = -12.0102223315977933882352385148661841260301942633996815127277_wp real(wp),parameter :: b3114 = 0.476601466242493239430442776862061899602963782003580209476163_wp real(wp),parameter :: b3115 = -29.0243011221036390525802623213654099596251221332470910692353_wp real(wp),parameter :: b3116 = 0.0_wp real(wp),parameter :: b3117 = 0.0_wp real(wp),parameter :: b3118 = 0.0_wp real(wp),parameter :: b3119 = 0.0_wp real(wp),parameter :: b3120 = 0.0_wp real(wp),parameter :: b3121 = 0.0_wp real(wp),parameter :: b3122 = 0.0_wp real(wp),parameter :: b3123 = 29.0243011221036390525802623213654099596251221332470910692353_wp real(wp),parameter :: b3124 = -0.476601466242493239430442776862061899602963782003580209476163_wp real(wp),parameter :: b3125 = 12.0102223315977933882352385148661841260301942633996815127277_wp real(wp),parameter :: b3126 = 0.0_wp real(wp),parameter :: b3127 = -3.32416449114093083106799552786572018336860092936986407160200_wp real(wp),parameter :: b3128 = -0.865570777116694254343770343821098281832847401233011859346737_wp real(wp),parameter :: b3129 = 0.254297048076270161384068506997153122141835626976703920846242_wp real(wp),parameter :: b3130 = -0.109993425580724703919462404865068340845119058295846426463652_wp real(wp),parameter :: b320 = -0.833333333333333333333333333333333333333333333333333333333333_wp real(wp),parameter :: b321 = 1.38888888888888888888888888888888888888888888888888888888889_wp real(wp),parameter :: b322 = 0.0_wp real(wp),parameter :: b323 = 0.0_wp real(wp),parameter :: b324 = -0.750000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b325 = 0.0_wp real(wp),parameter :: b326 = -0.492529543718026304422682049114021320200214681580657784719074_wp real(wp),parameter :: b327 = 0.0_wp real(wp),parameter :: b328 = 0.0_wp real(wp),parameter :: b329 = 0.0_wp real(wp),parameter :: b3210 = 0.0_wp real(wp),parameter :: b3211 = 0.0_wp real(wp),parameter :: b3212 = 0.0_wp real(wp),parameter :: b3213 = 0.0_wp real(wp),parameter :: b3214 = 0.0_wp real(wp),parameter :: b3215 = 0.0_wp real(wp),parameter :: b3216 = 0.0_wp real(wp),parameter :: b3217 = 0.0_wp real(wp),parameter :: b3218 = 0.0_wp real(wp),parameter :: b3219 = 0.0_wp real(wp),parameter :: b3220 = 0.0_wp real(wp),parameter :: b3221 = 0.0_wp real(wp),parameter :: b3222 = 0.0_wp real(wp),parameter :: b3223 = 0.0_wp real(wp),parameter :: b3224 = 0.0_wp real(wp),parameter :: b3225 = 0.0_wp real(wp),parameter :: b3226 = 0.0_wp real(wp),parameter :: b3227 = 0.0_wp real(wp),parameter :: b3228 = 0.0_wp real(wp),parameter :: b3229 = 0.0_wp real(wp),parameter :: b3230 = 0.492529543718026304422682049114021320200214681580657784719074_wp real(wp),parameter :: b3231 = 0.750000000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b330 = 0.111111111111111111111111111111111111111111111111111111111111_wp real(wp),parameter :: b331 = 0.0_wp real(wp),parameter :: b332 = -0.222222222222222222222222222222222222222222222222222222222222_wp real(wp),parameter :: b333 = 0.0_wp real(wp),parameter :: b334 = 0.0_wp real(wp),parameter :: b335 = 0.0_wp real(wp),parameter :: b336 = 0.0_wp real(wp),parameter :: b337 = 0.0_wp real(wp),parameter :: b338 = 0.0_wp real(wp),parameter :: b339 = 0.0_wp real(wp),parameter :: b3310 = 0.0_wp real(wp),parameter :: b3311 = 0.0_wp real(wp),parameter :: b3312 = 0.0_wp real(wp),parameter :: b3313 = 0.0_wp real(wp),parameter :: b3314 = 0.0_wp real(wp),parameter :: b3315 = 0.0_wp real(wp),parameter :: b3316 = 0.0_wp real(wp),parameter :: b3317 = 0.0_wp real(wp),parameter :: b3318 = 0.0_wp real(wp),parameter :: b3319 = 0.0_wp real(wp),parameter :: b3320 = 0.0_wp real(wp),parameter :: b3321 = 0.0_wp real(wp),parameter :: b3322 = 0.0_wp real(wp),parameter :: b3323 = 0.0_wp real(wp),parameter :: b3324 = 0.0_wp real(wp),parameter :: b3325 = 0.0_wp real(wp),parameter :: b3326 = 0.0_wp real(wp),parameter :: b3327 = 0.0_wp real(wp),parameter :: b3328 = 0.0_wp real(wp),parameter :: b3329 = 0.0_wp real(wp),parameter :: b3330 = 0.0_wp real(wp),parameter :: b3331 = 0.0_wp real(wp),parameter :: b3332 = 0.222222222222222222222222222222222222222222222222222222222222_wp real(wp),parameter :: b340 = 0.285835140388971558796088842163836414852927537894596466840753_wp real(wp),parameter :: b341 = 0.291666666666666666666666666666666666666666666666666666666667_wp real(wp),parameter :: b342 = 0.218750000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b343 = 0.0_wp real(wp),parameter :: b344 = 0.164062500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b345 = 0.0_wp real(wp),parameter :: b346 = 0.218194354945556658327188241581352107093288824322187941141516_wp real(wp),parameter :: b347 = 0.180392898478697766863635221946775437719620053641849228562435_wp real(wp),parameter :: b348 = 0.0_wp real(wp),parameter :: b349 = 0.205713839404845018859120755122929542277570094982808905393991_wp real(wp),parameter :: b3410 = 0.242715791581770239970282927959446515762745971386670541948576_wp real(wp),parameter :: b3411 = 0.246465780813629305833609291181891407799228103869305705137021_wp real(wp),parameter :: b3412 = -3.44991940790890824979834154601622662060370460614931644223924_wp real(wp),parameter :: b3413 = 0.228875562160036081760729060738458584294220372552740218459295_wp real(wp),parameter :: b3414 = 0.283290599702151415321527419056733335978436595493855789831434_wp real(wp),parameter :: b3415 = 3.21085125837766640960131490544236787005557320332238705967955_wp real(wp),parameter :: b3416 = -0.223538777364845699920233756214162507964125230083674032084065_wp real(wp),parameter :: b3417 = -0.707121157204419073518727286207487212130091231955206160635271_wp real(wp),parameter :: b3418 = 3.21123345150287080408174729202856500893260034443022374267639_wp real(wp),parameter :: b3419 = 1.40954348309669766030414474301123175769045945573548986335553_wp real(wp),parameter :: b3420 = -0.151362053443742613121602276742518111090963026203676055891793_wp real(wp),parameter :: b3421 = 0.372350574527014276454724080214619984397121028202148298716575_wp real(wp),parameter :: b3422 = 0.252978746406361336722199907762141285915775728129414319261111_wp real(wp),parameter :: b3423 = -3.21085125837766640960131490544236787005557320332238705967955_wp real(wp),parameter :: b3424 = -0.283290599702151415321527419056733335978436595493855789831434_wp real(wp),parameter :: b3425 = -0.228875562160036081760729060738458584294220372552740218459295_wp real(wp),parameter :: b3426 = -0.246465780813629305833609291181891407799228103869305705137021_wp real(wp),parameter :: b3427 = -0.242715791581770239970282927959446515762745971386670541948576_wp real(wp),parameter :: b3428 = -0.205713839404845018859120755122929542277570094982808905393991_wp real(wp),parameter :: b3429 = -0.180392898478697766863635221946775437719620053641849228562435_wp real(wp),parameter :: b3430 = -0.218194354945556658327188241581352107093288824322187941141516_wp real(wp),parameter :: b3431 = -0.164062500000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b3432 = -0.218750000000000000000000000000000000000000000000000000000000_wp real(wp),parameter :: b3433 = -0.291666666666666666666666666666666666666666666666666666666667_wp real(wp),dimension(me%n) :: f0,f1,f2,f3,f4,f5,f6,f7,f8,f9,f10,f11,f12,& f13,f14,f15,f16,f17,f18,f19,f20,f21,f22,f23,f24,& f25,f26,f27,f28,f29,f30,f31,f32,f33,f34 if (h==zero) then xf = x terr = zero return end if call me%f(t+a0*h, x,f0) call me%f(t+a1*h, x+h*(b10*f0),f1) call me%f(t+a2*h, x+h*(b20*f0 + b21 *f1),f2) call me%f(t+a3*h, x+h*(b30*f0 + b31 *f1 + b32*f2),f3) call me%f(t+a4*h, x+h*(b40*f0 + b41 *f1 + b42*f2 + b43*f3),f4) call me%f(t+a5*h, x+h*(b50*f0 + b51 *f1 + b52*f2 + b53*f3 + & b54*f4),f5) call me%f(t+a6*h, x+h*(b60*f0 + b61 *f1 + b62*f2 + b63*f3 + & b64*f4 + b65*f5),f6) call me%f(t+a7*h, x+h*(b70*f0 + b71 *f1 + b72*f2 + b73*f3 + & b74*f4 + b75*f5 + b76*f6),f7) call me%f(t+a8*h, x+h*(b80*f0 + b81 *f1 + b82*f2 + b83*f3 + & b84*f4 + b85*f5 + b86*f6 + b87*f7),f8) call me%f(t+a9*h, x+h*(b90*f0 + b91 *f1 + b92*f2 + b93*f3 + & b94*f4 + b95*f5 + b96*f6 + b97*f7 + & b98*f8),f9) call me%f(t+a10*h, x+h*(b100*f0 + b101*f1 + b102*f2 + b103*f3 + & b104*f4 + b105*f5 + b106*f6 + b107*f7 + & b108*f8 + b109*f9),f10) call me%f(t+a11*h, x+h*(b110*f0 + b111*f1 + b112*f2 + b113*f3 + & b114*f4 + b115*f5 + b116*f6 + b117*f7 + & b118*f8 + b119*f9 + b1110*f10),f11) call me%f(t+a12*h, x+h*(b120*f0 + b121*f1 + b122*f2 + b123*f3 + & b124*f4 + b125*f5 + b126*f6 + b127*f7 + & b128*f8 + b129*f9 + b1210*f10 + b1211*f11),f12) call me%f(t+a13*h, x+h*(b130*f0 + b131*f1 + b132*f2 + b133*f3 + & b134*f4 + b135*f5 + b136*f6 + b137*f7 + & b138*f8 + b139*f9 + b1310*f10 + b1311*f11 + & b1312*f12),f13) call me%f(t+a14*h, x+h*(b140*f0 + b141*f1 + b142*f2 + b143*f3 + & b144*f4 + b145*f5 + b146*f6 + b147*f7 + & b148*f8 + b149*f9 + b1410*f10 + b1411*f11 + & b1412*f12 + b1413*f13),f14) call me%f(t+a15*h, x+h*(b150*f0 + b151*f1 + b152*f2 + b153*f3 + & b154*f4 + b155*f5 + b156*f6 + b157*f7 + & b158*f8 + b159*f9 + b1510*f10 + b1511*f11 + & b1512*f12 + b1513*f13 + b1514*f14),f15) call me%f(t+a16*h, x+h*(b160*f0 + b161*f1 + b162*f2 + b163*f3 + & b164*f4 + b165*f5 + b166*f6 + b167*f7 + & b168*f8 + b169*f9 + b1610*f10 + b1611*f11 + & b1612*f12 + b1613*f13 + b1614*f14 + b1615*f15),f16) call me%f(t+a17*h, x+h*(b170*f0 + b171*f1 + b172*f2 + b173*f3 + & b174*f4 + b175*f5 + b176*f6 + b177*f7 + & b178*f8 + b179*f9 + b1710*f10 + b1711*f11 + & b1712*f12 + b1713*f13 + b1714*f14 + b1715*f15 + & b1716*f16),f17) call me%f(t+a18*h, x+h*(b180*f0 + b181*f1 + b182*f2 + b183*f3 + & b184*f4 + b185*f5 + b186*f6 + b187*f7 + & b188*f8 + b189*f9 + b1810*f10 + b1811*f11 + & b1812*f12 + b1813*f13 + b1814*f14 + b1815*f15 + & b1816*f16 + b1817*f17),f18) call me%f(t+a19*h, x+h*(b190*f0 + b191*f1 + b192*f2 + b193*f3 + & b194*f4 + b195*f5 + b196*f6 + b197*f7 + & b198*f8 + b199*f9 + b1910*f10 + b1911*f11 + & b1912*f12 + b1913*f13 + b1914*f14 + b1915*f15 + & b1916*f16 + b1917*f17 + b1918*f18),f19) call me%f(t+a20*h, x+h*(b200*f0 + b201*f1 + b202*f2 + b203*f3 + & b204*f4 + b205*f5 + b206*f6 + b207*f7 + & b208*f8 + b209*f9 + b2010*f10 + b2011*f11 + & b2012*f12 + b2013*f13 + b2014*f14 + b2015*f15 + & b2016*f16 + b2017*f17 + b2018*f18 + b2019*f19),f20) call me%f(t+a21*h, x+h*(b210*f0 + b211*f1 + b212*f2 + b213*f3 + & b214*f4 + b215*f5 + b216*f6 + b217*f7 + & b218*f8 + b219*f9 + b2110*f10 + b2111*f11 + & b2112*f12 + b2113*f13 + b2114*f14 + b2115*f15 + & b2116*f16 + b2117*f17 + b2118*f18 + b2119*f19 + & b2120*f20),f21) call me%f(t+a22*h, x+h*(b220*f0 + b221*f1 + b222*f2 + b223*f3 + & b224*f4 + b225*f5 + b226*f6 + b227*f7 + & b228*f8 + b229*f9 + b2210*f10 + b2211*f11 + & b2212*f12 + b2213*f13 + b2214*f14 + b2215*f15 + & b2216*f16 + b2217*f17 + b2218*f18 + b2219*f19 + & b2220*f20 + b2221*f21),f22) call me%f(t+a23*h, x+h*(b230*f0 + b231*f1 + b232*f2 + b233*f3 + & b234*f4 + b235*f5 + b236*f6 + b237*f7 + & b238*f8 + b239*f9 + b2310*f10 + b2311*f11 + & b2312*f12 + b2313*f13 + b2314*f14 + b2315*f15 + & b2316*f16 + b2317*f17 + b2318*f18 + b2319*f19 + & b2320*f20 + b2321*f21 + b2322*f22),f23) call me%f(t+a24*h, x+h*(b240*f0 + b241*f1 + b242*f2 + b243*f3 + & b244*f4 + b245*f5 + b246*f6 + b247*f7 + & b248*f8 + b249*f9 + b2410*f10 + b2411*f11 + & b2412*f12 + b2413*f13 + b2414*f14 + b2415*f15 + & b2416*f16 + b2417*f17 + b2418*f18 + b2419*f19 + & b2420*f20 + b2421*f21 + b2422*f22 + b2423*f23),f24) call me%f(t+a25*h, x+h*(b250*f0 + b251*f1 + b252*f2 + b253*f3 + & b254*f4 + b255*f5 + b256*f6 + b257*f7 + & b258*f8 + b259*f9 + b2510*f10 + b2511*f11 + & b2512*f12 + b2513*f13 + b2514*f14 + b2515*f15 + & b2516*f16 + b2517*f17 + b2518*f18 + b2519*f19 + & b2520*f20 + b2521*f21 + b2522*f22 + b2523*f23 + & b2524*f24),f25) call me%f(t+a26*h, x+h*(b260*f0 + b261*f1 + b262*f2 + b263*f3 + & b264*f4 + b265*f5 + b266*f6 + b267*f7 + & b268*f8 + b269*f9 + b2610*f10 + b2611*f11 + & b2612*f12 + b2613*f13 + b2614*f14 + b2615*f15 + & b2616*f16 + b2617*f17 + b2618*f18 + b2619*f19 + & b2620*f20 + b2621*f21 + b2622*f22 + b2623*f23 + & b2624*f24 + b2625*f25),f26) call me%f(t+a27*h, x+h*(b270*f0 + b271*f1 + b272*f2 + b273*f3 + & b274*f4 + b275*f5 + b276*f6 + b277*f7 + & b278*f8 + b279*f9 + b2710*f10 + b2711*f11 + & b2712*f12 + b2713*f13 + b2714*f14 + b2715*f15 + & b2716*f16 + b2717*f17 + b2718*f18 + b2719*f19 + & b2720*f20 + b2721*f21 + b2722*f22 + b2723*f23 + & b2724*f24 + b2725*f25 + b2726*f26),f27) call me%f(t+a28*h, x+h*(b280*f0 + b281*f1 + b282*f2 + b283*f3 + & b284*f4 + b285*f5 + b286*f6 + b287*f7 + & b288*f8 + b289*f9 + b2810*f10 + b2811*f11 + & b2812*f12 + b2813*f13 + b2814*f14 + b2815*f15 + & b2816*f16 + b2817*f17 + b2818*f18 + b2819*f19 + & b2820*f20 + b2821*f21 + b2822*f22 + b2823*f23 + & b2824*f24 + b2825*f25 + b2826*f26 + b2827*f27),f28) call me%f(t+a29*h, x+h*(b290*f0 + b291*f1 + b292*f2 + b293*f3 + & b294*f4 + b295*f5 + b296*f6 + b297*f7 + & b298*f8 + b299*f9 + b2910*f10 + b2911*f11 + & b2912*f12 + b2913*f13 + b2914*f14 + b2915*f15 + & b2916*f16 + b2917*f17 + b2918*f18 + b2919*f19 + & b2920*f20 + b2921*f21 + b2922*f22 + b2923*f23 + & b2924*f24 + b2925*f25 + b2926*f26 + b2927*f27 + & b2928*f28),f29) call me%f(t+a30*h, x+h*(b300*f0 + b301*f1 + b302*f2 + b303*f3 + & b304*f4 + b305*f5 + b306*f6 + b307*f7 + & b308*f8 + b309*f9 + b3010*f10 + b3011*f11 + & b3012*f12 + b3013*f13 + b3014*f14 + b3015*f15 + & b3016*f16 + b3017*f17 + b3018*f18 + b3019*f19 + & b3020*f20 + b3021*f21 + b3022*f22 + b3023*f23 + & b3024*f24 + b3025*f25 + b3026*f26 + b3027*f27 + & b3028*f28 + b3029*f29),f30) call me%f(t+a31*h, x+h*(b310*f0 + b311*f1 + b312*f2 + b313*f3 + & b314*f4 + b315*f5 + b316*f6 + b317*f7 + & b318*f8 + b319*f9 + b3110*f10 + b3111*f11 + & b3112*f12 + b3113*f13 + b3114*f14 + b3115*f15 + & b3116*f16 + b3117*f17 + b3118*f18 + b3119*f19 + & b3120*f20 + b3121*f21 + b3122*f22 + b3123*f23 + & b3124*f24 + b3125*f25 + b3126*f26 + b3127*f27 + & b3128*f28 + b3129*f29 + b3130*f30),f31) call me%f(t+a32*h, x+h*(b320*f0 + b321*f1 + b322*f2 + b323*f3 + & b324*f4 + b325*f5 + b326*f6 + b327*f7 + & b328*f8 + b329*f9 + b3210*f10 + b3211*f11 + & b3212*f12 + b3213*f13 + b3214*f14 + b3215*f15 + & b3216*f16 + b3217*f17 + b3218*f18 + b3219*f19 + & b3220*f20 + b3221*f21 + b3222*f22 + b3223*f23 + & b3224*f24 + b3225*f25 + b3226*f26 + b3227*f27 + & b3228*f28 + b3229*f29 + b3230*f30 + b3231*f31),f32) call me%f(t+a33*h, x+h*(b330*f0 + b331*f1 + b332*f2 + b333*f3 + & b334*f4 + b335*f5 + b336*f6 + b337*f7 + & b338*f8 + b339*f9 + b3310*f10 + b3311*f11 + & b3312*f12 + b3313*f13 + b3314*f14 + b3315*f15 + & b3316*f16 + b3317*f17 + b3318*f18 + b3319*f19 + & b3320*f20 + b3321*f21 + b3322*f22 + b3323*f23 + & b3324*f24 + b3325*f25 + b3326*f26 + b3327*f27 + & b3328*f28 + b3329*f29 + b3330*f30 + b3331*f31 + & b3332*f32),f33) call me%f(t+a34*h, x+h*(b340*f0 + b341*f1 + b342*f2 + b343*f3 + & b344*f4 + b345*f5 + b346*f6 + b347*f7 + & b348*f8 + b349*f9 + b3410*f10 + b3411*f11 + & b3412*f12 + b3413*f13 + b3414*f14 + b3415*f15 + & b3416*f16 + b3417*f17 + b3418*f18 + b3419*f19 + & b3420*f20 + b3421*f21 + b3422*f22 + b3423*f23 + & b3424*f24 + b3425*f25 + b3426*f26 + b3427*f27 + & b3428*f28 + b3429*f29 + b3430*f30 + b3431*f31 + & b3432*f32 + b3433*f33),f34) xf = x+h*( c0*f0 + & c1*f1 + & c2*f2 + & c3*f3 + & c4*f4 + & c5*f5 + & c6*f6 + & c7*f7 + & c8*f8 + & c9*f9 + & c10*f10 + & c11*f11 + & c12*f12 + & c13*f13 + & c14*f14 + & c15*f15 + & c16*f16 + & c17*f17 + & c18*f18 + & c19*f19 + & c20*f20 + & c21*f21 + & c22*f22 + & c23*f23 + & c24*f24 + & c25*f25 + & c26*f26 + & c27*f27 + & c28*f28 + & c29*f29 + & c30*f30 + & c31*f31 + & c32*f32 + & c33*f33 + & c34*f34 ) terr = (1.0_wp/1000.0_wp)*h*(f1-f33) end subroutine rkf1412