aboutsummaryrefslogtreecommitdiff
path: root/src/Derivatives_4_3_min_err_coeff.F90
diff options
context:
space:
mode:
Diffstat (limited to 'src/Derivatives_4_3_min_err_coeff.F90')
-rw-r--r--src/Derivatives_4_3_min_err_coeff.F90235
1 files changed, 235 insertions, 0 deletions
diff --git a/src/Derivatives_4_3_min_err_coeff.F90 b/src/Derivatives_4_3_min_err_coeff.F90
new file mode 100644
index 0000000..938c4bc
--- /dev/null
+++ b/src/Derivatives_4_3_min_err_coeff.F90
@@ -0,0 +1,235 @@
+#include "cctk.h"
+#include "cctk_Functions.h"
+#include "cctk_Parameters.h"
+
+subroutine deriv_gf_4_3_opt ( 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(2), save :: a
+ CCTK_REAL, dimension(7,5), save :: q
+ CCTK_REAL :: idel
+
+ CCTK_INT :: il, ir, jl, jr, kl, kr
+
+ logical, save :: first = .true.
+
+ if ( first ) then
+ a(1) = 2.0_wp/3.0_wp; a(2) = -1.0_wp/12.0_wp
+
+ q(1,1) = -2.09329763466349871588733_wp
+ q(2,1) = 4.0398572053206615302160_wp
+ q(3,1) = -3.0597858079809922953240_wp
+ q(4,1) = 1.37319053865399486354933_wp
+ q(5,1) = -0.25996430133016538255400_wp
+ q(6,1) = 0_wp
+ q(7,1) = 0_wp
+ q(1,2) = -0.31641585285940445272297_wp
+ q(2,2) = -0.53930788973980422327388_wp
+ q(3,2) = 0.98517732028644343383297_wp
+ q(4,2) = -0.05264665989297578146709_wp
+ q(5,2) = -0.113807251750624235013258_wp
+ q(6,2) = 0.039879767889849911803103_wp
+ q(7,2) = -0.0028794339334846531588787_wp
+ q(1,3) = 0.13026916185021164524452_wp
+ q(2,3) = -0.87966858995059249256890_wp
+ q(3,3) = 0.38609640961100070000134_wp
+ q(4,3) = 0.31358369072435588745988_wp
+ q(5,3) = 0.085318941913678384633511_wp
+ q(6,3) = -0.039046615792734640274641_wp
+ q(7,3) = 0.0034470016440805155042908_wp
+ q(1,4) = -0.01724512193824647912172_wp
+ q(2,4) = 0.16272288227127504381134_wp
+ q(3,4) = -0.81349810248648813029217_wp
+ q(4,4) = 0.13833269266479833215645_wp
+ q(5,4) = 0.59743854328548053399616_wp
+ q(6,4) = -0.066026434346299887619324_wp
+ q(7,4) = -0.0017244594505194129307249_wp
+ q(1,5) = -0.00883569468552192965061_wp
+ q(2,5) = 0.03056074759203203857284_wp
+ q(3,5) = 0.05021168274530854232278_wp
+ q(4,5) = -0.66307364652444929534068_wp
+ q(5,5) = 0.014878787464005191116088_wp
+ q(6,5) = 0.65882706381707471953820_wp
+ q(7,5) = -0.082568940408449266558615_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,:,:) ) * idel
+ dvar(2,:,:) = ( q(1,2) * var(1,:,:) + q(2,2) * var(2,:,:) + &
+ 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,:,:) ) * idel
+ dvar(3,:,:) = ( q(1,3) * var(1,:,:) + q(2,3) * var(2,:,:) + &
+ q(3,3) * var(3,:,:) + q(4,3) * var(4,:,:) + &
+ q(5,3) * var(5,:,:) + q(6,3) * var(6,:,:) + &
+ q(7,3) * var(7,:,:) ) * idel
+ dvar(4,:,:) = ( q(1,4) * var(1,:,:) + q(2,4) * var(2,:,:) + &
+ q(3,4) * var(3,:,:) + q(4,4) * var(4,:,:) + &
+ q(5,4) * var(5,:,:) + q(6,4) * var(6,:,:) + &
+ q(7,4) * var(7,:,:) ) * idel
+ dvar(5,:,:) = ( q(1,5) * var(1,:,:) + q(2,5) * var(2,:,:) + &
+ q(3,5) * var(3,:,:) + q(4,5) * var(4,:,:) + &
+ q(5,5) * var(5,:,:) + q(6,5) * var(6,:,:) + &
+ q(7,5) * var(7,:,:) ) * idel
+ il = 6
+ end if
+ if ( bb(2) == 0 ) then
+ ir = ni - gsize
+ else
+ 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(5,5) * var(ni-4,:,:) + q(6,5) * var(ni-5,:,:) + &
+ q(7,5) * var(ni-6,:,:) ) * idel
+ dvar(ni-3,:,:) = - ( q(1,4) * var(ni,:,:) + q(2,4) * var(ni-1,:,:) + &
+ q(3,4) * var(ni-2,:,:) + q(4,4) * var(ni-3,:,:) + &
+ q(5,4) * var(ni-4,:,:) + q(6,4) * var(ni-5,:,:) + &
+ q(7,4) * var(ni-6,:,:) ) * idel
+ dvar(ni-2,:,:) = - ( q(1,3) * var(ni,:,:) + q(2,3) * var(ni-1,:,:) + &
+ q(3,3) * var(ni-2,:,:) + q(4,3) * var(ni-3,:,:) + &
+ q(5,3) * var(ni-4,:,:) + q(6,3) * var(ni-5,:,:) + &
+ q(7,3) * var(ni-6,:,:) ) * idel
+ dvar(ni-1,:,:) = - ( q(1,2) * var(ni,:,:) + q(2,2) * var(ni-1,:,:) + &
+ 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,:,:) ) * idel
+ 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,:,:) ) * idel
+ ir = ni - 5
+ 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,:,:) ) ) * 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,:) ) * idel
+ dvar(:,2,:) = ( q(1,2) * var(:,1,:) + q(2,2) * var(:,2,:) + &
+ 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,:) ) * idel
+ dvar(:,3,:) = ( q(1,3) * var(:,1,:) + q(2,3) * var(:,2,:) + &
+ q(3,3) * var(:,3,:) + q(4,3) * var(:,4,:) + &
+ q(5,3) * var(:,5,:) + q(6,3) * var(:,6,:) + &
+ q(7,3) * var(:,7,:) ) * idel
+ dvar(:,4,:) = ( q(1,4) * var(:,1,:) + q(2,4) * var(:,2,:) + &
+ q(3,4) * var(:,3,:) + q(4,4) * var(:,4,:) + &
+ q(5,4) * var(:,5,:) + q(6,4) * var(:,6,:) + &
+ q(7,4) * var(:,7,:) ) * idel
+ dvar(:,5,:) = ( q(1,5) * var(:,1,:) + q(2,5) * var(:,2,:) + &
+ q(3,5) * var(:,3,:) + q(4,5) * var(:,4,:) + &
+ q(5,5) * var(:,5,:) + q(6,5) * var(:,6,:) + &
+ q(7,5) * var(:,7,:) ) * idel
+ jl = 6
+ end if
+ if ( bb(2) == 0 ) then
+ jr = nj - gsize
+ else
+ 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(5,5) * var(:,nj-4,:) + q(6,5) * var(:,nj-5,:) + &
+ q(7,5) * var(:,nj-6,:) ) * idel
+ dvar(:,nj-3,:) = - ( q(1,4) * var(:,nj,:) + q(2,4) * var(:,nj-1,:) + &
+ q(3,4) * var(:,nj-2,:) + q(4,4) * var(:,nj-3,:) + &
+ q(5,4) * var(:,nj-4,:) + q(6,4) * var(:,nj-5,:) + &
+ q(7,4) * var(:,nj-6,:) ) * idel
+ dvar(:,nj-2,:) = - ( q(1,3) * var(:,nj,:) + q(2,3) * var(:,nj-1,:) + &
+ q(3,3) * var(:,nj-2,:) + q(4,3) * var(:,nj-3,:) + &
+ q(5,3) * var(:,nj-4,:) + q(6,3) * var(:,nj-5,:) + &
+ q(7,3) * var(:,nj-6,:) ) * idel
+ dvar(:,nj-1,:) = - ( q(1,2) * var(:,nj,:) + q(2,2) * var(:,nj-1,:) + &
+ 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,:) ) * idel
+ 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,:) ) * idel
+ jr = nj - 5
+ 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,:) ) ) * 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) ) * idel
+ dvar(:,:,2) = ( q(1,2) * var(:,:,1) + q(2,2) * var(:,:,2) + &
+ 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) ) * idel
+ dvar(:,:,3) = ( q(1,3) * var(:,:,1) + q(2,3) * var(:,:,2) + &
+ q(3,3) * var(:,:,3) + q(4,3) * var(:,:,4) + &
+ q(5,3) * var(:,:,5) + q(6,3) * var(:,:,6) + &
+ q(7,3) * var(:,:,7) ) * idel
+ dvar(:,:,4) = ( q(1,4) * var(:,:,1) + q(2,4) * var(:,:,2) + &
+ q(3,4) * var(:,:,3) + q(4,4) * var(:,:,4) + &
+ q(5,4) * var(:,:,5) + q(6,4) * var(:,:,6) + &
+ q(7,4) * var(:,:,7) ) * idel
+ dvar(:,:,5) = ( q(1,5) * var(:,:,1) + q(2,5) * var(:,:,2) + &
+ q(3,5) * var(:,:,3) + q(4,5) * var(:,:,4) + &
+ q(5,5) * var(:,:,5) + q(6,5) * var(:,:,6) + &
+ q(7,5) * var(:,:,7) ) * idel
+ kl = 6
+ end if
+ if ( bb(2) == 0 ) then
+ kr = nk - gsize
+ else
+ 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(5,5) * var(:,:,nk-4) + q(6,5) * var(:,:,nk-5) + &
+ q(7,5) * var(:,:,nk-6) ) * idel
+ dvar(:,:,nk-3) = - ( q(1,4) * var(:,:,nk) + q(2,4) * var(:,:,nk-1) + &
+ q(3,4) * var(:,:,nk-2) + q(4,4) * var(:,:,nk-3) + &
+ q(5,4) * var(:,:,nk-4) + q(6,4) * var(:,:,nk-5) + &
+ q(7,4) * var(:,:,nk-6) ) * idel
+ dvar(:,:,nk-2) = - ( q(1,3) * var(:,:,nk) + q(2,3) * var(:,:,nk-1) + &
+ q(3,3) * var(:,:,nk-2) + q(4,3) * var(:,:,nk-3) + &
+ q(5,3) * var(:,:,nk-4) + q(6,3) * var(:,:,nk-5) + &
+ q(7,3) * var(:,:,nk-6) ) * idel
+ dvar(:,:,nk-1) = - ( q(1,2) * var(:,:,nk) + q(2,2) * var(:,:,nk-1) + &
+ 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) ) * idel
+ 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) ) * idel
+ kr = nk - 5
+ 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) ) ) * idel
+ end select direction
+end subroutine deriv_gf_4_3_opt