SCALE-RM
scale_atmos_diagnostic_cartesC.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
20  !-----------------------------------------------------------------------------
21  implicit none
22  private
23  !-----------------------------------------------------------------------------
24  !
25  !++ Public procedure
26  !
28 
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public parameters & variables
32  !
33  !-----------------------------------------------------------------------------
34  !
35  !++ Private procedure
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private parameters & variables
40  !
41  !-----------------------------------------------------------------------------
42 contains
43  !-----------------------------------------------------------------------------
48  KA, KS, KE, &
49  IA, IS, IE, &
50  JA, JS, JE, &
51  DENS, &
52  MOMZ, &
53  MOMX, &
54  MOMY, &
55  W, &
56  U, &
57  V )
58  use scale_prc_cartesc, only: &
59  prc_twod
60  use scale_comm_cartesc, only: &
61  comm_vars8, &
62  comm_wait
68  i_xyz, &
69  i_xyw
70  integer, intent(in) :: ka, ks, ke
71  integer, intent(in) :: ia, is, ie
72  integer, intent(in) :: ja, js, je
73 
74  real(rp), intent(in) :: dens(ka,ia,ja)
75  real(rp), intent(in) :: momz(ka,ia,ja)
76  real(rp), intent(in) :: momx(ka,ia,ja)
77  real(rp), intent(in) :: momy(ka,ia,ja)
78 
79  real(rp), intent(out) :: w (ka,ia,ja)
80  real(rp), intent(out) :: u (ka,ia,ja)
81  real(rp), intent(out) :: v (ka,ia,ja)
82 
83  real(rp) :: momws
84 
85  integer :: k, i, j
86  !---------------------------------------------------------------------------
87 
88  !$acc data copyin(DENS, MOMZ, MOMY, MOMX, GSQRT, J23G, J13G) copyout(W, U, V)
89 
90  ! Note: W(KS,:,:) is filled because the values at i=1 or j=1 are not calculated below.
91 !OCL XFILL
92  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
93  !$acc kernels
94  do j = js, je
95  do i = is, ie
96  do k = ks, ke-1
97  w(k,i,j) = 0.5_rp * ( momz(k-1,i,j)+momz(k,i,j) ) / dens(k,i,j)
98  enddo
99  enddo
100  enddo
101  !$acc end kernels
102  if ( prc_twod ) then
103 !OCL XFILL
104  !$omp parallel do OMP_SCHEDULE_ &
105  !$omp private(j,momws)
106  !$acc kernels
107  do j = max(2,js), je
108  ! at KS+1/2
109  momws = momz(ks,is,j) &
110  + j23g(ks,is,j,i_xyw) * ( momy(ks,is,j) + momy(ks,is,j-1) + momy(ks+1,is,j) + momy(ks+1,is,j-1) ) &
111  * 0.25_rp / gsqrt(ks,is,j,i_xyw)
112  ! at KS
113  ! momws at the surface is assumed to be zero
114  w(ks,is,j) = ( momws * 0.5_rp &
115  - j23g(ks,is,j,i_xyz) * ( momy(ks,is,j) + momy(ks,is,j-1) ) &
116  * 0.5_rp / gsqrt(ks,is,j,i_xyz) &
117  ) / dens(ks,is,j)
118  enddo
119  !$acc end kernels
120  else
121 !OCL XFILL
122  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
123  !$omp private(i,j,momws)
124  !$acc kernels
125  do j = max(2,js), je
126  do i = max(2,is), ie
127  ! at KS+1/2
128  momws = momz(ks,i,j) &
129  + ( j13g(ks,i,j,i_xyw) * ( momx(ks,i,j) + momx(ks,i-1,j) + momx(ks+1,i,j) + momx(ks+1,i-1,j) ) &
130  + j23g(ks,i,j,i_xyw) * ( momy(ks,i,j) + momy(ks,i,j-1) + momy(ks+1,i,j) + momy(ks+1,i,j-1) ) ) &
131  * 0.25_rp / gsqrt(ks,i,j,i_xyw)
132  ! at KS
133  ! momws at the surface is assumed to be zero
134  w(ks,i,j) = ( momws * 0.5_rp &
135  - ( j13g(ks,i,j,i_xyz) * ( momx(ks,i,j) + momx(ks,i-1,j) ) &
136  + j23g(ks,i,j,i_xyz) * ( momy(ks,i,j) + momy(ks,i,j-1) ) ) &
137  * 0.5_rp / gsqrt(ks,i,j,i_xyz) &
138  ) / dens(ks,i,j)
139  enddo
140  enddo
141  !$acc end kernels
142  end if
143 !OCL XFILL
144  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
145  !$acc kernels
146  do j = js, je
147  do i = is, ie
148  w(ke,i,j) = 0.5_rp * ( momz(ke-1,i,j) ) / dens(ke,i,j)
149  enddo
150  enddo
151  !$acc end kernels
152 
153  if ( prc_twod ) then
154 !OCL XFILL
155  !$omp parallel do private(j,k) OMP_SCHEDULE_
156  !$acc kernels
157  do j = js, je
158  do k = ks, ke
159  u(k,is,j) = momx(k,is,j) / dens(k,is,j)
160  enddo
161  enddo
162  !$acc end kernels
163  else
164 !OCL XFILL
165  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
166  !$acc kernels
167  do j = js, je
168  do i = max(2,is), ie
169  do k = ks, ke
170  u(k,i,j) = 0.5_rp * ( momx(k,i-1,j)+momx(k,i,j) ) / dens(k,i,j)
171  enddo
172  enddo
173  enddo
174  !$acc end kernels
175 !OCL XFILL
176  !$omp parallel do private(j,k) OMP_SCHEDULE_ collapse(2)
177  !$acc kernels
178  do j = js, je
179  do k = ks, ke
180  u(k,1,j) = momx(k,1,j) / dens(k,1,j)
181  enddo
182  enddo
183  !$acc end kernels
184  end if
185 
186  !OCL XFILL
187  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
188  !$acc kernels
189  do j = max(2,js), je
190  do i = is, ie
191  do k = ks, ke
192  v(k,i,j) = 0.5_rp * ( momy(k,i,j-1)+momy(k,i,j) ) / dens(k,i,j)
193  enddo
194  enddo
195  enddo
196  !$acc end kernels
197 !OCL XFILL
198  !$omp parallel do private(i,k) OMP_SCHEDULE_ collapse(2)
199  !$acc kernels
200  do i = is, ie
201  do k = ks, ke
202  v(k,i,1) = momy(k,i,1) / dens(k,i,1)
203  enddo
204  enddo
205  !$acc end kernels
206 
207  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
208  !$acc kernels
209  do j = js, je
210  do i = is, ie
211  !$acc loop seq
212  do k = 1, ks-1
213  w(k,i,j) = w(ks,i,j)
214  u(k,i,j) = u(ks,i,j)
215  v(k,i,j) = v(ks,i,j)
216  end do
217  !$acc loop seq
218  do k = ke+1, ka
219  w(k,i,j) = w(ke,i,j)
220  u(k,i,j) = u(ke,i,j)
221  v(k,i,j) = v(ke,i,j)
222  end do
223  enddo
224  enddo
225  !$acc end kernels
226 
227  call comm_vars8( w(:,:,:), 1 )
228  call comm_vars8( u(:,:,:), 2 )
229  call comm_vars8( v(:,:,:), 3 )
230  call comm_wait ( w(:,:,:), 1, .false. )
231  call comm_wait ( u(:,:,:), 2, .false. )
232  call comm_wait ( v(:,:,:), 3, .false. )
233 
234  !$acc end data
235 
236  return
237  end subroutine atmos_diagnostic_cartesc_get_vel
238 
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
scale_atmos_grid_cartesc_index::i_xyz
integer, public i_xyz
Definition: scale_atmos_grid_cartesC_index.F90:91
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_diagnostic_cartesc::atmos_diagnostic_cartesc_get_vel
subroutine, public atmos_diagnostic_cartesc_get_vel(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, MOMZ, MOMX, MOMY, W, U, V)
ATMOS_DIAGNOSTIC_CARTESC_get_vel W, U, V.
Definition: scale_atmos_diagnostic_cartesC.F90:58
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_gsqrt
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_gsqrt
transformation metrics from Z to Xi, {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:38
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j23g
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j23g
(2,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:40
scale_atmos_diagnostic_cartesc
module atmosphere / diagnostic / CartesianC
Definition: scale_atmos_diagnostic_cartesC.F90:12
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j13g
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j13g
(1,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:39
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:56
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:92