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 252 of file scale_filter.F90.

252  use scale_comm_cartesc, only: &
253  comm_vars8, &
254  comm_wait
255  use scale_prc, only: &
256  prc_abort
257  integer, intent(in) :: KA, KS, KE
258  integer, intent(in) :: IA, IS, IE
259  integer, intent(in) :: JA, JS, JE
260 
261  real(RP), intent(inout) :: data(KA,IA,JA)
262 
263  integer, intent(in) :: order
264  integer, intent(in) :: nite
265 
266  real(RP), intent(in), optional :: limiter_sign(KA,IA,JA)
267 
268  real(RP) :: work_data(IA,JA)
269  real(RP) :: work_sign(IA,JA)
270 
271  logical :: flag
272 
273  integer :: k, i, j
274 
275  if ( present(limiter_sign) ) then
276  flag = .true.
277  else
278  flag = .false.
279  end if
280 
281  !$acc data copy(data) create(work_data)
282  !$acc data copyin(limiter_sign) create(work_sign) if(flag)
283 
284  do k = ks, ke
285 
286  !$omp parallel do
287  !$acc kernels
288  do j = js, je
289  do i = is, ie
290  work_data(i,j) = data(k,i,j)
291  end do
292  end do
293  !$acc end kernels
294  if ( flag ) then
295  !$omp parallel do
296  !$acc kernels
297  do j = js, je
298  do i = is, ie
299  work_sign(i,j) = limiter_sign(k,i,j)
300  end do
301  end do
302  !$acc end kernels
303 
304  call filter_hyperdiff_2d( ia, is, ie, ja, js, je, &
305  work_data(:,:), order, nite, &
306  limiter_sign = work_sign )
307 
308  else
309 
310  call filter_hyperdiff_2d( ia, is, ie, ja, js, je, &
311  work_data(:,:), order, nite )
312  end if
313 
314 
315  !$omp parallel do
316  !$acc kernels
317  do j = 1, ja
318  do i = 1, ia
319  data(k,i,j) = work_data(i,j)
320  end do
321  end do
322  !$acc end kernels
323 
324  end do
325 
326  !$acc end data
327  !$acc end data
328 
329  return

References scale_prc::prc_abort().

Here is the call graph for this function:
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56