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  ! Note: W(KS,:,:) is filled because the values at i=1 or j=1 are not calculated below.
89 !OCL XFILL
90  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
91  do j = js, je
92  do i = is, ie
93  do k = ks, ke-1
94  w(k,i,j) = 0.5_rp * ( momz(k-1,i,j)+momz(k,i,j) ) / dens(k,i,j)
95  enddo
96  enddo
97  enddo
98  if ( prc_twod ) then
99 !OCL XFILL
100  !$omp parallel do OMP_SCHEDULE_ &
101  !$omp private(j,momws)
102  do j = max(2,js), je
103  ! at KS+1/2
104  momws = momz(ks,is,j) &
105  + 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) ) &
106  * 0.25_rp / gsqrt(ks,is,j,i_xyw)
107  ! at KS
108  ! momws at the surface is assumed to be zero
109  w(ks,is,j) = ( momws * 0.5_rp &
110  - j23g(ks,is,j,i_xyz) * ( momy(ks,is,j) + momy(ks,is,j-1) ) &
111  * 0.5_rp / gsqrt(ks,is,j,i_xyz) &
112  ) / dens(ks,is,j)
113  enddo
114  else
115 !OCL XFILL
116  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
117  !$omp private(i,j,momws)
118  do j = max(2,js), je
119  do i = max(2,is), ie
120  ! at KS+1/2
121  momws = momz(ks,i,j) &
122  + ( 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) ) &
123  + 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) ) ) &
124  * 0.25_rp / gsqrt(ks,i,j,i_xyw)
125  ! at KS
126  ! momws at the surface is assumed to be zero
127  w(ks,i,j) = ( momws * 0.5_rp &
128  - ( j13g(ks,i,j,i_xyz) * ( momx(ks,i,j) + momx(ks,i-1,j) ) &
129  + j23g(ks,i,j,i_xyz) * ( momy(ks,i,j) + momy(ks,i,j-1) ) ) &
130  * 0.5_rp / gsqrt(ks,i,j,i_xyz) &
131  ) / dens(ks,i,j)
132  enddo
133  enddo
134  end if
135 !OCL XFILL
136  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
137  do j = js, je
138  do i = is, ie
139  w(ke,i,j) = 0.5_rp * ( momz(ke-1,i,j) ) / dens(ke,i,j)
140  enddo
141  enddo
142 
143  if ( prc_twod ) then
144 !OCL XFILL
145  !$omp parallel do private(j,k) OMP_SCHEDULE_
146  do j = js, je
147  do k = ks, ke
148  u(k,is,j) = momx(k,is,j) / dens(k,is,j)
149  enddo
150  enddo
151  else
152 !OCL XFILL
153  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
154  do j = js, je
155  do i = max(2,is), ie
156  do k = ks, ke
157  u(k,i,j) = 0.5_rp * ( momx(k,i-1,j)+momx(k,i,j) ) / dens(k,i,j)
158  enddo
159  enddo
160  enddo
161 !OCL XFILL
162  !$omp parallel do private(j,k) OMP_SCHEDULE_ collapse(2)
163  do j = js, je
164  do k = ks, ke
165  u(k,1,j) = momx(k,1,j) / dens(k,1,j)
166  enddo
167  enddo
168  end if
169 
170  !OCL XFILL
171  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
172  do j = max(2,js), je
173  do i = is, ie
174  do k = ks, ke
175  v(k,i,j) = 0.5_rp * ( momy(k,i,j-1)+momy(k,i,j) ) / dens(k,i,j)
176  enddo
177  enddo
178  enddo
179 !OCL XFILL
180  !$omp parallel do private(i,k) OMP_SCHEDULE_ collapse(2)
181  do i = is, ie
182  do k = ks, ke
183  v(k,i,1) = momy(k,i,1) / dens(k,i,1)
184  enddo
185  enddo
186 
187  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
188  do j = js, je
189  do i = is, ie
190  w( 1:ks-1,i,j) = w(ks,i,j)
191  u( 1:ks-1,i,j) = u(ks,i,j)
192  v( 1:ks-1,i,j) = v(ks,i,j)
193  w(ke+1:ka, i,j) = w(ke,i,j)
194  u(ke+1:ka, i,j) = u(ke,i,j)
195  v(ke+1:ka, i,j) = v(ke,i,j)
196  enddo
197  enddo
198 
199  call comm_vars8( w(:,:,:), 1 )
200  call comm_vars8( u(:,:,:), 2 )
201  call comm_vars8( v(:,:,:), 3 )
202  call comm_wait ( w(:,:,:), 1, .false. )
203  call comm_wait ( u(:,:,:), 2, .false. )
204  call comm_wait ( v(:,:,:), 3, .false. )
205 
206  return
207  end subroutine atmos_diagnostic_cartesc_get_vel
208 
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:90
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:37
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:39
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:38
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:55
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:91