#include "cctk.h" #include "cctk_Functions.h" #include "cctk_Parameters.h" subroutine deriv_gf_4_3 ( var, ni, nj, nk, dir, bb, gsize, delta, dvar ) use All_Coeffs_mod implicit none DECLARE_CCTK_FUNCTIONS DECLARE_CCTK_PARAMETERS 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 call coeffs_1_4_3 ( a, q ) first = .false. end if idel = 1.0_wp / delta if (gsize < 2) call CCTK_WARN (0, "not enough ghostzones") direction: select case (dir) case (0) direction if ( bb(1) == 0 ) then il = 1 + gsize else !$omp parallel workshare 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,:,:) + & q(5,2) * var(5,:,:) + q(6,2) * var(6,:,:) ) * 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,:,:) ) * 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 !$omp end parallel workshare il = 6 end if if ( bb(2) == 0 ) then ir = ni - gsize else !$omp parallel workshare 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,:,:) ) * 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,:,:) ) * 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 !$omp end parallel workshare ir = ni - 5 end if if (il > ir+1) call CCTK_WARN (0, "domain too small") !$omp parallel workshare 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 !$omp end parallel workshare case (1) direction if ( zero_derivs_y /= 0 ) then dvar = zero else if ( bb(1) == 0 ) then jl = 1 + gsize else !$omp parallel workshare 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,:) + & q(5,2) * var(:,5,:) + q(6,2) * var(:,6,:) ) * 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,:) ) * 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 !$omp end parallel workshare jl = 6 end if if ( bb(2) == 0 ) then jr = nj - gsize else !$omp parallel workshare 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,:) ) * 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,:) ) * 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 !$omp end parallel workshare jr = nj - 5 end if if (jl > jr+1) call CCTK_WARN (0, "domain too small") !$omp parallel workshare 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 !$omp end parallel workshare end if case (2) direction if ( zero_derivs_z /= 0 ) then dvar = zero else if ( bb(1) == 0 ) then kl = 1 + gsize else !$omp parallel workshare 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) + & q(5,2) * var(:,:,5) + q(6,2) * var(:,:,6) ) * 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) ) * 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 !$omp end parallel workshare kl = 6 end if if ( bb(2) == 0 ) then kr = nk - gsize else !$omp parallel workshare 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) ) * 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) ) * 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 !$omp end parallel workshare kr = nk - 5 end if if (kl > kr+1) call CCTK_WARN (0, "domain too small") !$omp parallel workshare 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 !$omp end parallel workshare end if end select direction end subroutine deriv_gf_4_3