SCALE-RM
Functions/Subroutines
scale_atmos_diagnostic Module Reference

module atmosphere / diagnostic More...

Functions/Subroutines

subroutine, public atmos_diagnostic_get_therm_rhot (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, RHOT, Rtot, CVtot, CPtot, POTT, TEMP, PRES, EXNER)
 ATMOS_DIAGNOSTIC_get_therm_rhot potential temperature, temperature, pressure. More...
 
subroutine, public atmos_diagnostic_get_therm_rhoe (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, RHOE, Rtot, CPtot, CVtot, TEMP, POTT, PRES, EXNER)
 ATMOS_DIAGNOSTIC_get_therm_rhoe potential temperature, temperature, pressure. More...
 
subroutine, public atmos_diagnostic_get_phyd (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, CZ, FZ, PHYD, PHYDH)
 ATMOS_DIAGNOSTIC_get_phyd hydrostatic pressure. More...
 
subroutine, public atmos_diagnostic_get_n2 (KA, KS, KE, IA, IS, IE, JA, JS, JE, POTT, Rtot, CZ, FZ, F2H, N2)
 ATMOS_DIAGNOSTIC_get_n2 N^2. More...
 
subroutine, public atmos_diagnostic_get_potv (KA, KS, KE, IA, IS, IE, JA, JS, JE, POTT, Rtot, POTV)
 ATMOS_DIAGNOSTIC_get_potv virtual potential temperature. More...
 
subroutine, public atmos_diagnostic_get_teml (KA, KS, KE, IA, IS, IE, JA, JS, JE, TEMP, LHV, LHS, QC, QI, CPtot, TEML)
 ATMOS_DIAGNOSTIC_get_teml liqued water temperature. More...
 

Detailed Description

module atmosphere / diagnostic

Description
Calculate diagnostic variables
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_diagnostic_get_therm_rhot()

subroutine, public scale_atmos_diagnostic::atmos_diagnostic_get_therm_rhot ( 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(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  RHOT,
real(rp), dimension (ka,ia,ja), intent(in)  Rtot,
real(rp), dimension(ka,ia,ja), intent(in)  CVtot,
real(rp), dimension(ka,ia,ja), intent(in)  CPtot,
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)  EXNER 
)

ATMOS_DIAGNOSTIC_get_therm_rhot potential temperature, temperature, pressure.

Definition at line 59 of file scale_atmos_diagnostic.F90.

59  use scale_atmos_thermodyn, only: &
60  thermodyn_rhot2temp_pres => atmos_thermodyn_rhot2temp_pres
61  integer, intent(in) :: KA, KS, KE
62  integer, intent(in) :: IA, IS, IE
63  integer, intent(in) :: JA, JS, JE
64 
65  real(RP), intent(in) :: DENS (KA,IA,JA)
66  real(RP), intent(in) :: RHOT (KA,IA,JA)
67  real(RP), intent(in) :: Rtot (KA,IA,JA)
68  real(RP), intent(in) :: CVtot(KA,IA,JA)
69  real(RP), intent(in) :: CPtot(KA,IA,JA)
70 
71  real(RP), intent(out) :: POTT (KA,IA,JA)
72  real(RP), intent(out) :: TEMP (KA,IA,JA)
73  real(RP), intent(out) :: PRES (KA,IA,JA)
74  real(RP), intent(out) :: EXNER(KA,IA,JA)
75 
76  integer :: k, i, j
77  !---------------------------------------------------------------------------
78 
79  call thermodyn_rhot2temp_pres( ka, ks, ke, ia, is, ie, ja, js, je, &
80  dens(:,:,:), rhot(:,:,:), & ! (in)
81  rtot(:,:,:), cvtot(:,:,:), cptot(:,:,:), & ! (in)
82  temp(:,:,:), pres(:,:,:) ) ! (out)
83 
84 
85 !OCL XFILL
86  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
87  !$omp private(i,j,k) &
88  !$omp shared(POTT,EXNER,RHOT,DENS,TEMP) &
89  !$omp shared(KS,KE,IS,IE,JS,JE)
90  do j = js, je
91  do i = is, ie
92  do k = ks, ke
93  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
94  exner(k,i,j) = temp(k,i,j) / pott(k,i,j)
95  enddo
96  enddo
97  enddo
98 
99  return

