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_comm_cartesc, only: &
59  comm_vars8, &
60  comm_wait
66  i_xyz, &
67  i_xyw
68  integer, intent(in) :: KA, KS, KE
69  integer, intent(in) :: IA, IS, IE
70  integer, intent(in) :: JA, JS, JE
71 
72  real(RP), intent(in) :: DENS(ka,ia,ja)
73  real(RP), intent(in) :: MOMZ(ka,ia,ja)
74  real(RP), intent(in) :: MOMX(ka,ia,ja)
75  real(RP), intent(in) :: MOMY(ka,ia,ja)
76 
77  real(RP), intent(out) :: W (ka,ia,ja)
78  real(RP), intent(out) :: U (ka,ia,ja)
79  real(RP), intent(out) :: V (ka,ia,ja)
80 
81  real(RP) :: momws
82 
83  integer :: k, i, j
84  !---------------------------------------------------------------------------
85 
86  ! Note: W(KS,:,:) is filled because the values at i=1 or j=1 are not calculated below.
87 !OCL XFILL
88  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
89  do j = js, je
90  do i = is, ie
91  do k = ks, ke-1
92  w(k,i,j) = 0.5_rp * ( momz(k-1,i,j)+momz(k,i,j) ) / dens(k,i,j)
93  enddo
94  enddo
95  enddo
96 !OCL XFILL
97  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
98  !$omp private(i,j,momws)
99  do j = max(2,js), je
100  do i = max(2,is), ie
101  ! at KS+1/2
102  momws = momz(ks,i,j) &
103  + ( 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) ) &
104  + 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) ) ) &
105  * 0.25_rp / gsqrt(ks,i,j,i_xyw)
106  ! at KS
107  ! momws at the surface is assumed to be zero
108  w(ks,i,j) = momws * 0.5_rp &
109  - ( j13g(ks,i,j,i_xyz) * ( momx(ks,i,j) + momx(ks,i-1,j) ) &
110  + j23g(ks,i,j,i_xyz) * ( momy(ks,i,j) + momy(ks,i,j-1) ) ) &
111  * 0.5_rp / ( dens(ks,i,j) * gsqrt(ks,i,j,i_xyz) )
112  enddo
113  enddo
114 !OCL XFILL
115  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
116  do j = js, je
117  do i = is, ie
118  w(ke,i,j) = 0.5_rp * ( momz(ke-1,i,j) ) / dens(ke,i,j)
119  enddo
120  enddo
121 
122 !OCL XFILL
123  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
124  do j = js, je
125  do i = max(2,is), ie
126  do k = ks, ke
127  u(k,i,j) = 0.5_rp * ( momx(k,i-1,j)+momx(k,i,j) ) / dens(k,i,j)
128  enddo
129  enddo
130  enddo
131 !OCL XFILL
132  !$omp parallel do private(j,k) OMP_SCHEDULE_ collapse(2)
133  do j = js, je
134  do k = ks, ke
135  u(k,1,j) = momx(k,1,j) / dens(k,1,j)
136  enddo
137  enddo
138 
139  !OCL XFILL
140  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
141  do j = max(2,js), je
142  do i = is, ie
143  do k = ks, ke
144  v(k,i,j) = 0.5_rp * ( momy(k,i,j-1)+momy(k,i,j) ) / dens(k,i,j)
145  enddo
146  enddo
147  enddo
148 !OCL XFILL
149  !$omp parallel do private(i,k) OMP_SCHEDULE_ collapse(2)
150  do i = is, ie
151  do k = ks, ke
152  v(k,i,1) = momy(k,i,1) / dens(k,i,1)
153  enddo
154  enddo
155 
156  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
157  do j = js, je
158  do i = is, ie
159  w( 1:ks-1,i,j) = w(ks,i,j)
160  u( 1:ks-1,i,j) = u(ks,i,j)
161  v( 1:ks-1,i,j) = v(ks,i,j)
162  w(ke+1:ka, i,j) = w(ke,i,j)
163  u(ke+1:ka, i,j) = u(ke,i,j)
164  v(ke+1:ka, i,j) = v(ke,i,j)
165  enddo
166  enddo
167 
168  call comm_vars8( w(:,:,:), 1 )
169  call comm_vars8( u(:,:,:), 2 )
170  call comm_vars8( v(:,:,:), 3 )
171  call comm_wait ( w(:,:,:), 1, .false. )
172  call comm_wait ( u(:,:,:), 2, .false. )
173  call comm_wait ( v(:,:,:), 3, .false. )
174 
175  return
176  end subroutine atmos_diagnostic_cartesc_get_vel
177 
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_gsqrt
transformation metrics from Z to Xi, {G}^1/2
module Atmosphere Grid CartesianC metirc
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.
module COMMUNICATION
module atmosphere / grid / cartesC index
module atmosphere / diagnostic / CartesianC
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j23g
(2,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j13g
(1,3) element of Jacobian matrix * {G}^1/2
module profiler
Definition: scale_prof.F90:11
module PRECISION
module STDIO
Definition: scale_io.F90:10