SCALE-RM
Data Types | Functions/Subroutines
scale_filter Module Reference

module FILTER More...

Functions/Subroutines

subroutine filter_hyperdiff_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, data, order, nite, limiter_sign)
 Hyper diffusion filter 3D. More...
 

Detailed Description

module FILTER

Description
Horizontal filter
Author
Team SCALE

Function/Subroutine Documentation

◆ filter_hyperdiff_3d()

subroutine scale_filter::filter_hyperdiff_3d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), intent(inout)  data,
integer, intent(in)  order,
integer, intent(in)  nite,
real(rp), dimension(ka,ia,ja), intent(in), optional  limiter_sign 
)

Hyper diffusion filter 3D.

Definition at line 169 of file scale_filter.F90.

References scale_prc::prc_abort().

169  use scale_comm_cartesc, only: &
170  comm_vars8, &
171  comm_wait
172  use scale_prc, only: &
173  prc_abort
174  integer, intent(in) :: ka, ks, ke
175  integer, intent(in) :: ia, is, ie
176  integer, intent(in) :: ja, js, je
177 
178  real(RP), intent(inout) :: data(ka,ia,ja)
179 
180  integer, intent(in) :: order
181  integer, intent(in) :: nite
182 
183  real(RP), intent(in), optional :: limiter_sign(ka,ia,ja)
184 
185  real(RP) :: work_data(ia,ja)
186  real(RP) :: work_sign(ia,ja)
187 
188  logical :: flag
189 
190  integer :: k, i, j
191 
192  if ( present(limiter_sign) ) then
193  flag = .true.
194  else
195  flag = .false.
196  end if
197 
198  do k = ks, ke
199 
200  !$omp parallel do
201  do j = js, je
202  do i = is, ie
203  work_data(i,j) = data(k,i,j)
204  end do
205  end do
206  if ( flag ) then
207  !$omp parallel do
208  do j = js, je
209  do i = is, ie
210  work_sign(i,j) = limiter_sign(k,i,j)
211  end do
212  end do
213 
214  call filter_hyperdiff_2d( ia, is, ie, ja, js, je, &
215  work_data(:,:), order, nite, &
216  limiter_sign = work_sign )
217 
218  else
219 
220  call filter_hyperdiff_2d( ia, is, ie, ja, js, je, &
221  work_data(:,:), order, nite )
222  end if
223 
224 
225  !$omp parallel do
226  do j = 1, ja
227  do i = 1, ia
228  data(k,i,j) = work_data(i,j)
229  end do
230  end do
231 
232  end do
233 
234 
235  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
module COMMUNICATION
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
Here is the call graph for this function: