SCALE-RM
Functions/Subroutines
scale_atmos_diagnostic Module Reference

module ATMOSPHERE / Diagnostic More...

Functions/Subroutines

subroutine, public atmos_diagnostic_get (POTT, TEMP, PRES, PHYD, W, U, V, N2, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC)
 

Detailed Description

module ATMOSPHERE / Diagnostic

Description
Calculate diagnostic variables
Author
Team SCALE
History
  • 2017-02-24 (S.Nishizawa) [new]

Function/Subroutine Documentation

◆ atmos_diagnostic_get()

subroutine, public scale_atmos_diagnostic::atmos_diagnostic_get ( real(rp), dimension(ka,ia,ja), intent(out)  POTT,
real(rp), dimension(ka,ia,ja), intent(out)  TEMP,
real(rp), dimension(ka,ia,ja), intent(out)  PRES,
real(rp), dimension(ka,ia,ja), intent(out)  PHYD,
real(rp), dimension (ka,ia,ja), intent(out)  W,
real(rp), dimension (ka,ia,ja), intent(out)  U,
real(rp), dimension (ka,ia,ja), intent(out)  V,
real(rp), dimension (ka,ia,ja), intent(out)  N2,
real(rp), dimension(ka,ia,ja), intent(in)  DENS,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC 
)

Definition at line 64 of file scale_atmos_sub_diagnostic.F90.

References scale_const::const_grav, scale_gridtrans::gtrans_gsqrt, scale_gridtrans::gtrans_j13g, scale_gridtrans::gtrans_j23g, scale_gridtrans::i_xyw, scale_gridtrans::i_xyz, scale_grid_index::ia, scale_grid_index::ja, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_tracer::qa, scale_grid_real::real_cz, scale_grid_real::real_fz, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by mod_atmos_vars::atmos_vars_diagnostics().

