diff options
author | diener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a> | 2004-10-07 23:32:09 +0000 |
---|---|---|
committer | diener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a> | 2004-10-07 23:32:09 +0000 |
commit | 820cf085c7fba79fa34163a535e5d7a37f1a2697 (patch) | |
tree | 88d9a097ce6a59adcf1a6e3a2e1a7d9882f02e2f | |
parent | 52862a4d965f4a147fc2da64251a7981e9c08241 (diff) |
Added 8-4 derivatives. Note only Diff_gv has been updated to use it. Diff_gf
has to be modified in other ways as well.
git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@9 f69c4107-0314-4c4f-9ad4-17e986b73f4a
-rw-r--r-- | param.ccl | 2 | ||||
-rw-r--r-- | src/Derivatives_8_4.F90 | 457 | ||||
-rw-r--r-- | src/call_derivs.c | 13 | ||||
-rw-r--r-- | src/make.code.defn | 2 |
4 files changed, 472 insertions, 2 deletions
@@ -3,5 +3,5 @@ INT order "Order of accuracy" STEERABLE=always { - 2:6:2 :: "" + 2:8:2 :: "" } 2 diff --git a/src/Derivatives_8_4.F90 b/src/Derivatives_8_4.F90 new file mode 100644 index 0000000..fb312f1 --- /dev/null +++ b/src/Derivatives_8_4.F90 @@ -0,0 +1,457 @@ +#include "cctk.h" +#include "cctk_Functions.h" +#include "cctk_Parameters.h" + +subroutine deriv_gf_8_4 ( var, ni, nj, nk, dir, bb, gsize, delta, dvar ) + + implicit none + + DECLARE_CCTK_FUNCTIONS + DECLARE_CCTK_PARAMETERS + + CCTK_REAL, parameter :: zero = 0.0 + integer, parameter :: wp = kind(zero) + CCTK_INT, intent(IN) :: ni, nj, nk + CCTK_REAL, dimension(ni,nj,nk), intent(IN) :: var + CCTK_INT, intent(IN) :: dir + CCTK_INT, intent(IN) :: bb(2) + CCTK_INT, intent(IN) :: gsize + CCTK_REAL, intent(IN) :: delta + CCTK_REAL, dimension(ni,nj,nk), intent(OUT) :: dvar + + CCTK_REAL, dimension(4), save :: a + CCTK_REAL, dimension(12,8), save :: q + CCTK_REAL, parameter :: x1 = 0.541_wp, x2 = -0.0675_wp, x3 = 0.748_wp + CCTK_REAL :: idel + + CCTK_INT :: il, ir, jl, jr, kl, kr, j, k + + logical, save :: first = .true. + + if ( first ) then + a(1) = 4.0_wp/5.0_wp; a(2) = -1.0_wp/5.0_wp + a(3) = 4.0_wp/105.0_wp; a(4) = -1.0_wp/280.0_wp + + q(1,1) = -2540160.0_wp/1498139.0_wp + q(2,1) = 9.0_wp * ( 2257920.0_wp*x1 + 11289600.0_wp*x2 + & + 22579200.0_wp*x3 - 15849163.0_wp ) / 5992556.0_wp + q(3,1) = 3.0_wp * ( -33868800.0_wp*x1 - 162570240.0_wp*x2 - & + 304819200.0_wp*x3 + 235236677.0_wp ) / 5992556.0_wp + q(4,1) = ( 609638400.0_wp*x1 + 2743372800.0_wp*x2 + & + 4572288000.0_wp * x3 - 3577778591.0_wp ) / 17977668.0_wp + q(5,1) = 3.0_wp * ( -16934400*x1 - 67737600.0_wp*x2 - & + 84672000.0_wp*x3 + 67906303.0_wp ) / 1498139.0_wp + q(6,1) = 105.0_wp * ( 967680.0_wp*x1 + 2903040.0_wp*x2 - & + 305821.0_wp ) / 5992556.0_wp + q(7,1) = 49.0_wp * ( -1244160.0_wp*x1 + 18662400.0_wp * x3 - & + 13322233.0_wp ) / 17977668.0_wp + q(8,1) = 3.0_wp * ( -6773760.0_wp*x2 - 33868800.0_wp*x3 + & + 24839327.0_wp ) / 5992556.0_wp + q(9,1) = zero; q(10,1) = zero; q(11,1) = zero; q(12,1) = zero + + q(1,2) = 9.0_wp * ( -2257920.0_wp*x1 - 11289600.0_wp*x2 - & + 22579200.0_wp*x3 + 15849163.0_wp ) / 31004596.0_wp + q(2,2) = zero + q(3,2) = 3.0_wp * ( 7257600.0_wp*x1 + 33868800.0_wp*x2 + & + 60963840.0_wp*x3 - 47167457.0_wp ) / 2214614.0_wp + q(4,2) = 3.0_wp * ( -9676800.0_wp*x1 - 42336000.0_wp*x2 - & + 67737600.0_wp*x3 + 53224573.0_wp ) / 1107307.0_wp + q(5,2) = 7.0_wp * ( 55987200.0_wp*x1 + 217728000.0_wp*x2 + & + 261273600.0_wp*x3 - 211102099.0_wp ) / 13287684.0_wp + q(6,2) = 3.0_wp * ( -11612160.0_wp*x1 - 33868800.0_wp*x2 + & + 3884117.0_wp ) / 2214614.0_wp + q(7,2) = 150.0_wp * ( 24192.0_wp*x1 - 338688.0_wp*x3 + & + 240463.0_wp ) / 1107307.0_wp + q(8,2) = ( 152409600.0_wp*x2 + 731566080.0_wp*x3 - & + 536324953.0_wp ) / 46506894.0_wp + q(9,2) = zero; q(10,2) = zero; q(11,2) = zero; q(12,2) = zero + + q(1,3) = ( 33868800.0_wp*x1 + 162570240.0_wp*x2 + & + 304819200.0_wp*x3 - 235236677.0_wp ) / 1743924.0_wp + q(2,3) = ( -7257600.0_wp*x1 - 33868800.0_wp*x2 - & + 60963840.0_wp*x3 + 47167457.0_wp ) / 124566.0_wp + q(3,3) = zero + q(4,3) = ( 24192000.0_wp*x1 + 101606400.0_wp*x2 + & + 152409600.0_wp*x3 - 120219461.0_wp ) / 124566.0_wp + q(5,3) = ( -72576000.0_wp*x1 - 270950400.0_wp*x2 - & + 304819200.0_wp*x3 + 249289259.0_wp ) / 249132.0_wp + q(6,3) = 9.0_wp * ( 806400.0_wp*x1 + 2257920.0_wp*x2 - & + 290167.0_wp ) / 41522.0_wp + q(7,3) = 6.0_wp * ( -134400.0_wp*x1 + 1693440.0_wp*x3 - & + 1191611.0_wp ) / 20761.0_wp + q(8,3) = 5.0_wp * ( -2257920.0_wp*x2 - 10160640.0_wp*x3 + & + 7439833.0_wp ) / 290654.0_wp + q(9,3) = zero; q(10,3) = zero; q(11,3) = zero; q(12,3) = zero + + q(1,4) = ( -609638400.0_wp*x1 - 2743372800.0_wp*x2 - & + 4572288000.0_wp*x3 + 3577778591.0_wp ) / 109619916.0_wp + q(2,4) = 3.0_wp * ( 9676800.0_wp*x1 + 42336000.0_wp*x2 + & + 67737600.0_wp*x3 - 53224573.0_wp ) / 1304999.0_wp + q(3,4) = 3.0_wp * ( -24192000.0_wp*x1 - 101606400.0_wp*x2 - & + 152409600.0_wp*x3 + 120219461.0_wp ) / 2609998.0_wp + q(4,4) = zero + q(5,4) = 9.0_wp * ( 16128000.0_wp*x1 + 56448000.0_wp*x2 + & + 56448000.0_wp*x3 - 47206049.0_wp ) / 5219996.0_wp + q(6,4) = 3.0_wp * ( -19353600.0_wp*x1 - 50803200.0_wp*x2 + & + 7628371.0_wp ) / 2609998.0_wp + q(7,4) = 2.0_wp * ( 10886400.0_wp*x1 - 114307200.0_wp*x3 + & + 79048289.0_wp ) / 3914997.0_wp + q(8,4) = 75.0_wp * ( 1354752.0_wp*x2 + 5419008.0_wp*x3 - & + 3952831.0_wp ) / 18269986.0_wp + q(9,4) = zero; q(10,4) = zero; q(11,4) = zero; q(12,4) = zero + + q(1,5) = 3.0_wp * ( 16934400.0_wp*x1 + 67737600.0_wp*x2 + & + 84672000.0_wp*x3 - 67906303.0_wp ) / 2096689.0_wp + q(2,5) = 7.0_wp * ( -55987200.0_wp*x1 - 217728000.0_wp*x2 - & + 261273600.0_wp*x3 + 211102099.0_wp ) / 3594324.0_wp + q(3,5) = 3.0_wp * ( 72576000.0_wp*x1 + 270950400.0_wp*x2 + & + 304819200.0_wp*x3 - 249289259.0_wp ) / 1198108.0_wp + q(4,5) = 9.0_wp * ( -16128000.0_wp*x1 - 56448000.0_wp*x2 - & + 56448000.0_wp*x3 + 47206049.0_wp ) / 1198108.0_wp + q(5,5) = zero + q(6,5) = 105.0_wp * ( 414720.0_wp*x1 + 967680.0_wp*x2 - & + 165527.0_wp ) / 1198108.0_wp + q(7,5) = 15.0_wp * ( -967680.0_wp*x1 + 6773760.0_wp*x3 - & + 4472029.0_wp ) / 1198108.0_wp + q(8,5) = ( -304819200.0_wp*x2 - 914457600.0_wp*x3 + & + 657798011.0_wp ) / 25160268.0_wp + q(9,5) = -2592.0_wp/299527.0_wp; + q(10,5) = zero; q(11,5) = zero; q(12,5) = zero + + q(1,6) = 5.0_wp * ( -967680.0_wp*x1 - 2903040.0_wp*x2 + & + 305821.0_wp ) / 1237164.0_wp + q(2,6) = ( 11612160.0_wp*x1 + 33868800.0_wp*x2 - & + 3884117.0_wp ) / 618582.0_wp + q(3,6) = 9.0_wp * ( -806400.0_wp*x1 - 2257920.0_wp*x2 + & + 290167.0_wp ) / 206194.0_wp + q(4,6) = ( 19353600.0_wp*x1 + 50803200.0_wp*x2 - & + 7628371.0_wp ) / 618582.0_wp + q(5,6) = 35.0_wp * ( -414720.0_wp*x1 - 967680.0_wp*x2 + & + 165527.0_wp ) / 1237164.0_wp + q(6,6) = zero; q(7,6) = 80640.0_wp*x1 / 103097.0_wp + q(8,6) = 80640.0_wp*x2 / 103097.0_wp; q(9,6) = 3072.0_wp / 103097.0_wp + q(10,6) = -288.0_wp/103097.0_wp; q(11,6) = zero; q(12,6) = zero + + q(1,7) = 7.0_wp * ( 1244160.0_wp*x1 - 18662400.0_wp*x3 + & + 13322233.0_wp ) / 8041092.0_wp + q(2,7) = 150.0_wp * ( -24192.0_wp*x1 + 338688.0_wp*x3 - & + 240463.0_wp ) / 670091.0_wp + q(3,7) = 54.0_wp * ( 134400.0_wp*x1 - 1693440.0_wp*x3 + & + 1191611.0_wp ) / 670091.0_wp + q(4,7) = 2.0_wp * ( -10886400.0_wp*x1 + 114307200.0_wp*x3 - & + 79048289.0_wp ) / 2010273.0_wp + q(5,7) = 15.0_wp * ( 967680.0_wp*x1 - 6773760.0_wp*x3 + & + 4472029.0_wp ) / 2680364.0_wp + q(6,7) = -725760.0_wp*x1 / 670091.0_wp; q(7,7) = zero + q(8,7) = 725760.0_wp*x3 / 670091.0_wp; q(9,7) = -145152.0_wp/670091.0_wp + q(10,7) = 27648.0_wp/670091.0_wp; q(11,7) = -2592.0_wp/670091.0_wp + q(12,7) = zero + + q(1,8) = 3.0_wp * ( 6773760.0_wp*x2 + 33868800.0_wp*x3 - & + 24839327.0_wp ) / 20510956.0_wp + q(2,8) = ( -152409600.0_wp*x2 - 731566080.0_wp*x3 + & + 536324953.0_wp ) / 30766434.0_wp + q(3,8) = 45.0_wp * ( 2257920.0_wp*x2 + 10160640.0_wp*x3 - & + 7439833.0_wp ) / 10255478.0_wp + q(4,8) = 75.0_wp * ( -1354752.0_wp*x2 - 5419008.0_wp*x3 + & + 3952831.0_wp ) / 10255478.0_wp + q(5,8) = ( 304819200.0_wp*x2 + 914457600.0_wp*x3 - & + 657798011.0_wp ) / 61532868.0_wp + q(6,8) = -5080320.0_wp*x2 / 5127739.0_wp + q(7,8) = -5080320.0_wp*x3 / 5127739.0_wp; q(8,8) = zero + q(9,8) = 4064256.0_wp/5127739.0_wp; q(10,8) = -1016064.0_wp/5127739.0_wp + q(11,8) = 193536.0_wp/5127739.0_wp; q(12,8) = -18144.0_wp/5127739.0_wp + + first = .false. + end if + + idel = 1.0_wp / delta + + direction: select case (dir) + case (0) direction + if ( bb(1) == 0 ) then + il = 1 + gsize + else + dvar(1,:,:) = ( q(1,1) * var(1,:,:) + q(2,1) * var(2,:,:) + & + q(3,1) * var(3,:,:) + q(4,1) * var(4,:,:) + & + q(5,1) * var(5,:,:) + q(6,1) * var(6,:,:) + & + q(7,1) * var(7,:,:) + q(8,1) * var(8,:,:) ) * idel + dvar(2,:,:) = ( q(1,2) * var(1,:,:) + q(3,2) * var(3,:,:) + & + q(4,2) * var(4,:,:) + q(5,2) * var(5,:,:) + & + q(6,2) * var(6,:,:) + q(7,2) * var(7,:,:) + & + q(8,2) * var(8,:,:) ) * idel + dvar(3,:,:) = ( q(1,3) * var(1,:,:) + q(2,3) * var(2,:,:) + & + q(4,3) * var(4,:,:) + q(5,3) * var(5,:,:) + & + q(6,3) * var(6,:,:) + q(7,3) * var(7,:,:) + & + q(8,3) * var(8,:,:) ) * idel + dvar(4,:,:) = ( q(1,4) * var(1,:,:) + q(2,4) * var(2,:,:) + & + q(3,4) * var(3,:,:) + q(5,4) * var(5,:,:) + & + q(6,4) * var(6,:,:) + q(7,4) * var(7,:,:) + & + q(8,4) * var(8,:,:) ) * idel + dvar(5,:,:) = ( q(1,5) * var(1,:,:) + q(2,5) * var(2,:,:) + & + q(3,5) * var(3,:,:) + q(4,5) * var(4,:,:) + & + q(6,5) * var(6,:,:) + q(7,5) * var(7,:,:) + & + q(8,5) * var(8,:,:) + q(9,5) * var(9,:,:) ) * idel + dvar(6,:,:) = ( q(1,6) * var(1,:,:) + q(2,6) * var(2,:,:) + & + q(3,6) * var(3,:,:) + q(4,6) * var(4,:,:) + & + q(5,6) * var(5,:,:) + q(7,6) * var(7,:,:) + & + q(8,6) * var(8,:,:) + q(9,6) * var(9,:,:) + & + q(10,6) * var(10,:,:) ) * idel + dvar(7,:,:) = ( q(1,7) * var(1,:,:) + q(2,7) * var(2,:,:) + & + q(3,7) * var(3,:,:) + q(4,7) * var(4,:,:) + & + q(5,7) * var(5,:,:) + q(6,7) * var(6,:,:) + & + q(8,7) * var(8,:,:) + q(9,7) * var(9,:,:) + & + q(10,7) * var(10,:,:) + q(11,7) * var(11,:,:) ) * idel + dvar(8,:,:) = ( q(1,8) * var(1,:,:) + q(2,8) * var(2,:,:) + & + q(3,8) * var(3,:,:) + q(4,8) * var(4,:,:) + & + q(5,8) * var(5,:,:) + q(6,8) * var(6,:,:) + & + q(7,8) * var(7,:,:) + q(9,8) * var(9,:,:) + & + q(10,8) * var(10,:,:) + q(11,8) * var(11,:,:) + & + q(12,8) * var(12,:,:) ) * idel + il = 9 + end if + if ( bb(2) == 0 ) then + ir = ni - gsize + else + dvar(ni,:,:) = - ( q(1,1) * var(ni,:,:) + q(2,1) * var(ni-1,:,:) + & + q(3,1) * var(ni-2,:,:) + q(4,1) * var(ni-3,:,:) + & + q(5,1) * var(ni-4,:,:) + q(6,1) * var(ni-5,:,:) + & + q(7,1) * var(ni-6,:,:) + & + q(8,1) * var(ni-7,:,:) ) * idel + dvar(ni-1,:,:) = - ( q(1,2) * var(ni,:,:) + q(3,2) * var(ni-2,:,:) + & + q(4,2) * var(ni-3,:,:) + q(5,2) * var(ni-4,:,:) + & + q(6,2) * var(ni-5,:,:) + q(7,2) * var(ni-6,:,:) + & + q(8,2) * var(ni-7,:,:) ) * idel + dvar(ni-2,:,:) = - ( q(1,3) * var(ni,:,:) + q(2,3) * var(ni-1,:,:) + & + q(4,3) * var(ni-3,:,:) + q(5,3) * var(ni-4,:,:) + & + q(6,3) * var(ni-5,:,:) + q(7,3) * var(ni-6,:,:) + & + q(8,3) * var(ni-7,:,:) ) * idel + dvar(ni-3,:,:) = - ( q(1,4) * var(ni,:,:) + q(2,4) * var(ni-1,:,:) + & + q(3,4) * var(ni-2,:,:) + q(5,4) * var(ni-4,:,:) + & + q(6,4) * var(ni-5,:,:) + q(7,4) * var(ni-6,:,:) + & + q(8,4) * var(ni-7,:,:) ) * idel + dvar(ni-4,:,:) = - ( q(1,5) * var(ni,:,:) + q(2,5) * var(ni-1,:,:) + & + q(3,5) * var(ni-2,:,:) + q(4,5) * var(ni-3,:,:) + & + q(6,5) * var(ni-5,:,:) + q(7,5) * var(ni-6,:,:) + & + q(8,5) * var(ni-7,:,:) + & + q(9,5) * var(ni-8,:,:) ) * idel + dvar(ni-5,:,:) = - ( q(1,6) * var(ni,:,:) + q(2,6) * var(ni-1,:,:) + & + q(3,6) * var(ni-2,:,:) + q(4,6) * var(ni-3,:,:) + & + q(5,6) * var(ni-4,:,:) + q(7,6) * var(ni-6,:,:) + & + q(8,6) * var(ni-7,:,:) + q(9,6) * var(ni-8,:,:) + & + q(10,6) * var(ni-9,:,:) ) * idel + dvar(ni-6,:,:) = - ( q(1,7) * var(ni,:,:) + q(2,7) * var(ni-1,:,:) + & + q(3,7) * var(ni-2,:,:) + q(4,7) * var(ni-3,:,:) + & + q(5,7) * var(ni-4,:,:) + q(6,7) * var(ni-5,:,:) + & + q(8,7) * var(ni-7,:,:) + q(9,7) * var(ni-8,:,:) + & + q(10,7) * var(ni-9,:,:) + & + q(11,7) * var(ni-10,:,:) ) * idel + dvar(ni-7,:,:) = - ( q(1,8) * var(ni,:,:) + q(2,8) * var(ni-1,:,:) + & + q(3,8) * var(ni-2,:,:) + q(4,8) * var(ni-3,:,:) + & + q(5,8) * var(ni-4,:,:) + q(6,8) * var(ni-5,:,:) + & + q(7,8) * var(ni-6,:,:) + q(9,8) * var(ni-8,:,:) + & + q(10,8) * var(ni-9,:,:) + & + q(11,8) * var(ni-10,:,:) + & + q(12,8) * var(ni-11,:,:) ) * idel + ir = ni - 8 + end if + dvar(il:ir,:,:) = ( a(1) * ( var(il+1:ir+1,:,:) - & + var(il-1:ir-1,:,:) ) + & + a(2) * ( var(il+2:ir+2,:,:) - & + var(il-2:ir-2,:,:) ) + & + a(3) * ( var(il+3:ir+3,:,:) - & + var(il-3:ir-3,:,:) ) + & + a(4) * ( var(il+4:ir+4,:,:) - & + var(il-4:ir-4,:,:) ) ) * idel + case (1) direction + if ( bb(1) == 0 ) then + jl = 1 + gsize + else + dvar(:,1,:) = ( q(1,1) * var(:,1,:) + q(2,1) * var(:,2,:) + & + q(3,1) * var(:,3,:) + q(4,1) * var(:,4,:) + & + q(5,1) * var(:,5,:) + q(6,1) * var(:,6,:) + & + q(7,1) * var(:,7,:) + q(8,1) * var(:,8,:) ) * idel + dvar(:,2,:) = ( q(1,2) * var(:,1,:) + q(3,2) * var(:,3,:) + & + q(4,2) * var(:,4,:) + q(5,2) * var(:,5,:) + & + q(6,2) * var(:,6,:) + q(7,2) * var(:,7,:) + & + q(8,2) * var(:,8,:) ) * idel + dvar(:,3,:) = ( q(1,3) * var(:,1,:) + q(2,3) * var(:,2,:) + & + q(4,3) * var(:,4,:) + q(5,3) * var(:,5,:) + & + q(6,3) * var(:,6,:) + q(7,3) * var(:,7,:) + & + q(8,3) * var(:,8,:) ) * idel + dvar(:,4,:) = ( q(1,4) * var(:,1,:) + q(2,4) * var(:,2,:) + & + q(3,4) * var(:,3,:) + q(5,4) * var(:,5,:) + & + q(6,4) * var(:,6,:) + q(7,4) * var(:,7,:) + & + q(8,4) * var(:,8,:) ) * idel + dvar(:,5,:) = ( q(1,5) * var(:,1,:) + q(2,5) * var(:,2,:) + & + q(3,5) * var(:,3,:) + q(4,5) * var(:,4,:) + & + q(6,5) * var(:,6,:) + q(7,5) * var(:,7,:) + & + q(8,5) * var(:,8,:) + q(9,5) * var(:,9,:) ) * idel + dvar(:,6,:) = ( q(1,6) * var(:,1,:) + q(2,6) * var(:,2,:) + & + q(3,6) * var(:,3,:) + q(4,6) * var(:,4,:) + & + q(5,6) * var(:,5,:) + q(7,6) * var(:,7,:) + & + q(8,6) * var(:,8,:) + q(9,6) * var(:,9,:) + & + q(10,6) * var(:,10,:) ) * idel + dvar(:,7,:) = ( q(1,7) * var(:,1,:) + q(2,7) * var(:,2,:) + & + q(3,7) * var(:,3,:) + q(4,7) * var(:,4,:) + & + q(5,7) * var(:,5,:) + q(6,7) * var(:,6,:) + & + q(8,7) * var(:,8,:) + q(9,7) * var(:,9,:) + & + q(10,7) * var(:,10,:) + q(11,7) * var(:,11,:) ) * idel + dvar(:,8,:) = ( q(1,8) * var(:,1,:) + q(2,8) * var(:,2,:) + & + q(3,8) * var(:,3,:) + q(4,8) * var(:,4,:) + & + q(5,8) * var(:,5,:) + q(6,8) * var(:,6,:) + & + q(7,8) * var(:,7,:) + q(9,8) * var(:,9,:) + & + q(10,8) * var(:,10,:) + q(11,8) * var(:,11,:) + & + q(12,8) * var(:,12,:) ) * idel + jl = 9 + end if + if ( bb(2) == 0 ) then + jr = nj - gsize + else + dvar(:,nj,:) = - ( q(1,1) * var(:,nj,:) + q(2,1) * var(:,nj-1,:) + & + q(3,1) * var(:,nj-2,:) + q(4,1) * var(:,nj-3,:) + & + q(5,1) * var(:,nj-4,:) + q(6,1) * var(:,nj-5,:) + & + q(7,1) * var(:,nj-6,:) + & + q(8,1) * var(:,nj-7,:) ) * idel + dvar(:,nj-1,:) = - ( q(1,2) * var(:,nj,:) + q(3,2) * var(:,nj-2,:) + & + q(4,2) * var(:,nj-3,:) + q(5,2) * var(:,nj-4,:) + & + q(6,2) * var(:,nj-5,:) + q(7,2) * var(:,nj-6,:) + & + q(8,2) * var(:,nj-7,:) ) * idel + dvar(:,nj-2,:) = - ( q(1,3) * var(:,nj,:) + q(2,3) * var(:,nj-1,:) + & + q(4,3) * var(:,nj-3,:) + q(5,3) * var(:,nj-4,:) + & + q(6,3) * var(:,nj-5,:) + q(7,3) * var(:,nj-6,:) + & + q(8,3) * var(:,nj-7,:) ) * idel + dvar(:,nj-3,:) = - ( q(1,4) * var(:,nj,:) + q(2,4) * var(:,nj-1,:) + & + q(3,4) * var(:,nj-2,:) + q(5,4) * var(:,nj-4,:) + & + q(6,4) * var(:,nj-5,:) + q(7,4) * var(:,nj-6,:) + & + q(8,4) * var(:,nj-7,:) ) * idel + dvar(:,nj-4,:) = - ( q(1,5) * var(:,nj,:) + q(2,5) * var(:,nj-1,:) + & + q(3,5) * var(:,nj-2,:) + q(4,5) * var(:,nj-3,:) + & + q(6,5) * var(:,nj-5,:) + q(7,5) * var(:,nj-6,:) + & + q(8,5) * var(:,nj-7,:) + & + q(9,5) * var(:,nj-8,:) ) * idel + dvar(:,nj-5,:) = - ( q(1,6) * var(:,nj,:) + q(2,6) * var(:,nj-1,:) + & + q(3,6) * var(:,nj-2,:) + q(4,6) * var(:,nj-3,:) + & + q(5,6) * var(:,nj-4,:) + q(7,6) * var(:,nj-6,:) + & + q(8,6) * var(:,nj-7,:) + q(9,6) * var(:,nj-8,:) + & + q(10,6) * var(:,nj-9,:) ) * idel + dvar(:,nj-6,:) = - ( q(1,7) * var(:,nj,:) + q(2,7) * var(:,nj-1,:) + & + q(3,7) * var(:,nj-2,:) + q(4,7) * var(:,nj-3,:) + & + q(5,7) * var(:,nj-4,:) + q(6,7) * var(:,nj-5,:) + & + q(8,7) * var(:,nj-7,:) + q(9,7) * var(:,nj-8,:) + & + q(10,7) * var(:,nj-9,:) + & + q(11,7) * var(:,nj-10,:) ) * idel + dvar(:,nj-7,:) = - ( q(1,8) * var(:,nj,:) + q(2,8) * var(:,nj-1,:) + & + q(3,8) * var(:,nj-2,:) + q(4,8) * var(:,nj-3,:) + & + q(5,8) * var(:,nj-4,:) + q(6,8) * var(:,nj-5,:) + & + q(7,8) * var(:,nj-6,:) + q(9,8) * var(:,nj-8,:) + & + q(10,8) * var(:,nj-9,:) + & + q(11,8) * var(:,nj-10,:) + & + q(12,8) * var(:,nj-11,:) ) * idel + jr = nj - 8 + end if + dvar(:,jl:jr,:) = ( a(1) * ( var(:,jl+1:jr+1,:) - & + var(:,jl-1:jr-1,:) ) + & + a(2) * ( var(:,jl+2:jr+2,:) - & + var(:,jl-2:jr-2,:) ) + & + a(3) * ( var(:,jl+3:jr+3,:) - & + var(:,jl-3:jr-3,:) ) + & + a(4) * ( var(:,jl+4:jr+4,:) - & + var(:,jl-4:jr-4,:) ) ) * idel + case (2) direction + if ( bb(1) == 0 ) then + kl = 1 + gsize + else + dvar(:,:,1) = ( q(1,1) * var(:,:,1) + q(2,1) * var(:,:,2) + & + q(3,1) * var(:,:,3) + q(4,1) * var(:,:,4) + & + q(5,1) * var(:,:,5) + q(6,1) * var(:,:,6) + & + q(7,1) * var(:,:,7) + q(8,1) * var(:,:,8) ) * idel + dvar(:,:,2) = ( q(1,2) * var(:,:,1) + q(3,2) * var(:,:,3) + & + q(4,2) * var(:,:,4) + q(5,2) * var(:,:,5) + & + q(6,2) * var(:,:,6) + q(7,2) * var(:,:,7) + & + q(8,2) * var(:,:,8) ) * idel + dvar(:,:,3) = ( q(1,3) * var(:,:,1) + q(2,3) * var(:,:,2) + & + q(4,3) * var(:,:,4) + q(5,3) * var(:,:,5) + & + q(6,3) * var(:,:,6) + q(7,3) * var(:,:,7) + & + q(8,3) * var(:,:,8) ) * idel + dvar(:,:,4) = ( q(1,4) * var(:,:,1) + q(2,4) * var(:,:,2) + & + q(3,4) * var(:,:,3) + q(5,4) * var(:,:,5) + & + q(6,4) * var(:,:,6) + q(7,4) * var(:,:,7) + & + q(8,4) * var(:,:,8) ) * idel + dvar(:,:,5) = ( q(1,5) * var(:,:,1) + q(2,5) * var(:,:,2) + & + q(3,5) * var(:,:,3) + q(4,5) * var(:,:,4) + & + q(6,5) * var(:,:,6) + q(7,5) * var(:,:,7) + & + q(8,5) * var(:,:,8) + q(9,5) * var(:,:,9) ) * idel + dvar(:,:,6) = ( q(1,6) * var(:,:,1) + q(2,6) * var(:,:,2) + & + q(3,6) * var(:,:,3) + q(4,6) * var(:,:,4) + & + q(5,6) * var(:,:,5) + q(7,6) * var(:,:,7) + & + q(8,6) * var(:,:,8) + q(9,6) * var(:,:,9) + & + q(10,6) * var(:,:,10) ) * idel + dvar(:,:,7) = ( q(1,7) * var(:,:,1) + q(2,7) * var(:,:,2) + & + q(3,7) * var(:,:,3) + q(4,7) * var(:,:,4) + & + q(5,7) * var(:,:,5) + q(6,7) * var(:,:,6) + & + q(8,7) * var(:,:,8) + q(9,7) * var(:,:,9) + & + q(10,7) * var(:,:,10) + q(11,7) * var(:,:,11) ) * idel + dvar(:,:,8) = ( q(1,8) * var(:,:,1) + q(2,8) * var(:,:,2) + & + q(3,8) * var(:,:,3) + q(4,8) * var(:,:,4) + & + q(5,8) * var(:,:,5) + q(6,8) * var(:,:,6) + & + q(7,8) * var(:,:,7) + q(9,8) * var(:,:,9) + & + q(10,8) * var(:,:,10) + q(11,8) * var(:,:,11) + & + q(12,8) * var(:,:,12) ) * idel + kl = 9 + end if + if ( bb(2) == 0 ) then + kr = nk - gsize + else + dvar(:,:,nk) = - ( q(1,1) * var(:,:,nk) + q(2,1) * var(:,:,nk-1) + & + q(3,1) * var(:,:,nk-2) + q(4,1) * var(:,:,nk-3) + & + q(5,1) * var(:,:,nk-4) + q(6,1) * var(:,:,nk-5) + & + q(7,1) * var(:,:,nk-6) + & + q(8,1) * var(:,:,nk-7) ) * idel + dvar(:,:,nk-1) = - ( q(1,2) * var(:,:,nk) + q(3,2) * var(:,:,nk-2) + & + q(4,2) * var(:,:,nk-3) + q(5,2) * var(:,:,nk-4) + & + q(6,2) * var(:,:,nk-5) + q(7,2) * var(:,:,nk-6) + & + q(8,2) * var(:,:,nk-7) ) * idel + dvar(:,:,nk-2) = - ( q(1,3) * var(:,:,nk) + q(2,3) * var(:,:,nk-1) + & + q(4,3) * var(:,:,nk-3) + q(5,3) * var(:,:,nk-4) + & + q(6,3) * var(:,:,nk-5) + q(7,3) * var(:,:,nk-6) + & + q(8,3) * var(:,:,nk-7) ) * idel + dvar(:,:,nk-3) = - ( q(1,4) * var(:,:,nk) + q(2,4) * var(:,:,nk-1) + & + q(3,4) * var(:,:,nk-2) + q(5,4) * var(:,:,nk-4) + & + q(6,4) * var(:,:,nk-5) + q(7,4) * var(:,:,nk-6) + & + q(8,4) * var(:,:,nk-7) ) * idel + dvar(:,:,nk-4) = - ( q(1,5) * var(:,:,nk) + q(2,5) * var(:,:,nk-1) + & + q(3,5) * var(:,:,nk-2) + q(4,5) * var(:,:,nk-3) + & + q(6,5) * var(:,:,nk-5) + q(7,5) * var(:,:,nk-6) + & + q(8,5) * var(:,:,nk-7) + & + q(9,5) * var(:,:,nk-8) ) * idel + dvar(:,:,nk-5) = - ( q(1,6) * var(:,:,nk) + q(2,6) * var(:,:,nk-1) + & + q(3,6) * var(:,:,nk-2) + q(4,6) * var(:,:,nk-3) + & + q(5,6) * var(:,:,nk-4) + q(7,6) * var(:,:,nk-6) + & + q(8,6) * var(:,:,nk-7) + q(9,6) * var(:,:,nk-8) + & + q(10,6) * var(:,:,nk-9) ) * idel + dvar(:,:,nk-6) = - ( q(1,7) * var(:,:,nk) + q(2,7) * var(:,:,nk-1) + & + q(3,7) * var(:,:,nk-2) + q(4,7) * var(:,:,nk-3) + & + q(5,7) * var(:,:,nk-4) + q(6,7) * var(:,:,nk-5) + & + q(8,7) * var(:,:,nk-7) + q(9,7) * var(:,:,nk-8) + & + q(10,7) * var(:,:,nk-9) + & + q(11,7) * var(:,:,nk-10) ) * idel + dvar(:,:,nk-7) = - ( q(1,8) * var(:,:,nk) + q(2,8) * var(:,:,nk-1) + & + q(3,8) * var(:,:,nk-2) + q(4,8) * var(:,:,nk-3) + & + q(5,8) * var(:,:,nk-4) + q(6,8) * var(:,:,nk-5) + & + q(7,8) * var(:,:,nk-6) + q(9,8) * var(:,:,nk-8) + & + q(10,8) * var(:,:,nk-9) + & + q(11,8) * var(:,:,nk-10) + & + q(12,8) * var(:,:,nk-11) ) * idel + kr = nk - 8 + end if + dvar(:,:,kl:kr) = ( a(1) * ( var(:,:,kl+1:kr+1) - & + var(:,:,kl-1:kr-1) ) + & + a(2) * ( var(:,:,kl+2:kr+2) - & + var(:,:,kl-2:kr-2) ) + & + a(3) * ( var(:,:,kl+3:kr+3) - & + var(:,:,kl-3:kr-3) ) + & + a(4) * ( var(:,:,kl+4:kr+4) - & + var(:,:,kl-4:kr-4) ) ) * idel + end select direction +end subroutine deriv_gf_8_4 diff --git a/src/call_derivs.c b/src/call_derivs.c index 30c280a..e8c051a 100644 --- a/src/call_derivs.c +++ b/src/call_derivs.c @@ -42,6 +42,15 @@ void DiffGv ( const CCTK_POINTER cctkGH_, const CCTK_INT dir, const CCTK_INT *gsize, const CCTK_REAL *delta, CCTK_REAL *dvar); + void CCTK_FCALL CCTK_FNAME(deriv_gf_8_4)(const CCTK_REAL *var, + const CCTK_INT *ni, + const CCTK_INT *nj, + const CCTK_INT *nk, + const CCTK_INT *dir, + const CCTK_INT *bb, + const CCTK_INT *gsize, + const CCTK_REAL *delta, + CCTK_REAL *dvar); ni = cctk_lsh[0]; nj = cctk_lsh[1]; nk = cctk_lsh[2]; @@ -82,6 +91,10 @@ void DiffGv ( const CCTK_POINTER cctkGH_, const CCTK_INT dir, CCTK_FNAME(deriv_gf_6_3)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar); break; } + case 8: { + CCTK_FNAME(deriv_gf_8_4)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar); + break; + } default: assert(0); } diff --git a/src/make.code.defn b/src/make.code.defn index 70a8873..a5db2d8 100644 --- a/src/make.code.defn +++ b/src/make.code.defn @@ -2,7 +2,7 @@ # $Header$ # Source files in this directory -SRCS = Derivatives_2_1.F90 Derivatives_4_2.F90 Derivatives_6_3.F90 call_derivs.c call_derivs_name.c +SRCS = Derivatives_2_1.F90 Derivatives_4_2.F90 Derivatives_6_3.F90 Derivatives_8_4.F90 call_derivs.c call_derivs_name.c # Subdirectories containing source files SUBDIRS = |