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, 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.

Referenced by mod_atmos_vars::atmos_vars_calc_diagnostics().

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
module atmosphere / thermodyn
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.

References scale_const::const_pre00.

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
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
module CONSTANT
Definition: scale_const.F90:11

◆ 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.

References scale_const::const_grav.

Referenced by mod_atmos_vars::atmos_vars_calc_diagnostics().

157  use scale_const, only: &
158  grav => const_grav
159  implicit none
160 
161  integer, intent(in) :: ka, ks, ke
162  integer, intent(in) :: ia, is, ie
163  integer, intent(in) :: ja, js, je
164  real(RP), intent(in) :: dens ( ka,ia,ja)
165  real(RP), intent(in) :: pres ( ka,ia,ja)
166  real(RP), intent(in) :: cz ( ka,ia,ja)
167  real(RP), intent(in) :: fz (0:ka,ia,ja)
168  real(RP), intent(out) :: phyd ( ka,ia,ja)
169  real(RP), intent(out) :: phydh(0:ka,ia,ja)
170 
171  integer :: k, i, j
172  !---------------------------------------------------------------------------
173 
174  !$omp parallel do default(none) OMP_SCHEDULE_ &
175  !$omp private(i,j,k) &
176  !$omp shared(PHYD,PHYDH,DENS,PRES,CZ,FZ,GRAV) &
177  !$omp shared(KS,KE,IS,IE,JS,JE)
178  do j = js, je
179  do i = is, ie
180  phydh(ke,i,j) = pres(ke,i,j) - dens(ke,i,j) * grav * ( fz(ke,i,j) - cz(ke,i,j) )
181  do k = ke, ks, -1
182  phydh(k-1,i,j) = phydh(k,i,j) + dens(k,i,j) * grav * ( fz(k,i,j) - fz(k-1,i,j) )
183 ! PHYD (k ,i,j) = 0.5_RP * ( PHYDH(k,i,j) + PHYDH(k-1,i,j) )
184  phyd(k ,i,j) = sqrt(phydh(k,i,j)) * sqrt(phydh(k-1,i,j))
185  end do
186  enddo
187  enddo
188 
189  return
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
module CONSTANT
Definition: scale_const.F90:11
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 (ka,ia,ja), intent(out)  N2 
)

ATMOS_DIAGNOSTIC_get_n2 N^2.

Definition at line 204 of file scale_atmos_diagnostic.F90.

References scale_const::const_grav.

Referenced by mod_atmos_vars::atmos_vars_get_diagnostic_3d().

204  use scale_const, only: &
205  grav => const_grav
206  integer, intent(in) :: ka, ks, ke
207  integer, intent(in) :: ia, is, ie
208  integer, intent(in) :: ja, js, je
209 
210  real(RP), intent(in) :: pott(ka,ia,ja)
211  real(RP), intent(in) :: rtot(ka,ia,ja)
212 
213  real(RP), intent(in) :: cz(ka,ia,ja)
214 
215  real(RP), intent(out) :: n2 (ka,ia,ja)
216 
217  real(RP) :: rpt(ka)
218 
219  integer :: k, i, j
220  !---------------------------------------------------------------------------
221 
222  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
223  !$omp private(i,j,k) &
224  !$omp private(RPT) &
225  !$omp shared(N2,POTT,Rtot,CZ,GRAV) &
226  !$omp shared(KS,KE,IS,IE,JS,JE)
227  do j = js, je
228  do i = is, ie
229  do k = ks, ke
230  rpt(k) = rtot(k,i,j) * pott(k,i,j)
231  end do
232 
233  n2(ks,i,j) = grav * ( rpt(ks+1) - rpt(ks) ) / ( ( cz(ks+1,i,j) - cz(ks,i,j) ) * rpt(ks) )
234  do k = ks+1,ke-1
235  n2(k,i,j) = grav * ( rpt(k+1) - rpt(k-1) ) / ( ( cz(k+1,i,j) - cz(k-1,i,j) ) * rpt(k ) )
236  end do
237  n2(ke,i,j) = grav * ( rpt(ke) - rpt(ke-1) ) / ( ( cz(ke,i,j) - cz(ke-1,i,j) ) * rpt(ke) )
238  end do
239  end do
240 
241  return
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
module CONSTANT
Definition: scale_const.F90:11
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 255 of file scale_atmos_diagnostic.F90.

References scale_const::const_rdry.

Referenced by mod_atmos_vars::atmos_vars_get_diagnostic_3d().

255  use scale_const, only: &
256  rdry => const_rdry
257  integer, intent(in) :: ka, ks, ke
258  integer, intent(in) :: ia, is, ie
259  integer, intent(in) :: ja, js, je
260 
261  real(RP), intent(in) :: pott(ka,ia,ja)
262  real(RP), intent(in) :: rtot(ka,ia,ja)
263 
264  real(RP), intent(out) :: potv(ka,ia,ja)
265 
266  integer :: k, i, j
267  !---------------------------------------------------------------------------
268 
269 !OCL XFILL
270  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
271  !$omp private(i,j,k) &
272  !$omp shared(POTV,POTT,Rtot,Rdry) &
273  !$omp shared(KS,KE,IS,IE,JS,JE)
274  do j = js, je
275  do i = is, ie
276  do k = ks, ke
277  potv(k,i,j) = pott(k,i,j) * rtot(k,i,j) / rdry
278  end do
279  end do
280  end do
281 
282  return
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
module CONSTANT
Definition: scale_const.F90:11
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 298 of file scale_atmos_diagnostic.F90.

Referenced by mod_atmos_vars::atmos_vars_get_diagnostic_3d().

298  integer, intent(in) :: ka, ks, ke
299  integer, intent(in) :: ia, is, ie
300  integer, intent(in) :: ja, js, je
301 
302  real(RP), intent(in) :: temp (ka,ia,ja)
303  real(RP), intent(in) :: lhv (ka,ia,ja)
304  real(RP), intent(in) :: lhs (ka,ia,ja)
305  real(RP), intent(in) :: qc (ka,ia,ja)
306  real(RP), intent(in) :: qi (ka,ia,ja)
307  real(RP), intent(in) :: cptot(ka,ia,ja)
308 
309  real(RP), intent(out) :: teml(ka,ia,ja)
310 
311  integer :: k, i, j
312  !---------------------------------------------------------------------------
313 
314 !OCL XFILL
315  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
316  !$omp private(i,j,k) &
317  !$omp shared(TEML,TEMP,LHV,LHS,QC,QI,CPtot) &
318  !$omp shared(KS,KE,IS,IE,JS,JE)
319  do j = js, je
320  do i = is, ie
321  do k = ks, ke
322  teml(k,i,j) = temp(k,i,j) &
323  - ( lhv(k,i,j) * qc(k,i,j) + lhs(k,i,j) * qi(k,i,j) ) / cptot(k,i,j)
324  end do
325  end do
326  end do
327 
328  return
Here is the caller graph for this function: