SCALE-RM
scale_atmos_diagnostic.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  !
32  public :: atmos_diagnostic_get_n2
33 
34  !-----------------------------------------------------------------------------
35  !
36  !++ Public parameters & variables
37  !
38  !-----------------------------------------------------------------------------
39  !
40  !++ Private procedure
41  !
42  !-----------------------------------------------------------------------------
43  !
44  !++ Private parameters & variables
45  !
46  !-----------------------------------------------------------------------------
47 contains
48  !-----------------------------------------------------------------------------
53  KA, KS, KE, &
54  IA, IS, IE, &
55  JA, JS, JE, &
56  DENS, RHOT, &
57  Rtot, CVtot, CPtot, &
58  POTT, TEMP, PRES, EXNER )
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
100  end subroutine atmos_diagnostic_get_therm_rhot
101 
102  !-----------------------------------------------------------------------------
106  subroutine atmos_diagnostic_get_therm_rhoe( &
107  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
108  DENS, RHOE, Rtot, CPtot, CVtot, & ! (in)
109  temp, pott, pres, exner ) ! (out)
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
143  end subroutine atmos_diagnostic_get_therm_rhoe
144 
145  !-----------------------------------------------------------------------------
149  subroutine atmos_diagnostic_get_phyd( &
150  KA, KS, KE, &
151  IA, IS, IE, &
152  JA, JS, JE, &
153  DENS, PRES, &
154  CZ, FZ, &
155  PHYD, &
156  PHYDH )
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
206  end subroutine atmos_diagnostic_get_phyd
207 
208  !-----------------------------------------------------------------------------
212  subroutine atmos_diagnostic_get_n2( &
213  KA, KS, KE, &
214  IA, IS, IE, &
215  JA, JS, JE, &
216  POTT, Rtot, &
217  CZ, FZ, F2H, &
218  N2 )
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
264  end subroutine atmos_diagnostic_get_n2
265 
266  !-----------------------------------------------------------------------------
270  subroutine atmos_diagnostic_get_potv( &
271  KA, KS, KE, &
272  IA, IS, IE, &
273  JA, JS, JE, &
274  POTT, &
275  Rtot, &
276  POTV )
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
305  end subroutine atmos_diagnostic_get_potv
306 
307  !-----------------------------------------------------------------------------
311  subroutine atmos_diagnostic_get_teml( &
312  KA, KS, KE, &
313  IA, IS, IE, &
314  JA, JS, JE, &
315  TEMP, &
316  LHV, LHS, &
317  QC, QI, &
318  CPtot, &
319  TEML )
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
351  end subroutine atmos_diagnostic_get_teml
352 
353 end module scale_atmos_diagnostic
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
scale_atmos_diagnostic::atmos_diagnostic_get_teml
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.
Definition: scale_atmos_diagnostic.F90:320
scale_atmos_diagnostic::atmos_diagnostic_get_potv
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.
Definition: scale_atmos_diagnostic.F90:277
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_atmos_diagnostic
module atmosphere / diagnostic
Definition: scale_atmos_diagnostic.F90:12
scale_atmos_diagnostic::atmos_diagnostic_get_phyd
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.
Definition: scale_atmos_diagnostic.F90:157
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_diagnostic::atmos_diagnostic_get_therm_rhot
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.
Definition: scale_atmos_diagnostic.F90:59
scale_atmos_diagnostic::atmos_diagnostic_get_n2
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.
Definition: scale_atmos_diagnostic.F90:219
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_atmos_diagnostic::atmos_diagnostic_get_therm_rhoe
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.
Definition: scale_atmos_diagnostic.F90:110
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88