Referenced by mod_atmos_vars::atmos_vars_calc_diagnostics().

Here is the caller graph for this function:

◆ atmos_diagnostic_get_therm_rhoe()

subroutine, public scale_atmos_diagnostic::atmos_diagnostic_get_therm_rhoe ( 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(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  RHOE,
real(rp), dimension (ka,ia,ja), intent(in)  Rtot,
real(rp), dimension(ka,ia,ja), intent(in)  CPtot,
real(rp), dimension(ka,ia,ja), intent(in)  CVtot,
real(rp), dimension (ka,ia,ja), intent(out)  TEMP,
real(rp), dimension (ka,ia,ja), intent(out)  POTT,
real(rp), dimension (ka,ia,ja), intent(out)  PRES,
real(rp), dimension(ka,ia,ja), intent(out)  EXNER 
)

ATMOS_DIAGNOSTIC_get_therm_rhoe potential temperature, temperature, pressure.

Definition at line 110 of file scale_atmos_diagnostic.F90.

110  use scale_const, only: &
111  pre00 => const_pre00
112  implicit none
113  integer, intent(in) :: KA, KS, KE
114  integer, intent(in) :: IA, IS, IE
115  integer, intent(in) :: JA, JS, JE
116 
117  real(RP), intent(in) :: DENS (KA,IA,JA)
118  real(RP), intent(in) :: RHOE (KA,IA,JA)
119  real(RP), intent(in) :: Rtot (KA,IA,JA)
120  real(RP), intent(in) :: CPtot(KA,IA,JA)
121  real(RP), intent(in) :: CVtot(KA,IA,JA)
122 
123  real(RP), intent(out) :: POTT (KA,IA,JA)
124  real(RP), intent(out) :: TEMP (KA,IA,JA)
125  real(RP), intent(out) :: PRES (KA,IA,JA)
126  real(RP), intent(out) :: EXNER(KA,IA,JA)
127 
128  integer :: i, j, k
129 
130  !$omp parallel do
131  do j = js, je
132  do i = is, ie
133  do k = ks, ke
134  temp(k,i,j) = rhoe(k,i,j) / ( cvtot(k,i,j) * dens(k,i,j) )
135  pres(k,i,j) = dens(k,i,j) * rtot(k,i,j) * temp(k,i,j)
136  exner(k,i,j) = ( pres(k,i,j) / pre00 )**( rtot(k,i,j) / cptot(k,i,j) )
137  pott(k,i,j) = temp(k,i,j) / exner(k,i,j)
138  end do
139  end do
140  end do
141 
142  return

References scale_const::const_pre00.

◆ atmos_diagnostic_get_phyd()

subroutine, public scale_atmos_diagnostic::atmos_diagnostic_get_phyd ( 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(in)  DENS,
real(rp), dimension ( ka,ia,ja), intent(in)  PRES,
real(rp), dimension ( ka,ia,ja), intent(in)  CZ,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
real(rp), dimension ( ka,ia,ja), intent(out)  PHYD,
real(rp), dimension(0:ka,ia,ja), intent(out)  PHYDH 
)

ATMOS_DIAGNOSTIC_get_phyd hydrostatic pressure.

Definition at line 157 of file scale_atmos_diagnostic.F90.

157  use scale_const, only: &
158  eps => const_eps, &
159  grav => const_grav
160  implicit none
161 
162  integer, intent(in) :: KA, KS, KE
163  integer, intent(in) :: IA, IS, IE
164  integer, intent(in) :: JA, JS, JE
165  real(RP), intent(in) :: DENS ( KA,IA,JA)
166  real(RP), intent(in) :: PRES ( KA,IA,JA)
167  real(RP), intent(in) :: CZ ( KA,IA,JA)
168  real(RP), intent(in) :: FZ (0:KA,IA,JA)
169  real(RP), intent(out) :: PHYD ( KA,IA,JA)
170  real(RP), intent(out) :: PHYDH(0:KA,IA,JA)
171 
172  real(RP) :: diff
173 
174  integer :: k, i, j
175  !---------------------------------------------------------------------------
176 
177  !$omp parallel do default(none) OMP_SCHEDULE_ &
178  !$omp private(i,j,k) &
179  !$omp private(diff) &
180  !$omp shared(PHYD,PHYDH,DENS,PRES,CZ,FZ,GRAV,EPS) &
181  !$omp shared(KS,KE,IS,IE,JS,JE)
182  do j = js, je
183  do i = is, ie
184  phydh(ke,i,j) = pres(ke,i,j) - dens(ke,i,j) * grav * ( fz(ke,i,j) - cz(ke,i,j) )
185  do k = ke-1, ks-1, -1
186  phydh(k,i,j) = phydh(k+1,i,j) + dens(k+1,i,j) * grav * ( fz(k+1,i,j) - fz(k,i,j) )
187  end do
188  phyd(ke,i,j) = pres(ke,i,j)
189  diff = 0.0_rp
190  do k = ke-1, ks, -1
191  phyd(k,i,j) = phyd(k+1,i,j) + ( dens(k+1,i,j) + dens(k,i,j) ) * grav * ( cz(k+1,i,j) - cz(k,i,j) ) * 0.5_rp
192  diff = diff + ( pres(k,i,j) - phyd(k,i,j) ) * ( fz(k,i,j) - fz(k-1,i,j) )
193  end do
194  diff = diff / ( fz(ke,i,j) - fz(ks-1,i,j) )
195  diff = max( diff, eps - phydh(ke,i,j) )
196  do k = ks-1, ke
197  phydh(k,i,j) = phydh(k,i,j) + diff
198  end do
199  do k = ks, ke
200  phyd(k,i,j) = phyd(k,i,j) + diff
201  end do
202  enddo
203  enddo
204 
205  return

References scale_const::const_eps, and scale_const::const_grav.

Referenced by mod_atmos_vars::atmos_vars_calc_diagnostics().

Here is the caller graph for this function:

◆ atmos_diagnostic_get_n2()

subroutine, public scale_atmos_diagnostic::atmos_diagnostic_get_n2 ( 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(in)  POTT,
real(rp), dimension(ka,ia,ja), intent(in)  Rtot,
real(rp), dimension ( ka,ia,ja), intent(in)  CZ,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
real(rp), dimension(ka,2,ia,ja), intent(in)  F2H,
real(rp), dimension (ka,ia,ja), intent(out)  N2 
)

ATMOS_DIAGNOSTIC_get_n2 N^2.

Definition at line 219 of file scale_atmos_diagnostic.F90.

219  use scale_const, only: &
220  grav => const_grav
221  integer, intent(in) :: KA, KS, KE
222  integer, intent(in) :: IA, IS, IE
223  integer, intent(in) :: JA, JS, JE
224 
225  real(RP), intent(in) :: POTT(KA,IA,JA)
226  real(RP), intent(in) :: Rtot(KA,IA,JA)
227 
228  real(RP), intent(in) :: CZ ( KA,IA,JA)
229  real(RP), intent(in) :: FZ (0:KA,IA,JA)
230  real(RP), intent(in) :: F2H(KA,2,IA,JA)
231 
232  real(RP), intent(out) :: N2 (KA,IA,JA)
233 
234  real(RP) :: RPT(KA)
235  real(RP) :: RPT_h(KA)
236 
237  integer :: k, i, j
238  !---------------------------------------------------------------------------
239 
240  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
241  !$omp private(i,j,k) &
242  !$omp private(RPT,RPT_h) &
243  !$omp shared(N2,POTT,Rtot,CZ,FZ,F2H,GRAV) &
244  !$omp shared(KS,KE,IS,IE,JS,JE)
245  do j = js, je
246  do i = is, ie
247  do k = ks, ke
248  rpt(k) = rtot(k,i,j) * pott(k,i,j)
249  end do
250 
251  do k = ks, ke-1
252  rpt_h(k) = f2h(k,1,i,j) * rpt(k+1) + f2h(k,2,i,j) * rpt(k)
253  end do
254 
255  n2(ks,i,j) = grav * ( rpt_h(ks) - rpt(ks) ) / ( ( fz(ks,i,j) - cz(ks,i,j) ) * rpt(ks) )
256  do k = ks+1,ke-1
257  n2(k,i,j) = grav * ( rpt_h(k) - rpt_h(k-1) ) / ( ( fz(k,i,j) - fz(k-1,i,j) ) * rpt(k) )
258  end do
259  n2(ke,i,j) = grav * ( rpt(ke) - rpt_h(ke-1) ) / ( ( cz(ke,i,j) - fz(ke-1,i,j) ) * rpt(ke) )
260  end do
261  end do
262 
263  return

References scale_const::const_grav.

Referenced by mod_atmos_vars::atmos_vars_get_diagnostic_3d().

Here is the caller graph for this function:

◆ atmos_diagnostic_get_potv()

subroutine, public scale_atmos_diagnostic::atmos_diagnostic_get_potv ( 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(in)  POTT,
real(rp), dimension(ka,ia,ja), intent(in)  Rtot,
real(rp), dimension(ka,ia,ja), intent(out)  POTV 
)

ATMOS_DIAGNOSTIC_get_potv virtual potential temperature.

Definition at line 277 of file scale_atmos_diagnostic.F90.

277  use scale_const, only: &
278  rdry => const_rdry
279  integer, intent(in) :: KA, KS, KE
280  integer, intent(in) :: IA, IS, IE
281  integer, intent(in) :: JA, JS, JE
282 
283  real(RP), intent(in) :: POTT(KA,IA,JA)
284  real(RP), intent(in) :: Rtot(KA,IA,JA)
285 
286  real(RP), intent(out) :: POTV(KA,IA,JA)
287 
288  integer :: k, i, j
289  !---------------------------------------------------------------------------
290 
291 !OCL XFILL
292  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
293  !$omp private(i,j,k) &
294  !$omp shared(POTV,POTT,Rtot,Rdry) &
295  !$omp shared(KS,KE,IS,IE,JS,JE)
296  do j = js, je
297  do i = is, ie
298  do k = ks, ke
299  potv(k,i,j) = pott(k,i,j) * rtot(k,i,j) / rdry
300  end do
301  end do
302  end do
303 
304  return

References scale_const::const_rdry.

Referenced by mod_atmos_vars::atmos_vars_get_diagnostic_3d().

Here is the caller graph for this function:

◆ atmos_diagnostic_get_teml()

subroutine, public scale_atmos_diagnostic::atmos_diagnostic_get_teml ( 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(in)  TEMP,
real(rp), dimension (ka,ia,ja), intent(in)  LHV,
real(rp), dimension (ka,ia,ja), intent(in)  LHS,
real(rp), dimension (ka,ia,ja), intent(in)  QC,
real(rp), dimension (ka,ia,ja), intent(in)  QI,
real(rp), dimension(ka,ia,ja), intent(in)  CPtot,
real(rp), dimension(ka,ia,ja), intent(out)  TEML 
)

ATMOS_DIAGNOSTIC_get_teml liqued water temperature.

Definition at line 320 of file scale_atmos_diagnostic.F90.

320  integer, intent(in) :: KA, KS, KE
321  integer, intent(in) :: IA, IS, IE
322  integer, intent(in) :: JA, JS, JE
323 
324  real(RP), intent(in) :: TEMP (KA,IA,JA)
325  real(RP), intent(in) :: LHV (KA,IA,JA)
326  real(RP), intent(in) :: LHS (KA,IA,JA)
327  real(RP), intent(in) :: QC (KA,IA,JA)
328  real(RP), intent(in) :: QI (KA,IA,JA)
329  real(RP), intent(in) :: CPtot(KA,IA,JA)
330 
331  real(RP), intent(out) :: TEML(KA,IA,JA)
332 
333  integer :: k, i, j
334  !---------------------------------------------------------------------------
335 
336 !OCL XFILL
337  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
338  !$omp private(i,j,k) &
339  !$omp shared(TEML,TEMP,LHV,LHS,QC,QI,CPtot) &
340  !$omp shared(KS,KE,IS,IE,JS,JE)
341  do j = js, je
342  do i = is, ie
343  do k = ks, ke
344  teml(k,i,j) = temp(k,i,j) &
345  - ( lhv(k,i,j) * qc(k,i,j) + lhs(k,i,j) * qi(k,i,j) ) / cptot(k,i,j)
346  end do
347  end do
348  end do
349 
350  return

Referenced by mod_atmos_vars::atmos_vars_get_diagnostic_3d().

Here is the caller graph for this function:
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_thermodyn
module atmosphere / thermodyn
Definition: scale_atmos_thermodyn.F90:11
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88