aboutsummaryrefslogtreecommitdiff
path: root/src/Derivatives_4_2.F90
diff options
context:
space:
mode:
authordiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2006-04-10 15:10:43 +0000
committerdiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2006-04-10 15:10:43 +0000
commit0c913f2c81ce7f1b10c89e38299badf23747763e (patch)
tree2aefb05f70fa447791d4523b910b750ee53656da /src/Derivatives_4_2.F90
parent3aaf16cfccc2ad7b57344ddc37fd45e69f6b7ee5 (diff)
First attempt at upwinding SBP operators. Don't use yet. They will change
soon to be more flexible. git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@62 f69c4107-0314-4c4f-9ad4-17e986b73f4a
Diffstat (limited to 'src/Derivatives_4_2.F90')
-rw-r--r--src/Derivatives_4_2.F90292
1 files changed, 292 insertions, 0 deletions
diff --git a/src/Derivatives_4_2.F90 b/src/Derivatives_4_2.F90
index 7fc40b1..85e7ab3 100644
--- a/src/Derivatives_4_2.F90
+++ b/src/Derivatives_4_2.F90
@@ -166,3 +166,295 @@ subroutine deriv_gf_4_2 ( var, ni, nj, nk, dir, bb, gsize, delta, dvar )
end if
end select direction
end subroutine deriv_gf_4_2
+
+subroutine up_deriv_gf_4_2 ( var, ni, nj, nk, dir, bb, gsize, delta, up, 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_INT, intent(IN) :: up
+ CCTK_REAL, dimension(ni,nj,nk), intent(OUT) :: dvar
+
+ CCTK_REAL, dimension(-2:2), save :: a1, a2
+ CCTK_REAL, dimension(-2:2) :: a
+ CCTK_REAL, dimension(6,4), save :: q1, q2
+ CCTK_REAL, dimension(6,4) :: q
+ CCTK_REAL :: idel
+
+ CCTK_INT :: il, ir, jl, jr, kl, kr
+
+ logical, save :: first = .true.
+
+ if ( first ) then
+ a1(-2) = 0.1666666666666666666666666666666666666667_wp
+ a1(-1) = -1.000000000000000000000000000000000000000_wp
+ a1(0) = 0.5000000000000000000000000000000000000000_wp
+ a1(1) = 0.3333333333333333333333333333333333333333_wp
+ a1(2) = 0_wp
+
+ q1(1,1) = -1.176470588235294117647058823529411764706_wp
+ q1(2,1) = 1.264705882352941176470588235294117647059_wp
+ q1(3,1) = 0_wp
+ q1(4,1) = -0.08823529411764705882352941176470588235294_wp
+ q1(5,1) = 0_wp
+ q1(6,1) = 0_wp
+ q1(1,2) = -0.6355932203389830508474576271186440677966_wp
+ q1(2,2) = 0.3389830508474576271186440677966101694915_wp
+ q1(3,2) = 0.2288135593220338983050847457627118644068_wp
+ q1(4,2) = 0.06779661016949152542372881355932203389831_wp
+ q1(5,2) = 0_wp
+ q1(6,2) = 0_wp
+ q1(1,3) = 0.1860465116279069767441860465116279069767_wp
+ q1(2,3) = -1.058139534883720930232558139534883720930_wp
+ q1(3,3) = 0.5581395348837209302325581395348837209302_wp
+ q1(4,3) = 0.3139534883720930232558139534883720930233_wp
+ q1(5,3) = 0_wp
+ q1(6,3) = 0_wp
+ q1(1,4) = 0.03061224489795918367346938775510204081633_wp
+ q1(2,4) = 0.08163265306122448979591836734693877551020_wp
+ q1(3,4) = -0.9285714285714285714285714285714285714286_wp
+ q1(4,4) = 0.4897959183673469387755102040816326530612_wp
+ q1(5,4) = 0.3265306122448979591836734693877551020408_wp
+ q1(6,4) = 0_wp
+
+ a2(-2) = 0_wp
+ a2(-1) = -0.3333333333333333333333333333333333333333_wp
+ a2(0) = -0.5000000000000000000000000000000000000000_wp
+ a2(1) = 1.000000000000000000000000000000000000000_wp
+ a2(2) = -0.1666666666666666666666666666666666666667_wp
+
+ q2(1,1) = -1.647058823529411764705882352941176470588_wp
+ q2(2,1) = 2.205882352941176470588235294117647058824_wp
+ q2(3,1) = -0.4705882352941176470588235294117647058824_wp
+ q2(4,1) = -0.08823529411764705882352941176470588235294_wp
+ q2(5,1) = 0_wp
+ q2(6,1) = 0_wp
+ q2(1,2) = -0.3644067796610169491525423728813559322034_wp
+ q2(2,2) = -0.3389830508474576271186440677966101694915_wp
+ q2(3,2) = 0.7711864406779661016949152542372881355932_wp
+ q2(4,2) = -0.06779661016949152542372881355932203389831_wp
+ q2(5,2) = 0_wp
+ q2(6,2) = 0_wp
+ q2(1,3) = 0_wp
+ q2(2,3) = -0.3139534883720930232558139534883720930233_wp
+ q2(3,3) = -0.5581395348837209302325581395348837209302_wp
+ q2(4,3) = 1.058139534883720930232558139534883720930_wp
+ q2(5,3) = -0.1860465116279069767441860465116279069767_wp
+ q2(6,3) = 0_wp
+ q2(1,4) = 0.03061224489795918367346938775510204081633_wp
+ q2(2,4) = -0.08163265306122448979591836734693877551020_wp
+ q2(3,4) = -0.2755102040816326530612244897959183673469_wp
+ q2(4,4) = -0.4897959183673469387755102040816326530612_wp
+ q2(5,4) = 0.9795918367346938775510204081632653061224_wp
+ q2(6,4) = -0.1632653061224489795918367346938775510204_wp
+
+ first = .false.
+ end if
+
+ idel = 1.0_wp / delta
+
+ upwinding: select case (up)
+ case (-1) upwinding
+ a = a1
+ case (1) upwinding
+ a = a2
+ end select upwinding
+
+ direction: select case (dir)
+ case (0) direction
+ if ( bb(1) == 0 ) then
+ il = 1 + gsize
+ else
+ upwinding_xl: select case (up)
+ case (-1) upwinding_xl
+ q = q1
+ case (1) upwinding_xl
+ q = q2
+ end select upwinding_xl
+
+ dvar(1,:,:) = ( q(1,1) * var(1,:,:) + q(2,1) * var(2,:,:) + &
+ q(3,1) * var(3,:,:) + q(4,1) * var(4,:,:) ) * idel
+ dvar(2,:,:) = ( q(1,2) * var(1,:,:) + q(2,2) * var(2,:,:) + &
+ q(3,2) * var(3,:,:) + q(4,2) * var(4,:,:) ) * 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,:,:) ) * 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,:,:) ) * idel
+ il = 5
+ end if
+ if ( bb(2) == 0 ) then
+ ir = ni - gsize
+ else
+ upwinding_xr: select case (up)
+ case (-1) upwinding_xr
+ q = q2
+ case (1) upwinding_xr
+ q = q1
+ end select upwinding_xr
+
+ 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,:,:) ) * 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,:,:) ) * 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,:,:) ) * 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,:,:) ) * idel
+ ir = ni - 4
+ end if
+ dvar(il:ir,:,:) = ( a(-2) * var(il-2:ir-2,:,:) + &
+ a(-1) * var(il-1:ir-1,:,:) + &
+ a(0) * var(il:ir,:,:) + &
+ a(1) * var(il+1:ir+1,:,:) + &
+ a(2) * var(il+2:ir+2,:,:) ) * idel
+ case (1) direction
+ if ( zero_derivs_y /= 0 ) then
+ dvar = zero
+ else
+ if ( bb(1) == 0 ) then
+ jl = 1 + gsize
+ else
+ upwinding_yl: select case (up)
+ case (-1) upwinding_yl
+ q = q1
+ case (1) upwinding_yl
+ q = q2
+ end select upwinding_yl
+
+ dvar(:,1,:) = ( q(1,1) * var(:,1,:) + q(2,1) * var(:,2,:) + &
+ q(3,1) * var(:,3,:) + q(4,1) * var(:,4,:) ) * idel
+ dvar(:,2,:) = ( q(1,2) * var(:,1,:) + q(2,2) * var(:,2,:) + &
+ q(3,2) * var(:,3,:) + q(4,2) * var(:,4,:) ) * 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,:) ) * 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,:) ) * idel
+ jl = 5
+ end if
+ if ( bb(2) == 0 ) then
+ jr = nj - gsize
+ else
+ upwinding_yr: select case (up)
+ case (-1) upwinding_yr
+ q = q2
+ case (1) upwinding_yr
+ q = q1
+ end select upwinding_yr
+
+ 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,:) ) * 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,:) ) * 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,:) ) * 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,:) ) * idel
+ jr = nj - 4
+ end if
+ dvar(:,jl:jr,:) = ( a(-2) * var(:,jl-2:jr-2,:) + &
+ a(-1) * var(:,jl-1:jr-1,:) + &
+ a(0) * var(:,jl:jr,:) + &
+ a(1) * var(:,jl+1:jr+1,:) + &
+ a(2) * var(:,jl+2:jr+2,:) ) * idel
+ end if
+ case (2) direction
+ if ( zero_derivs_z /= 0 ) then
+ dvar = zero
+ else
+ if ( bb(1) == 0 ) then
+ kl = 1 + gsize
+ else
+ upwinding_zl: select case (up)
+ case (-1) upwinding_zl
+ q = q1
+ case (1) upwinding_zl
+ q = q2
+ end select upwinding_zl
+
+ dvar(:,:,1) = ( q(1,1) * var(:,:,1) + q(2,1) * var(:,:,2) + &
+ q(3,1) * var(:,:,3) + q(4,1) * var(:,:,4) ) * idel
+ dvar(:,:,2) = ( q(1,2) * var(:,:,1) + q(2,2) * var(:,:,2) + &
+ q(3,2) * var(:,:,3) + q(4,2) * var(:,:,4) ) * 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) ) * 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) ) * idel
+ kl = 5
+ end if
+ if ( bb(2) == 0 ) then
+ kr = nk - gsize
+ else
+ upwinding_zr: select case (up)
+ case (-1) upwinding_zr
+ q = q2
+ case (1) upwinding_zr
+ q = q1
+ end select upwinding_zr
+
+ 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) ) * 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) ) * 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) ) * 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) ) * idel
+ kl = nk - 4
+ end if
+ dvar(:,:,kl:kr) = ( a(-2) * var(:,:,kl-2:kr-2) + &
+ a(-1) * var(:,:,kl-1:kr-1) + &
+ a(0) * var(:,:,kl:kr) + &
+ a(1) * var(:,:,kl+1:kr+1) + &
+ a(2) * var(:,:,kl+2:kr+2) ) * idel
+ end if
+ end select direction
+end subroutine up_deriv_gf_4_2