64  use scale_comm, only: &
65  comm_vars8, &
66  comm_wait
67  use scale_const, only: &
68  grav => const_grav
69  use scale_atmos_thermodyn, only: &
70  thermodyn_qd => atmos_thermodyn_qd, &
71  thermodyn_r => atmos_thermodyn_r, &
72  thermodyn_temp_pres => atmos_thermodyn_temp_pres
73  use scale_grid_real, only: &
74  cz => real_cz, &
75  fz => real_fz
76  use scale_gridtrans, only: &
77  gsqrt => gtrans_gsqrt, &
78  j13g => gtrans_j13g, &
79  j23g => gtrans_j23g, &
80  i_xyw, &
81  i_xyz
82  implicit none
83 
84  real(RP), intent(out) :: POTT(KA,IA,JA)
85  real(RP), intent(out) :: TEMP(KA,IA,JA)
86  real(RP), intent(out) :: PRES(KA,IA,JA)
87  real(RP), intent(out) :: PHYD(KA,IA,JA)
88  real(RP), intent(out) :: W (KA,IA,JA)
89  real(RP), intent(out) :: U (KA,IA,JA)
90  real(RP), intent(out) :: V (KA,IA,JA)
91  real(RP), intent(out) :: N2 (KA,IA,JA)
92  real(RP), intent(in) :: DENS(KA,IA,JA)
93  real(RP), intent(in) :: MOMZ(KA,IA,JA)
94  real(RP), intent(in) :: MOMX(KA,IA,JA)
95  real(RP), intent(in) :: MOMY(KA,IA,JA)
96  real(RP), intent(in) :: RHOT(KA,IA,JA)
97  real(RP), intent(in) :: QTRC(KA,IA,JA,QA)
98 
99  real(RP) :: ph(KA)
100  real(RP) :: RPT(KA)
101  real(RP) :: q(QA)
102  real(RP) :: qdry, Rtot
103  real(RP) :: momws
104 
105  integer :: k, i, j
106  integer :: iq
107  !---------------------------------------------------------------------------
108 
109  call thermodyn_temp_pres( temp(:,:,:), & ! [OUT]
110  pres(:,:,:), & ! [OUT]
111  dens(:,:,:), & ! [IN]
112  rhot(:,:,:), & ! [IN]
113  qtrc(:,:,:,:), & ! [IN]
114  tracer_cv(:), & ! [IN]
115  tracer_r(:), & ! [IN]
116  tracer_mass(:) ) ! [IN]
117 
118  !$omp parallel do private(i,j,k,ph) OMP_SCHEDULE_ collapse(2)
119  do j = 1, ja
120  do i = 1, ia
121  ph(ke) = pres(ke,i,j) - dens(ke,i,j) * grav * ( fz(ke,i,j) - cz(ke,i,j) )
122  do k = ke, ks, -1
123  ph(k-1) = ph(k) + dens(k,i,j) * grav * ( fz(k,i,j) - fz(k-1,i,j) )
124  phyd(k,i,j) = ( ph(k) + ph(k-1) ) * 0.5_rp
125  end do
126  enddo
127  enddo
128 
129 
130 !OCL XFILL
131  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
132  do j = 1, ja
133  do i = 1, ia
134  do k = ks+1, ke-1
135  w(k,i,j) = 0.5_rp * ( momz(k-1,i,j)+momz(k,i,j) ) / dens(k,i,j)
136  enddo
137  enddo
138  enddo
139 !OCL XFILL
140  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
141  do j = 2, ja
142  do i = 2, ia
143 ! W(KS,i,j) = 0.5_RP * ( MOMZ(KS,i,j) ) / DENS(KS,i,j)
144 
145  ! at KS+1/2
146  momws = momz(ks,i,j) &
147  + ( 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) ) &
148  + 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) ) ) &
149  * 0.25_rp / gsqrt(ks,i,j,i_xyw)
150  ! at KS
151  ! momws at the surface is assumed to be zero
152  w(ks,i,j) = momws * 0.5_rp &
153  - ( j13g(ks,i,j,i_xyz) * ( momx(ks,i,j) + momx(ks,i-1,j) ) &
154  + j23g(ks,i,j,i_xyz) * ( momy(ks,i,j) + momy(ks,i,j-1) ) ) &
155  * 0.5_rp / ( dens(ks,i,j) * gsqrt(ks,i,j,i_xyz) )
156  enddo
157  enddo
158 !OCL XFILL
159  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
160  do j = 1, ja
161  do i = 1, ia
162  w(ke,i,j) = 0.5_rp * ( momz(ke-1,i,j) ) / dens(ke,i,j)
163  enddo
164  enddo
165 
166 !OCL XFILL
167  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
168  do j = 1, ja
169  do i = 2, ia
170  do k = ks, ke
171  u(k,i,j) = 0.5_rp * ( momx(k,i-1,j)+momx(k,i,j) ) / dens(k,i,j)
172  enddo
173  enddo
174  enddo
175 !OCL XFILL
176  !$omp parallel do private(j,k) OMP_SCHEDULE_ collapse(2)
177  do j = 1, ja
178  do k = ks, ke
179  u(k,1,j) = momx(k,1,j) / dens(k,1,j)
180  enddo
181  enddo
182 
183 !OCL XFILL
184  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
185  do j = 2, ja
186  do i = 1, ia
187  do k = ks, ke
188  v(k,i,j) = 0.5_rp * ( momy(k,i,j-1)+momy(k,i,j) ) / dens(k,i,j)
189  enddo
190  enddo
191  enddo
192 !OCL XFILL
193  !$omp parallel do private(i,k) OMP_SCHEDULE_ collapse(2)
194  do i = 1, ia
195  do k = ks, ke
196  v(k,i,1) = momy(k,i,1) / dens(k,i,1)
197  enddo
198  enddo
199 
200  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
201  do j = 1, ja
202  do i = 1, ia
203  w( 1:ks-1,i,j) = w(ks,i,j)
204  u( 1:ks-1,i,j) = u(ks,i,j)
205  v( 1:ks-1,i,j) = v(ks,i,j)
206  w(ke+1:ka, i,j) = w(ke,i,j)
207  u(ke+1:ka, i,j) = u(ke,i,j)
208  v(ke+1:ka, i,j) = v(ke,i,j)
209  enddo
210  enddo
211 
212  call comm_vars8( w(:,:,:), 1 )
213  call comm_vars8( u(:,:,:), 2 )
214  call comm_vars8( v(:,:,:), 3 )
215  call comm_wait ( w(:,:,:), 1, .false. )
216  call comm_wait ( u(:,:,:), 2, .false. )
217  call comm_wait ( v(:,:,:), 3, .false. )
218 
219 !OCL XFILL
220  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
221  do j = 1, ja
222  do i = 1, ia
223  do k = ks, ke
224  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
225  enddo
226  enddo
227  enddo
228 
229  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
230  !$omp private(i,j,q,rpt,qdry,rtot)
231  do j = 1, ja
232  do i = 1, ia
233  do k = ks, ke
234  do iq = 1, qa
235  q(iq) = qtrc(k,i,j,iq)
236  end do
237  call thermodyn_qd( qdry, q(:), tracer_mass(:) )
238  call thermodyn_r ( rtot, q(:), tracer_r(:), qdry )
239  rpt(k) = rtot * pott(k,i,j)
240  end do
241 
242  n2(ks,i,j) = grav * ( rpt(ks+1) - rpt(ks) ) / ( ( cz(ks+1,i,j) - cz(ks,i,j) ) * rpt(ks) )
243  do k = ks+1,ke-1
244  n2(k,i,j) = grav * ( rpt(k+1) - rpt(k-1) ) / ( ( cz(k+1,i,j) - cz(k-1,i,j) ) * rpt(k) )
245  end do
246  n2(ke,i,j) = grav * ( rpt(ke) - rpt(ke-1) ) / ( ( cz(ke,i,j) - cz(ke-1,i,j) ) * rpt(ke) )
247  end do
248  end do
249 
250  return
real(rp), dimension(:,:,:,:), allocatable, public gtrans_j23g
(2,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
module GRIDTRANS
module GRID (real space)
integer, public i_xyw
module COMMUNICATION
Definition: scale_comm.F90:23
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
real(rp), dimension(:,:,:,:), allocatable, public gtrans_j13g
(1,3) element of Jacobian matrix * {G}^1/2
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(:,:,:,:), allocatable, public gtrans_gsqrt
transformation metrics from Z to Xi, {G}^1/2
module ATMOSPHERE / Thermodynamics
integer, public i_xyz
Here is the caller graph for this function: