From fb5fddf002e842058fb2504bc84f09fbb63a8f13 Mon Sep 17 00:00:00 2001 From: diener Date: Thu, 26 May 2005 23:26:04 +0000 Subject: Dissipation operator for the optimized 6-5 restricted full norm operator. Don't use yet. Still being tested. git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@35 f69c4107-0314-4c4f-9ad4-17e986b73f4a --- src/Derivatives_6_5_min_err_coeff.F90 | 134 +++--- src/Dissipation_6_5_min_err_coeff.F90 | 869 ++++++++++++++++++++++++++++++++++ src/dissipation.c | 8 +- src/make.code.defn | 1 + 4 files changed, 941 insertions(+), 71 deletions(-) create mode 100644 src/Dissipation_6_5_min_err_coeff.F90 diff --git a/src/Derivatives_6_5_min_err_coeff.F90 b/src/Derivatives_6_5_min_err_coeff.F90 index 5bc1a2f..cad827a 100644 --- a/src/Derivatives_6_5_min_err_coeff.F90 +++ b/src/Derivatives_6_5_min_err_coeff.F90 @@ -30,76 +30,76 @@ subroutine deriv_gf_6_5_opt ( var, ni, nj, nk, dir, bb, gsize, delta, dvar ) if ( first ) then a(1) = 3.0_wp/4.0_wp; a(2) = -3.0_wp/20.0_wp; a(3) = 1.0_wp/60.0_wp - q(1,1) = -2.465354921105726344726324628220306540433_wp - q(2,1) = 6.092129526634358068357947769321839242595_wp - q(3,1) = -7.730323816585895170894869423304598106488_wp - q(4,1) = 6.973765088781193561193159231072797475318_wp - q(5,1) = -3.980323816585895170894869423304598106488_wp - q(6,1) = 1.292129526634358068357947769321839242595_wp - q(7,1) = -0.1820215877723930113929912948869732070992_wp + q(1,1) = -2.465354921110524023660777656111276003457_wp + q(2,1) = 6.092129526663144141964665936667656020742_wp + q(3,1) = -7.730323816657860354911664841669140051855_wp + q(4,1) = 6.973765088877147139882219788892186735807_wp + q(5,1) = -3.980323816657860354911664841669140051855_wp + q(6,1) = 1.292129526663144141964665936667656020742_wp + q(7,1) = -0.1820215877771906903274443227779426701237_wp q(8,1) = 0_wp q(9,1) = 0_wp q(10,1) = 0_wp - q(1,2) = -0.2234725650933373070210703921399196916276_wp - q(2,2) = -0.9329308120213103664703924027998434010527_wp - q(3,2) = 1.586820596323380726926380762297803979729_wp - q(4,2) = -0.3647002337446073060093577955438077228163_wp - q(5,2) = -0.2666957787028035557343386019032589667016_wp - q(6,2) = 0.3112949049503615973044042012025239373575_wp - q(7,2) = -0.1404504214992091753318473921133396942483_wp - q(8,2) = 0.03488568515844113812839530849342152163114_wp - q(9,2) = -0.004964021892185918282918047104381964571176_wp - q(10,2) = 0.0002126465212701664907443596108020023005085_wp - q(1,3) = 0.1582216736972761112130174221319685578703_wp - q(2,3) = -1.137049297945701138153781830545857828001_wp - q(3,3) = 1.212364522772435724261769516755293203384_wp - q(4,3) = -0.9562288727032211388608555534303520293292_wp - q(5,3) = 1.066548057100650347439560790660954981316_wp - q(6,3) = -0.3478788549812286788557041220360802499887_wp - q(7,3) = -0.03133923299767926953262771814410417142710_wp - q(8,3) = 0.04098845957876514538943748816651847063573_wp - q(9,3) = -0.005963188640428636310134738586972978148865_wp - q(10,3) = 0.0003367341191315334093187450286320436894509_wp - q(1,4) = 0.02915734643197329316878786487755959617414_wp - q(2,4) = -0.1169665090568385435258745132657230868726_wp - q(3,4) = -0.1112219090396860235587194333712503100976_wp - q(4,4) = -0.7924486264107931688501411569472360917480_wp - q(5,4) = 1.266650705052489006015656810787656809460_wp - q(6,4) = -0.2899273291614735652270415777078748545257_wp - q(7,4) = 0.002515684289746880247060269725567464781645_wp - q(8,4) = 0.01329713961045538259784924004687749664184_wp - q(9,4) = -0.001124464397192287407208734867324966241272_wp - q(10,4) = 0.00006796268131902653963123072174794242752627_wp - q(1,5) = -0.04582150000219473783416432201491386259525_wp - q(2,5) = 0.2240986548763910403084011643939472909730_wp - q(3,5) = -0.3246718492669968647456547557775778272008_wp - q(4,5) = -0.3929792922479188555634359057335494798829_wp - q(5,5) = 0.1166355819597526695402713402151963096924_wp - q(6,5) = 0.3449626905274072306980547367786542540315_wp - q(7,5) = 0.1430419813692381466377067274780162668130_wp - q(8,5) = -0.07764802500390609299919055841918594590836_wp - q(9,5) = 0.01332439335682098023215359404071263206533_wp - q(10,5) = -0.0009426355685935162741420209612996379880357_wp - q(1,6) = 0.003172814427534452820786948997021353279521_wp - q(2,6) = 0.00001061461750788952405498488693810144188215_wp - q(3,6) = -0.08747763621208721666864461344012255496503_wp - q(4,6) = 0.3975827328151903159828221701608106653957_wp - q(5,6) = -1.148835072887093214621630866015514532074_wp - q(6,6) = 0.3583006652024372613320453990260111557237_wp - q(7,6) = 0.5647665153537601810077289284875161241444_wp - q(8,6) = -0.09698196886045636272181593036009187114508_wp - q(9,6) = 0.008843905090944639404474662278460308660277_wp - q(10,6) = 0.0006174304522620539401783159789712495382110_wp - q(1,7) = -0.008639107520047315872218762260033429623639_wp - q(2,7) = 0.04722773941989840822874396849234369996233_wp - q(3,7) = -0.1008747534523833098882185784719780553451_wp - q(4,7) = 0.08043834912116582462744784937371906752091_wp - q(5,7) = 0.1295138677847054989778048608621896918103_wp - q(6,7) = -0.7909424167745257943256055480973493995702_wp - q(7,7) = 0.03807866849747788018644730375057882400112_wp - q(8,7) = 0.7367055699547762148521594995917960242549_wp - q(9,7) = -0.1480235854664835893774029568937961320792_wp - q(10,7) = 0.01651566843541618259084236365252970906860_wp + q(1,2) = -0.2234725650784319828746535134412736890421_wp + q(2,2) = -0.9329308121107134563129925525068570679651_wp + q(3,2) = 1.586820596545839371759081303802027231274_wp + q(4,2) = -0.3647002340377160216914505558624668821400_wp + q(5,2) = -0.2666957784872806143914117440166232718819_wp + q(6,2) = 0.3112949048634705032101261273629794071371_wp + q(7,2) = -0.1404504214762266650000768489896480092493_wp + q(8,2) = 0.03488568514730479833596013512958238764128_wp + q(9,2) = -0.004964021886392518344179263072091597647654_wp + q(10,2) = 0.0002126465201465853095969115943714918742904_wp + q(1,3) = 0.1582216737061633151406179477554921935333_wp + q(2,3) = -1.137049298003377811733609086574457439398_wp + q(3,3) = 1.212364522932578587741649981040340946798_wp + q(4,3) = -0.9562288729513894906148167047868730813830_wp + q(5,3) = 1.066548057336766350478498057851678826640_wp + q(6,3) = -0.3478788551267041838265477441805600110467_wp + q(7,3) = -0.03133923293520187620333693909408071632123_wp + q(8,3) = 0.04098845955755862691072597869183962277781_wp + q(9,3) = -0.005963188634687155197078928402509551508436_wp + q(10,3) = 0.0003367341182936373038974376991292099082999_wp + q(1,4) = 0.02915734641890708196910927068736798144670_wp + q(2,4) = -0.1169665089768926152768236581512624861308_wp + q(3,4) = -0.1112219092451476301503253995474190870412_wp + q(4,4) = -0.7924486261248032107393766820001361351677_wp + q(5,4) = 1.266650704820613624987450232358951199911_wp + q(6,4) = -0.2899273290506621673153239836530375587273_wp + q(7,4) = 0.002515684257201926199329020583484434062150_wp + q(8,4) = 0.01329713961871764653006682056620518602804_wp + q(9,4) = -0.001124464399630667352932212208930962568134_wp + q(10,4) = 0.00006796268169601114882659136477742818715059_wp + q(1,5) = -0.04582150000326981674750984653096293434777_wp + q(2,5) = 0.2240986548857151482718685516611524323427_wp + q(3,5) = -0.3246718493011818141660859125588209338018_wp + q(4,5) = -0.3929792921782506986152017485694441380503_wp + q(5,5) = 0.1166355818729375628072830916953646214341_wp + q(6,5) = 0.3449626905957060254933930895775644438105_wp + q(7,5) = 0.1430419813354607083034935179267283951745_wp + q(8,5) = -0.07764802499372607792980458731991885121073_wp + q(9,5) = 0.01332439335504217034559288889042994978834_wp + q(10,5) = -0.0009426355684332077630290447720929851395193_wp + q(1,6) = 0.003172814452954821196677290327889903944225_wp + q(2,6) = 0.00001061446045061551877105554145609103530766_wp + q(3,6) = -0.08747763580209736614983637747947172321794_wp + q(4,6) = 0.3975827322299876034907453299884380895682_wp + q(5,6) = -1.148835072393422871630425744497391344782_wp + q(6,6) = 0.3583006649535242306065761818925080902380_wp + q(7,6) = 0.5647665154270147564019144982190032455071_wp + q(8,6) = -0.09698196887272109736153117076061707705561_wp + q(9,6) = 0.008843905091972988427261446924164441884143_wp + q(10,6) = 0.0006174304523363194998474898440202828786385_wp + q(1,7) = -0.008639107540858839028043929986084287776394_wp + q(2,7) = 0.04722773954485212324714352753530343274219_wp + q(3,7) = -0.1008747537650261142294540111407681552350_wp + q(4,7) = 0.08043834953845218736895768965086958762389_wp + q(5,7) = 0.1295138674713300902982857323205417604553_wp + q(6,7) = -0.7909424166489541737614153656634872155367_wp + q(7,7) = 0.03807866847647628589685997987877954466259_wp + q(8,7) = 0.7367055699548196242687865288427927434250_wp + q(9,7) = -0.1480235854665196220062411065981933720158_wp + q(10,7) = 0.01651566843542843794512095516024596165494_wp first = .false. end if diff --git a/src/Dissipation_6_5_min_err_coeff.F90 b/src/Dissipation_6_5_min_err_coeff.F90 new file mode 100644 index 0000000..a88d3ad --- /dev/null +++ b/src/Dissipation_6_5_min_err_coeff.F90 @@ -0,0 +1,869 @@ +! $Header$ + +#include "cctk.h" +#include "cctk_Parameters.h" +#include "cctk_Functions.h" + + +subroutine dissipation_6_5_opt (var, lsh, gsh, lbnd, bb, gsize, & + delta, epsdisl, dfl, rhs) + + implicit none + + DECLARE_CCTK_PARAMETERS + DECLARE_CCTK_FUNCTIONS + + CCTK_INT, dimension(3), intent(in) :: lsh, gsh, lbnd + CCTK_INT :: ni, nj, nk + CCTK_REAL, dimension(lsh(1),lsh(2),lsh(3)), intent(in) :: var + CCTK_REAL, dimension(lsh(1),lsh(2),lsh(3)), intent(inout) :: rhs + CCTK_INT, dimension(6), intent(in) :: bb + CCTK_INT, dimension(3), intent(in) :: gsize + CCTK_REAL, dimension(3), intent(in) :: delta + CCTK_REAL, intent(in) :: epsdisl, dfl + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + integer :: i, j, k, center + CCTK_REAL, dimension(:,:), allocatable :: a, d, b, h, tmp + CCTK_REAL :: idel + + CCTK_INT :: il, ir, jl, jr, kl, kr + + ni = lsh(1); nj = lsh(2); nk = lsh(3) + + allocate ( a(ni,ni), d(ni,ni), b(ni,ni), h(ni,ni), tmp(ni,ni) ) + a = zero; d = zero; b = zero; h = zero; tmp = zero + + call set_bmatrix ( b, bb(1:2), lsh(1), gsh(1), lbnd(1), delta(1), dfl ) + + call set_hmatrix ( h, bb(1:2) ) + + if ( any ( bb(1:2) == 0 ) ) then + + call set_dmatrix ( d, bb(1:2) ) + + a = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) + + else + + if ( lsh(1) /= gsh(1) ) then + call CCTK_WARN(0, 'Internal error. Bounding box inconsistent with grid' ) + end if + + center = gsh(1) / 2 + ir = center + gsize(1) + + call set_dmatrix_half ( d(1:ir,1:ir), 1 ) + + tmp(1:ir,1:ir) = -transpose ( & + matmul ( h(1:ir,1:ir), & + matmul ( & + transpose ( d(1:ir,1:ir) ), & + matmul ( b(1:ir,1:ir), & + d(1:ir,1:ir) ) ) ) ) + + a(1:ir,1:center) = tmp(1:ir,1:center) + + il = center + 1 - gsize(1) + d = zero + call set_dmatrix_half ( d(il:ni,il:ni), 2 ) + tmp(il:ni,il:ni) = -transpose ( & + matmul ( h(il:ni,il:ni), & + matmul ( & + transpose ( d(il:ni,il:ni) ), & + matmul ( b(il:ni,il:ni), & + d(il:ni,il:ni) ) ) ) ) + + a(il:ni,center+1:ni) = tmp(il:ni,center+1:ni) + + end if + + idel = epsdisl / 64.0_wp + + if ( bb(1) == 0 ) then + + il = 1 + gsize(1) + + else + + rhs(1,:,:) = rhs(1,:,:) + & + ( a(1,1) * var(1,:,:) + a(2,1) * var(2,:,:) + & + a(3,1) * var(3,:,:) + a(4,1) * var(4,:,:) ) * idel + rhs(2,:,:) = rhs(2,:,:) + & + ( a(1,2) * var(1,:,:) + a(2,2) * var(2,:,:) + & + a(3,2) * var(3,:,:) + a(4,2) * var(4,:,:) + & + a(5,2) * var(5,:,:) + a(6,2) * var(6,:,:) + & + a(7,2) * var(7,:,:) + a(8,2) * var(8,:,:) + & + a(9,2) * var(9,:,:) + a(10,2) * var(10,:,:) ) * idel + rhs(3,:,:) = rhs(3,:,:) + & + ( a(1,3) * var(1,:,:) + a(2,3) * var(2,:,:) + & + a(3,3) * var(3,:,:) + a(4,3) * var(4,:,:) + & + a(5,3) * var(5,:,:) + a(6,3) * var(6,:,:) + & + a(7,3) * var(7,:,:) + a(8,3) * var(8,:,:) + & + a(9,3) * var(9,:,:) + a(10,3) * var(10,:,:) ) * idel + rhs(4,:,:) = rhs(4,:,:) + & + ( a(1,4) * var(1,:,:) + a(2,4) * var(2,:,:) + & + a(3,4) * var(3,:,:) + a(4,4) * var(4,:,:) + & + a(5,4) * var(5,:,:) + a(6,4) * var(6,:,:) + & + a(7,4) * var(7,:,:) + a(8,4) * var(8,:,:) + & + a(9,4) * var(9,:,:) + a(10,4) * var(10,:,:) ) * idel + rhs(5,:,:) = rhs(5,:,:) + & + ( a(1,5) * var(1,:,:) + a(2,5) * var(2,:,:) + & + a(3,5) * var(3,:,:) + a(4,5) * var(4,:,:) + & + a(5,5) * var(5,:,:) + a(6,5) * var(6,:,:) + & + a(7,5) * var(7,:,:) + a(8,5) * var(8,:,:) + & + a(9,5) * var(9,:,:) + a(10,5) * var(10,:,:) ) * idel + rhs(6,:,:) = rhs(6,:,:) + & + ( a(1,6) * var(1,:,:) + a(2,6) * var(2,:,:) + & + a(3,6) * var(3,:,:) + a(4,6) * var(4,:,:) + & + a(5,6) * var(5,:,:) + a(6,6) * var(6,:,:) + & + a(7,6) * var(7,:,:) + a(8,6) * var(8,:,:) + & + a(9,6) * var(9,:,:) + a(10,6) * var(10,:,:) ) * idel + rhs(7,:,:) = rhs(7,:,:) + & + ( a(1,7) * var(1,:,:) + a(2,7) * var(2,:,:) + & + a(3,7) * var(3,:,:) + a(4,7) * var(4,:,:) + & + a(5,7) * var(5,:,:) + a(6,7) * var(6,:,:) + & + a(7,7) * var(7,:,:) + a(8,7) * var(8,:,:) + & + a(9,7) * var(9,:,:) + a(10,7) * var(10,:,:) ) * idel + + il = 8 + + end if + + if ( bb(2) == 0 ) then + + ir = ni - gsize(1) + + else + + rhs(ni-6,:,:) = rhs(ni-6,:,:) + & + ( a(ni-9,ni-6) * var(ni-9,:,:) + & + a(ni-8,ni-6) * var(ni-8,:,:) + & + a(ni-7,ni-6) * var(ni-7,:,:) + & + a(ni-6,ni-6) * var(ni-6,:,:) + & + a(ni-5,ni-6) * var(ni-5,:,:) + & + a(ni-4,ni-6) * var(ni-4,:,:) + & + a(ni-3,ni-6) * var(ni-3,:,:) + & + a(ni-2,ni-6) * var(ni-2,:,:) + & + a(ni-1,ni-6) * var(ni-1,:,:) + & + a(ni,ni-6) * var(ni,:,:) ) * idel + rhs(ni-5,:,:) = rhs(ni-5,:,:) + & + ( a(ni-9,ni-5) * var(ni-9,:,:) + & + a(ni-8,ni-5) * var(ni-8,:,:) + & + a(ni-7,ni-5) * var(ni-7,:,:) + & + a(ni-6,ni-5) * var(ni-6,:,:) + & + a(ni-5,ni-5) * var(ni-5,:,:) + & + a(ni-4,ni-5) * var(ni-4,:,:) + & + a(ni-3,ni-5) * var(ni-3,:,:) + & + a(ni-2,ni-5) * var(ni-2,:,:) + & + a(ni-1,ni-5) * var(ni-1,:,:) + & + a(ni,ni-5) * var(ni,:,:) ) * idel + rhs(ni-4,:,:) = rhs(ni-4,:,:) + & + ( a(ni-9,ni-4) * var(ni-9,:,:) + & + a(ni-8,ni-4) * var(ni-8,:,:) + & + a(ni-7,ni-4) * var(ni-7,:,:) + & + a(ni-6,ni-4) * var(ni-6,:,:) + & + a(ni-5,ni-4) * var(ni-5,:,:) + & + a(ni-4,ni-4) * var(ni-4,:,:) + & + a(ni-3,ni-4) * var(ni-3,:,:) + & + a(ni-2,ni-4) * var(ni-2,:,:) + & + a(ni-1,ni-4) * var(ni-1,:,:) + & + a(ni,ni-4) * var(ni,:,:) ) * idel + rhs(ni-3,:,:) = rhs(ni-3,:,:) + & + ( a(ni-9,ni-3) * var(ni-9,:,:) + & + a(ni-8,ni-3) * var(ni-8,:,:) + & + a(ni-7,ni-3) * var(ni-7,:,:) + & + a(ni-6,ni-3) * var(ni-6,:,:) + & + a(ni-5,ni-3) * var(ni-5,:,:) + & + a(ni-4,ni-3) * var(ni-4,:,:) + & + a(ni-3,ni-3) * var(ni-3,:,:) + & + a(ni-2,ni-3) * var(ni-2,:,:) + & + a(ni-1,ni-3) * var(ni-1,:,:) + & + a(ni,ni-3) * var(ni,:,:) ) * idel + rhs(ni-2,:,:) = rhs(ni-2,:,:) + & + ( a(ni-9,ni-2) * var(ni-9,:,:) + & + a(ni-8,ni-2) * var(ni-8,:,:) + & + a(ni-7,ni-2) * var(ni-7,:,:) + & + a(ni-6,ni-2) * var(ni-6,:,:) + & + a(ni-5,ni-2) * var(ni-5,:,:) + & + a(ni-4,ni-2) * var(ni-4,:,:) + & + a(ni-3,ni-2) * var(ni-3,:,:) + & + a(ni-2,ni-2) * var(ni-2,:,:) + & + a(ni-1,ni-2) * var(ni-1,:,:) + & + a(ni,ni-2) * var(ni,:,:) ) * idel + rhs(ni-1,:,:) = rhs(ni-1,:,:) + & + ( a(ni-9,ni-1) * var(ni-9,:,:) + & + a(ni-8,ni-1) * var(ni-8,:,:) + & + a(ni-7,ni-1) * var(ni-7,:,:) + & + a(ni-6,ni-1) * var(ni-6,:,:) + & + a(ni-5,ni-1) * var(ni-5,:,:) + & + a(ni-4,ni-1) * var(ni-4,:,:) + & + a(ni-3,ni-1) * var(ni-3,:,:) + & + a(ni-2,ni-1) * var(ni-2,:,:) + & + a(ni-1,ni-1) * var(ni-1,:,:) + & + a(ni,ni-1) * var(ni,:,:) ) * idel + rhs(ni,:,:) = rhs(ni,:,:) + & + ( a(ni-3,ni) * var(ni-3,:,:) + & + a(ni-2,ni) * var(ni-2,:,:) + & + a(ni-1,ni) * var(ni-1,:,:) + & + a(ni,ni) * var(ni,:,:) ) * idel + + ir = ni - 7 + + end if + + do i = il, ir + + rhs(i,:,:) = rhs(i,:,:) + & + ( a(i-3,i) * var(i-3,:,:) + & + a(i-2,i) * var(i-2,:,:) + & + a(i-1,i) * var(i-1,:,:) + & + a(i,i) * var(i,:,:) + & + a(i+1,i) * var(i+1,:,:) + & + a(i+2,i) * var(i+2,:,:) + & + a(i+3,i) * var(i+3,:,:) ) * idel + + end do + + deallocate ( a, d, b, h, tmp ) + + allocate ( a(nj,nj), d(nj,nj), b(nj,nj), h(nj,nj), tmp(nj,nj) ) + a = zero; d = zero; b = zero; h = zero; + + call set_bmatrix ( b, bb(3:4), lsh(2), gsh(2), lbnd(2), delta(2), dfl ) + + call set_hmatrix ( h, bb(3:4) ) + + if ( any ( bb(3:4) == 0 ) ) then + + call set_dmatrix ( d, bb(3:4) ) + + a = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) + + else + + if ( lsh(2) /= gsh(2) ) then + call CCTK_WARN(0, 'Internal error. Bounding box inconsistent with grid' ) + end if + + center = gsh(2) / 2 + jr = center + gsize(2) + + call set_dmatrix_half ( d(1:jr,1:jr), 1 ) + + tmp(1:jr,1:jr) = -transpose ( & + matmul ( h(1:jr,1:jr), & + matmul ( & + transpose ( d(1:jr,1:jr) ), & + matmul ( b(1:jr,1:jr), & + d(1:jr,1:jr) ) ) ) ) + + a(1:jr,1:center) = tmp(1:jr,1:center) + + jl = center + 1 - gsize(1) + d = zero + call set_dmatrix_half ( d(jl:nj,jl:nj), 2 ) + tmp(jl:nj,jl:nj) = -transpose ( & + matmul ( h(jl:nj,jl:nj), & + matmul ( & + transpose ( d(jl:nj,jl:nj) ), & + matmul ( b(jl:nj,jl:nj), & + d(jl:nj,jl:nj) ) ) ) ) + + a(jl:nj,center+1:nj) = tmp(jl:nj,center+1:nj) + + end if + + idel = epsdisl / 64.0_wp + + if ( bb(3) == 0 ) then + + jl = 1 + gsize(2) + + else + + rhs(:,1,:) = rhs(:,1,:) + & + ( a(1,1) * var(:,1,:) + a(2,1) * var(:,2,:) + & + a(3,1) * var(:,3,:) + a(4,1) * var(:,4,:) ) * idel + rhs(:,2,:) = rhs(:,2,:) + & + ( a(1,2) * var(:,1,:) + a(2,2) * var(:,2,:) + & + a(3,2) * var(:,3,:) + a(4,2) * var(:,4,:) + & + a(5,2) * var(:,5,:) + a(6,2) * var(:,6,:) + & + a(7,2) * var(:,7,:) + a(8,2) * var(:,8,:) + & + a(9,2) * var(:,9,:) + a(10,2) * var(:,10,:) ) * idel + rhs(:,3,:) = rhs(:,3,:) + & + ( a(1,3) * var(:,1,:) + a(2,3) * var(:,2,:) + & + a(3,3) * var(:,3,:) + a(4,3) * var(:,4,:) + & + a(5,3) * var(:,5,:) + a(6,3) * var(:,6,:) + & + a(7,3) * var(:,7,:) + a(8,3) * var(:,8,:) + & + a(9,3) * var(:,9,:) + a(10,3) * var(:,10,:) ) * idel + rhs(:,4,:) = rhs(:,4,:) + & + ( a(1,4) * var(:,1,:) + a(2,4) * var(:,2,:) + & + a(3,4) * var(:,3,:) + a(4,4) * var(:,4,:) + & + a(5,4) * var(:,5,:) + a(6,4) * var(:,6,:) + & + a(7,4) * var(:,7,:) + a(8,4) * var(:,8,:) + & + a(9,4) * var(:,9,:) + a(10,4) * var(:,10,:) ) * idel + rhs(:,5,:) = rhs(:,5,:) + & + ( a(1,5) * var(:,1,:) + a(2,5) * var(:,2,:) + & + a(3,5) * var(:,3,:) + a(4,5) * var(:,4,:) + & + a(5,5) * var(:,5,:) + a(6,5) * var(:,6,:) + & + a(7,5) * var(:,7,:) + a(8,5) * var(:,8,:) + & + a(9,5) * var(:,9,:) + a(10,5) * var(:,10,:) ) * idel + rhs(:,6,:) = rhs(:,6,:) + & + ( a(1,6) * var(:,1,:) + a(2,6) * var(:,2,:) + & + a(3,6) * var(:,3,:) + a(4,6) * var(:,4,:) + & + a(5,6) * var(:,5,:) + a(6,6) * var(:,6,:) + & + a(7,6) * var(:,7,:) + a(8,6) * var(:,8,:) + & + a(9,6) * var(:,9,:) + a(10,6) * var(:,10,:) ) * idel + rhs(:,7,:) = rhs(:,7,:) + & + ( a(1,7) * var(:,1,:) + a(2,7) * var(:,2,:) + & + a(3,7) * var(:,3,:) + a(4,7) * var(:,4,:) + & + a(5,7) * var(:,5,:) + a(6,7) * var(:,6,:) + & + a(7,7) * var(:,7,:) + a(8,7) * var(:,8,:) + & + a(9,7) * var(:,9,:) + a(10,7) * var(:,10,:) ) * idel + + jl = 8 + + end if + + if ( bb(4) == 0 ) then + + jr = nj - gsize(2) + + else + + rhs(:,nj-6,:) = rhs(:,nj-6,:) + & + ( a(nj-9,nj-6) * var(:,nj-9,:) + & + a(nj-8,nj-6) * var(:,nj-8,:) + & + a(nj-7,nj-6) * var(:,nj-7,:) + & + a(nj-6,nj-6) * var(:,nj-6,:) + & + a(nj-5,nj-6) * var(:,nj-5,:) + & + a(nj-4,nj-6) * var(:,nj-4,:) + & + a(nj-3,nj-6) * var(:,nj-3,:) + & + a(nj-2,nj-6) * var(:,nj-2,:) + & + a(nj-1,nj-6) * var(:,nj-1,:) + & + a(nj,nj-6) * var(:,nj,:) ) * idel + rhs(:,nj-5,:) = rhs(:,nj-5,:) + & + ( a(nj-9,nj-5) * var(:,nj-9,:) + & + a(nj-8,nj-5) * var(:,nj-8,:) + & + a(nj-7,nj-5) * var(:,nj-7,:) + & + a(nj-6,nj-5) * var(:,nj-6,:) + & + a(nj-5,nj-5) * var(:,nj-5,:) + & + a(nj-4,nj-5) * var(:,nj-4,:) + & + a(nj-3,nj-5) * var(:,nj-3,:) + & + a(nj-2,nj-5) * var(:,nj-2,:) + & + a(nj-1,nj-5) * var(:,nj-1,:) + & + a(nj,nj-5) * var(:,nj,:) ) * idel + rhs(:,nj-4,:) = rhs(:,nj-4,:) + & + ( a(nj-9,nj-4) * var(:,nj-9,:) + & + a(nj-8,nj-4) * var(:,nj-8,:) + & + a(nj-7,nj-4) * var(:,nj-7,:) + & + a(nj-6,nj-4) * var(:,nj-6,:) + & + a(nj-5,nj-4) * var(:,nj-5,:) + & + a(nj-4,nj-4) * var(:,nj-4,:) + & + a(nj-3,nj-4) * var(:,nj-3,:) + & + a(nj-2,nj-4) * var(:,nj-2,:) + & + a(nj-1,nj-4) * var(:,nj-1,:) + & + a(nj,nj-4) * var(:,nj,:) ) * idel + rhs(:,nj-3,:) = rhs(:,nj-3,:) + & + ( a(nj-9,nj-3) * var(:,nj-9,:) + & + a(nj-8,nj-3) * var(:,nj-8,:) + & + a(nj-7,nj-3) * var(:,nj-7,:) + & + a(nj-6,nj-3) * var(:,nj-6,:) + & + a(nj-5,nj-3) * var(:,nj-5,:) + & + a(nj-4,nj-3) * var(:,nj-4,:) + & + a(nj-3,nj-3) * var(:,nj-3,:) + & + a(nj-2,nj-3) * var(:,nj-2,:) + & + a(nj-1,nj-3) * var(:,nj-1,:) + & + a(nj,nj-3) * var(:,nj,:) ) * idel + rhs(:,nj-2,:) = rhs(:,nj-2,:) + & + ( a(nj-9,nj-2) * var(:,nj-9,:) + & + a(nj-8,nj-2) * var(:,nj-8,:) + & + a(nj-7,nj-2) * var(:,nj-7,:) + & + a(nj-6,nj-2) * var(:,nj-6,:) + & + a(nj-5,nj-2) * var(:,nj-5,:) + & + a(nj-4,nj-2) * var(:,nj-4,:) + & + a(nj-3,nj-2) * var(:,nj-3,:) + & + a(nj-2,nj-2) * var(:,nj-2,:) + & + a(nj-1,nj-2) * var(:,nj-1,:) + & + a(nj,nj-2) * var(:,nj,:) ) * idel + rhs(:,nj-1,:) = rhs(:,nj-1,:) + & + ( a(nj-9,nj-1) * var(:,nj-9,:) + & + a(nj-8,nj-1) * var(:,nj-8,:) + & + a(nj-7,nj-1) * var(:,nj-7,:) + & + a(nj-6,nj-1) * var(:,nj-6,:) + & + a(nj-5,nj-1) * var(:,nj-5,:) + & + a(nj-4,nj-1) * var(:,nj-4,:) + & + a(nj-3,nj-1) * var(:,nj-3,:) + & + a(nj-2,nj-1) * var(:,nj-2,:) + & + a(nj-1,nj-1) * var(:,nj-1,:) + & + a(nj,nj-1) * var(:,nj,:) ) * idel + rhs(:,nj,:) = rhs(:,nj,:) + & + ( a(nj-3,nj) * var(:,nj-3,:) + & + a(nj-2,nj) * var(:,nj-2,:) + & + a(nj-1,nj) * var(:,nj-1,:) + & + a(nj,nj) * var(:,nj,:) ) * idel + + jr = nj - 7 + + end if + + do j = jl, jr + + rhs(:,j,:) = rhs(:,j,:) + & + ( a(j-3,j) * var(:,j-3,:) + & + a(j-2,j) * var(:,j-2,:) + & + a(j-1,j) * var(:,j-1,:) + & + a(j,j) * var(:,j,:) + & + a(j+1,j) * var(:,j+1,:) + & + a(j+2,j) * var(:,j+2,:) + & + a(j+3,j) * var(:,j+3,:) ) * idel + + end do + + deallocate ( a, d, b, h, tmp ) + + allocate ( a(nk,nk), d(nk,nk), b(nk,nk), h(nk,nk), tmp(nk,nk) ) + a = zero; d = zero; b = zero; h = zero; tmp = zero + + call set_bmatrix ( b, bb(5:6), lsh(3), gsh(3), lbnd(3), delta(3), dfl ) + + call set_hmatrix ( h, bb(5:6) ) + + if ( any ( bb(5:6) == 0 ) ) then + + call set_dmatrix ( d, bb(5:6) ) + + a = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) + + else + + if ( lsh(3) /= gsh(3) ) then + call CCTK_WARN(0, 'Internal error. Bounding box inconsistent with grid' ) + end if + + center = gsh(3) / 2 + kr = center + gsize(3) + + call set_dmatrix_half ( d(1:kr,1:kr), 1 ) + + tmp(1:kr,1:kr) = -transpose ( & + matmul ( h(1:kr,1:kr), & + matmul ( & + transpose ( d(1:kr,1:kr) ), & + matmul ( b(1:kr,1:kr), & + d(1:kr,1:kr) ) ) ) ) + + a(1:kr,1:center) = tmp(1:kr,1:center) + + kl = center + 1 - gsize(3) + d = zero + call set_dmatrix_half ( d(kl:nk,kl:nk), 2 ) + tmp(kl:nk,kl:nk) = -transpose ( & + matmul ( h(kl:nk,kl:nk), & + matmul ( & + transpose ( d(kl:nk,kl:nk) ), & + matmul ( b(kl:nk,kl:nk), & + d(kl:nk,kl:nk) ) ) ) ) + + a(kl:nk,center+1:nk) = tmp(kl:nk,center+1:nk) + + end if + + idel = epsdisl / 64.0_wp + + if ( bb(5) == 0 ) then + + kl = 1 + gsize(3) + + else + + rhs(:,:,1) = rhs(:,:,1) + & + ( a(1,1) * var(:,:,1) + a(2,1) * var(:,:,2) + & + a(3,1) * var(:,:,3) + a(4,1) * var(:,:,4) ) * idel + rhs(:,:,2) = rhs(:,:,2) + & + ( a(1,2) * var(:,:,1) + a(2,2) * var(:,:,2) + & + a(3,2) * var(:,:,3) + a(4,2) * var(:,:,4) + & + a(5,2) * var(:,:,5) + a(6,2) * var(:,:,6) + & + a(7,2) * var(:,:,7) + a(8,2) * var(:,:,8) + & + a(9,2) * var(:,:,9) + a(10,2) * var(:,:,10) ) * idel + rhs(:,:,3) = rhs(:,:,3) + & + ( a(1,3) * var(:,:,1) + a(2,3) * var(:,:,2) + & + a(3,3) * var(:,:,3) + a(4,3) * var(:,:,4) + & + a(5,3) * var(:,:,5) + a(6,3) * var(:,:,6) + & + a(7,3) * var(:,:,7) + a(8,3) * var(:,:,8) + & + a(9,3) * var(:,:,9) + a(10,3) * var(:,:,10) ) * idel + rhs(:,:,4) = rhs(:,:,4) + & + ( a(1,4) * var(:,:,1) + a(2,4) * var(:,:,2) + & + a(3,4) * var(:,:,3) + a(4,4) * var(:,:,4) + & + a(5,4) * var(:,:,5) + a(6,4) * var(:,:,6) + & + a(7,4) * var(:,:,7) + a(8,4) * var(:,:,8) + & + a(9,4) * var(:,:,9) + a(10,4) * var(:,:,10) ) * idel + rhs(:,:,5) = rhs(:,:,5) + & + ( a(1,5) * var(:,:,1) + a(2,5) * var(:,:,2) + & + a(3,5) * var(:,:,3) + a(4,5) * var(:,:,4) + & + a(5,5) * var(:,:,5) + a(6,5) * var(:,:,6) + & + a(7,5) * var(:,:,7) + a(8,5) * var(:,:,8) + & + a(9,5) * var(:,:,9) + a(10,5) * var(:,:,10) ) * idel + rhs(:,:,6) = rhs(:,:,6) + & + ( a(1,6) * var(:,:,1) + a(2,6) * var(:,:,2) + & + a(3,6) * var(:,:,3) + a(4,6) * var(:,:,4) + & + a(5,6) * var(:,:,5) + a(6,6) * var(:,:,6) + & + a(7,6) * var(:,:,7) + a(8,6) * var(:,:,8) + & + a(9,6) * var(:,:,9) + a(10,6) * var(:,:,10) ) * idel + rhs(:,:,7) = rhs(:,:,7) + & + ( a(1,7) * var(:,:,1) + a(2,7) * var(:,:,2) + & + a(3,7) * var(:,:,3) + a(4,7) * var(:,:,4) + & + a(5,7) * var(:,:,5) + a(6,7) * var(:,:,6) + & + a(7,7) * var(:,:,7) + a(8,7) * var(:,:,8) + & + a(9,7) * var(:,:,9) + a(10,7) * var(:,:,10) ) * idel + + kl = 8 + + end if + + if ( bb(6) == 0 ) then + + kr = nk - gsize(3) + + else + + rhs(:,:,nk-6) = rhs(:,:,nk-6) + & + ( a(nk-9,nk-6) * var(:,:,nk-9) + & + a(nk-8,nk-6) * var(:,:,nk-8) + & + a(nk-7,nk-6) * var(:,:,nk-7) + & + a(nk-6,nk-6) * var(:,:,nk-6) + & + a(nk-5,nk-6) * var(:,:,nk-5) + & + a(nk-4,nk-6) * var(:,:,nk-4) + & + a(nk-3,nk-6) * var(:,:,nk-3) + & + a(nk-2,nk-6) * var(:,:,nk-2) + & + a(nk-1,nk-6) * var(:,:,nk-1) + & + a(nk,nk-6) * var(:,:,nk) ) * idel + rhs(:,:,nk-5) = rhs(:,:,nk-5) + & + ( a(nk-9,nk-5) * var(:,:,nk-9) + & + a(nk-8,nk-5) * var(:,:,nk-8) + & + a(nk-7,nk-5) * var(:,:,nk-7) + & + a(nk-6,nk-5) * var(:,:,nk-6) + & + a(nk-5,nk-5) * var(:,:,nk-5) + & + a(nk-4,nk-5) * var(:,:,nk-4) + & + a(nk-3,nk-5) * var(:,:,nk-3) + & + a(nk-2,nk-5) * var(:,:,nk-2) + & + a(nk-1,nk-5) * var(:,:,nk-1) + & + a(nk,nk-5) * var(:,:,nk) ) * idel + rhs(:,:,nk-4) = rhs(:,:,nk-4) + & + ( a(nk-9,nk-4) * var(:,:,nk-9) + & + a(nk-8,nk-4) * var(:,:,nk-8) + & + a(nk-7,nk-4) * var(:,:,nk-7) + & + a(nk-6,nk-4) * var(:,:,nk-6) + & + a(nk-5,nk-4) * var(:,:,nk-5) + & + a(nk-4,nk-4) * var(:,:,nk-4) + & + a(nk-3,nk-4) * var(:,:,nk-3) + & + a(nk-2,nk-4) * var(:,:,nk-2) + & + a(nk-1,nk-4) * var(:,:,nk-1) + & + a(nk,nk-4) * var(:,:,nk) ) * idel + rhs(:,:,nk-3) = rhs(:,:,nk-3) + & + ( a(nk-9,nk-3) * var(:,:,nk-9) + & + a(nk-8,nk-3) * var(:,:,nk-8) + & + a(nk-7,nk-3) * var(:,:,nk-7) + & + a(nk-6,nk-3) * var(:,:,nk-6) + & + a(nk-5,nk-3) * var(:,:,nk-5) + & + a(nk-4,nk-3) * var(:,:,nk-4) + & + a(nk-3,nk-3) * var(:,:,nk-3) + & + a(nk-2,nk-3) * var(:,:,nk-2) + & + a(nk-1,nk-3) * var(:,:,nk-1) + & + a(nk,nk-3) * var(:,:,nk) ) * idel + rhs(:,:,nk-2) = rhs(:,:,nk-2) + & + ( a(nk-9,nk-2) * var(:,:,nk-9) + & + a(nk-8,nk-2) * var(:,:,nk-8) + & + a(nk-7,nk-2) * var(:,:,nk-7) + & + a(nk-6,nk-2) * var(:,:,nk-6) + & + a(nk-5,nk-2) * var(:,:,nk-5) + & + a(nk-4,nk-2) * var(:,:,nk-4) + & + a(nk-3,nk-2) * var(:,:,nk-3) + & + a(nk-2,nk-2) * var(:,:,nk-2) + & + a(nk-1,nk-2) * var(:,:,nk-1) + & + a(nk,nk-2) * var(:,:,nk) ) * idel + rhs(:,:,nk-1) = rhs(:,:,nk-1) + & + ( a(nk-9,nk-1) * var(:,:,nk-9) + & + a(nk-8,nk-1) * var(:,:,nk-8) + & + a(nk-7,nk-1) * var(:,:,nk-7) + & + a(nk-6,nk-1) * var(:,:,nk-6) + & + a(nk-5,nk-1) * var(:,:,nk-5) + & + a(nk-4,nk-1) * var(:,:,nk-4) + & + a(nk-3,nk-1) * var(:,:,nk-3) + & + a(nk-2,nk-1) * var(:,:,nk-2) + & + a(nk-1,nk-1) * var(:,:,nk-1) + & + a(nk,nk-1) * var(:,:,nk) ) * idel + rhs(:,:,nk) = rhs(:,:,nk) + & + ( a(nk-3,nk) * var(:,:,nk-3) + & + a(nk-2,nk) * var(:,:,nk-2) + & + a(nk-1,nk) * var(:,:,nk-1) + & + a(nk,nk) * var(:,:,nk) ) * idel + + kr = nk - 7 + + end if + + do k = kl, kr + + rhs(:,:,k) = rhs(:,:,k) + & + ( a(k-3,k) * var(:,:,k-3) + & + a(k-2,k) * var(:,:,k-2) + & + a(k-1,k) * var(:,:,k-1) + & + a(k,k) * var(:,:,k) + & + a(k+1,k) * var(:,:,k+1) + & + a(k+2,k) * var(:,:,k+2) + & + a(k+3,k) * var(:,:,k+3) ) * idel + + end do + + deallocate ( a, d, b, h, tmp ) + +contains + + subroutine set_dmatrix ( d, bb ) + + implicit none + + CCTK_REAL, dimension(:,:), intent(out) :: d + CCTK_INT, dimension(2), intent(in) :: bb + CCTK_INT :: n + CCTK_REAL, dimension(4), save :: ac = (/ -1.0_wp, 3.0_wp, -3.0_wp, 1.0_wp /) + CCTK_INT :: i + + n = size(d,1) + + if ( any(bb /= 0) ) then + + if ( bb(1) /= 0 ) then + + d(1,1:4) = ac + d(2,1:4) = ac + + do i = 3, n-1 + d(i,i-2:i+1) = ac + end do + + d(n,n-2:n) = ac(1:3) + + else + + d(1,1:3) = ac(2:4) + + do i = 2, n-2 + d(i,i-1:i+2) = ac + end do + + d(n-1,n-3:n) = ac + d(n,n-3:n) = ac + + end if + + else + + d(1,1:2) = ac(3:4) + d(2,1:3) = ac(2:4) + + do i = 3, n-1 + d(i,i-2:i+1) = ac + end do + + d(n,n-2:n) = ac(1:3) + + end if + + end subroutine set_dmatrix + + + subroutine set_dmatrix_half ( d, part ) + + implicit none + + CCTK_REAL, dimension(:,:), intent(out) :: d + CCTK_INT, intent(in) :: part + CCTK_INT :: n + CCTK_REAL, dimension(4), save :: ac = (/ -1.0_wp, 3.0_wp, -3.0_wp, 1.0_wp /) + CCTK_INT :: i + + n = size(d,1) + + if ( part == 1 ) then + + d(1,1:4) = ac + d(2,1:4) = ac + + do i = 3, n-1 + d(i,i-2:i+1) = ac + end do + + d(n,n-2:n) = ac(1:3) + + else if ( part == 2 ) then + + d(1,1:3) = ac(2:4) + + do i = 2, n-2 + d(i,i-1:i+2) = ac + end do + + d(n-1,n-3:n) = ac + d(n,n-3:n) = ac + + end if + + end subroutine set_dmatrix_half + + + subroutine set_bmatrix ( b, bb, lsh, gsh, lbnd, h, dfl ) + + implicit none + + CCTK_REAL, dimension(:,:), intent(out) :: b + CCTK_INT, dimension(2), intent(in) :: bb + CCTK_INT, intent(in) :: lsh, gsh, lbnd + CCTK_REAL, intent(in) :: h, dfl + CCTK_INT :: n, nwt + CCTK_INT :: i + CCTK_REAL :: xp, xmax + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + n = size(b,1) + xp = ( gsh - 1 ) * dfl - 1 + + do i = 1, n + b(i,i) = bf ( real(i+lbnd-1,wp), xp, real(gsh-1,wp), h ) + end do + + end subroutine set_bmatrix + + + subroutine set_hmatrix ( h, bb ) + + implicit none + + CCTK_REAL, dimension(:,:), intent(out) :: h + CCTK_INT, dimension(2), intent(in) :: bb + CCTK_INT :: n, i + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + n = size(h,1) + + do i = 1, n + h(i,i) = 1.0_wp + end do + + if ( bb(1) /= 0 ) then + h(1:7,1:7) = sigma ( ) + end if + + if ( bb(2) /= 0 ) then + h(n:n-6:-1,n:n-6:-1) = sigma ( ) + end if + + end subroutine set_hmatrix + + + function bf ( x, xp, xmax, h ) + + implicit none + + CCTK_REAL :: bf + CCTK_REAL, intent(in) :: x, xp, xmax, h + CCTK_REAL :: h2, h21, xp2, xp3 + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + if ( x < xp ) then + + h2 = h**2 + h21 = 1.0_wp - h2 + xp2 = 1.0_wp / xp**2 + xp3 = xp2 / xp + bf = ( -2.0_wp * h21 * xp3 * x + 3.0_wp * h21 * xp2 ) * x**2 + h2 + + else if ( x > xmax - xp ) then + + h2 = h**2 + h21 = 1.0_wp - h2 + xp2 = 1.0_wp / xp**2 + xp3 = xp2 / xp + bf = ( 2.0_wp * h21 * xp3 * ( x - xmax ) + & + 3.0_wp * h21 * xp2 ) * ( x - xmax )**2 + h2 + + else + + bf = 1.0_wp + + end if + + end function bf + + + function sigma ( ) + + implicit none + + CCTK_REAL, dimension(7,7) :: sigma + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + sigma(1,1) = 4.930709842221048047321555312222552006914_wp + sigma(2,1) = 0_wp + sigma(3,1) = 0_wp + sigma(4,1) = 0_wp + sigma(5,1) = 0_wp + sigma(6,1) = 0_wp + sigma(7,1) = 0_wp + sigma(1,2) = 0_wp + sigma(2,2) = 1.499128158967098993672808967791430710830_wp + sigma(3,2) = 0.8363320637711490458766332755492551638364_wp + sigma(4,2) = -0.2634194189384563273844251938374368097488_wp + sigma(5,2) = -0.1281142262970477369698649073120970227895_wp + sigma(6,2) = -0.1830121923043950334684235233648902467425_wp + sigma(7,2) = 0.01275879120879511857581469566228951245742_wp + sigma(1,3) = 0_wp + sigma(2,3) = 0.8363320637711490458766332755492551638364_wp + sigma(3,3) = 0.8785380112113613457715293301885586949116_wp + sigma(4,3) = 0.1454039706555408804118493444561963110723_wp + sigma(5,3) = -0.03346859376328961536059718572568704595897_wp + sigma(6,3) = -0.1759548942026651677201193466207997400242_wp + sigma(7,3) = 0.02020404709761823823384626194775259449799_wp + sigma(1,4) = 0_wp + sigma(2,4) = -0.2634194189384563273844251938374368097488_wp + sigma(3,4) = 0.1454039706555408804118493444561963110723_wp + sigma(4,4) = 0.5690707861055695191099807759535994201984_wp + sigma(5,4) = 0.3374169937858825026860519890688362591356_wp + sigma(6,4) = -0.03076801586199402080957339555604653302670_wp + sigma(7,4) = 0.004077760901760668929595481886645691229036_wp + sigma(1,5) = 0_wp + sigma(2,5) = -0.1281142262970477369698649073120970227895_wp + sigma(3,5) = -0.03346859376328961536059718572568704595897_wp + sigma(4,5) = 0.3374169937858825026860519890688362591356_wp + sigma(5,5) = 0.5001980842834785426891480539161939116986_wp + sigma(6,5) = 0.2904403943485980286998891564955850119600_wp + sigma(7,5) = -0.05655813410599246578174268632557910837116_wp + sigma(1,6) = 0_wp + sigma(2,6) = -0.1830121923043950334684235233648902467425_wp + sigma(3,6) = -0.1759548942026651677201193466207997400242_wp + sigma(4,6) = -0.03076801586199402080957339555604653302670_wp + sigma(5,6) = 0.2904403943485980286998891564955850119600_wp + sigma(6,6) = 0.8640467497799918355533313312208192675134_wp + sigma(7,6) = 0.03704582714017916999084939064121697271831_wp + sigma(1,7) = 0_wp + sigma(2,7) = 0.01275879120879511857581469566228951245742_wp + sigma(3,7) = 0.02020404709761823823384626194775259449799_wp + sigma(4,7) = 0.004077760901760668929595481886645691229036_wp + sigma(5,7) = -0.05655813410599246578174268632557910837116_wp + sigma(6,7) = 0.03704582714017916999084939064121697271831_wp + sigma(7,7) = 0.9909401061257062767072573096147576992963_wp + + end function sigma + +end subroutine dissipation_6_5_opt diff --git a/src/dissipation.c b/src/dissipation.c index 59daa23..899bcde 100644 --- a/src/dissipation.c +++ b/src/dissipation.c @@ -17,7 +17,7 @@ void CCTK_FCALL CCTK_FNAME(dissipation_4_3_opt) (const CCTK_REAL *var, const CCTK_REAL *epsdis, const CCTK_REAL *diss_fraction, CCTK_REAL *rhs); -void CCTK_FCALL CCTK_FNAME(dissipation_6_5) (const CCTK_REAL *var, +void CCTK_FCALL CCTK_FNAME(dissipation_6_5_opt) (const CCTK_REAL *var, const CCTK_INT *ni, const CCTK_INT *nj, const CCTK_INT *nk, @@ -141,9 +141,9 @@ apply (int const varindex, char const * const optstring, void * const arg) break; } case 6: { - CCTK_FNAME(dissipation_6_5) - (varptr, &cctk_lsh[0], &cctk_lsh[1], &cctk_lsh[2], - bbox, gsize, dx, &epsdis, rhsptr); + CCTK_FNAME(dissipation_6_5_opt) + (varptr, cctk_lsh, cctk_gsh, cctk_lbnd, + bbox, gsize, dx, &epsdis, &diss_fraction, rhsptr); break; } default: diff --git a/src/make.code.defn b/src/make.code.defn index 171c48f..c113cfe 100644 --- a/src/make.code.defn +++ b/src/make.code.defn @@ -20,6 +20,7 @@ SRCS = call_derivs.c \ Dissipation_4_3_min_err_coeff.F90 \ Dissipation_8_4.F90 \ Dissipation_6_5.F90 \ + Dissipation_6_5_min_err_coeff.F90 \ get_coeffs.c \ Coefficients_2_1.F90 \ Coefficients_4_2.F90 \ -- cgit v1.2.3