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  !$acc data copyin(DENS, RHOT, Rtot, CVtot, CPtot) copyout(POTT, TEMP, PRES, EXNER)
80 
81  call thermodyn_rhot2temp_pres( ka, ks, ke, ia, is, ie, ja, js, je, &
82  dens(:,:,:), rhot(:,:,:), & ! (in)
83  rtot(:,:,:), cvtot(:,:,:), cptot(:,:,:), & ! (in)
84  temp(:,:,:), pres(:,:,:) ) ! (out)
85 
86 
87 !OCL XFILL
88  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
89  !$omp private(i,j,k) &
90  !$omp shared(POTT,EXNER,RHOT,DENS,TEMP) &
91  !$omp shared(KS,KE,IS,IE,JS,JE)
92  !$acc kernels
93  do j = js, je
94  do i = is, ie
95  do k = ks, ke
96  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
97  exner(k,i,j) = temp(k,i,j) / pott(k,i,j)
98  enddo
99  enddo
100  enddo
101  !$acc end kernels
102 
103  !$acc end data
104 
105  return
106  end subroutine atmos_diagnostic_get_therm_rhot
107 
108  !-----------------------------------------------------------------------------
112  subroutine atmos_diagnostic_get_therm_rhoe( &
113  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
114  DENS, RHOE, Rtot, CPtot, CVtot, & ! (in)
115  temp, pott, pres, exner ) ! (out)
116  use scale_const, only: &
117  pre00 => const_pre00
118  implicit none
119  integer, intent(in) :: ka, ks, ke
120  integer, intent(in) :: ia, is, ie
121  integer, intent(in) :: ja, js, je
122 
123  real(rp), intent(in) :: dens (ka,ia,ja)
124  real(rp), intent(in) :: rhoe (ka,ia,ja)
125  real(rp), intent(in) :: rtot (ka,ia,ja)
126  real(rp), intent(in) :: cptot(ka,ia,ja)
127  real(rp), intent(in) :: cvtot(ka,ia,ja)
128 
129  real(rp), intent(out) :: pott (ka,ia,ja)
130  real(rp), intent(out) :: temp (ka,ia,ja)
131  real(rp), intent(out) :: pres (ka,ia,ja)
132  real(rp), intent(out) :: exner(ka,ia,ja)
133 
134  integer :: i, j, k
135 
136  !$omp parallel do
137  !$acc kernels
138  do j = js, je
139  do i = is, ie
140  do k = ks, ke
141  temp(k,i,j) = rhoe(k,i,j) / ( cvtot(k,i,j) * dens(k,i,j) )
142  pres(k,i,j) = dens(k,i,j) * rtot(k,i,j) * temp(k,i,j)
143  exner(k,i,j) = ( pres(k,i,j) / pre00 )**( rtot(k,i,j) / cptot(k,i,j) )
144  pott(k,i,j) = temp(k,i,j) / exner(k,i,j)
145  end do
146  end do
147  end do
148  !$acc end kernels
149 
150  return
151  end subroutine atmos_diagnostic_get_therm_rhoe
152 
153  !-----------------------------------------------------------------------------
157  subroutine atmos_diagnostic_get_phyd( &
158  KA, KS, KE, &
159  IA, IS, IE, &
160  JA, JS, JE, &
161  DENS, PRES, &
162  CZ, FZ, &
163  PHYD, &
164  PHYDH )
165  use scale_const, only: &
166  eps => const_eps, &
167  grav => const_grav
168  implicit none
169 
170  integer, intent(in) :: ka, ks, ke
171  integer, intent(in) :: ia, is, ie
172  integer, intent(in) :: ja, js, je
173  real(rp), intent(in) :: dens ( ka,ia,ja)
174  real(rp), intent(in) :: pres ( ka,ia,ja)
175  real(rp), intent(in) :: cz ( ka,ia,ja)
176  real(rp), intent(in) :: fz (0:ka,ia,ja)
177  real(rp), intent(out) :: phyd ( ka,ia,ja)
178  real(rp), intent(out) :: phydh(0:ka,ia,ja)
179 
180  real(rp) :: diff
181 
182  integer :: k, i, j
183  !---------------------------------------------------------------------------
184 
185  !$omp parallel do default(none) OMP_SCHEDULE_ &
186  !$omp private(i,j,k) &
187  !$omp private(diff) &
188  !$omp shared(PHYD,PHYDH,DENS,PRES,CZ,FZ,GRAV,EPS) &
189  !$omp shared(KS,KE,IS,IE,JS,JE)
190  !$acc kernels copyin(DENS, CZ, FZ) copyout(PHYD, PHYDH)
191  do j = js, je
192  do i = is, ie
193  phydh(ke,i,j) = pres(ke,i,j) - dens(ke,i,j) * grav * ( fz(ke,i,j) - cz(ke,i,j) )
194  do k = ke-1, ks-1, -1
195  phydh(k,i,j) = phydh(k+1,i,j) + dens(k+1,i,j) * grav * ( fz(k+1,i,j) - fz(k,i,j) )
196  end do
197  phyd(ke,i,j) = pres(ke,i,j)
198  do k = ke-1, ks, -1
199  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
200  end do
201  diff = 0.0_rp
202  do k = ks, ke-1
203  diff = diff + ( pres(k,i,j) - phyd(k,i,j) ) * ( fz(k,i,j) - fz(k-1,i,j) )
204  end do
205  diff = diff / ( fz(ke,i,j) - fz(ks-1,i,j) )
206  diff = max( diff, eps - phydh(ke,i,j) )
207  do k = ks-1, ke
208  phydh(k,i,j) = phydh(k,i,j) + diff
209  end do
210  do k = ks, ke
211  phyd(k,i,j) = phyd(k,i,j) + diff
212  end do
213  enddo
214  enddo
215  !$acc end kernels
216 
217  return
218  end subroutine atmos_diagnostic_get_phyd
219 
220  !-----------------------------------------------------------------------------
224  subroutine atmos_diagnostic_get_n2( &
225  KA, KS, KE, &
226  IA, IS, IE, &
227  JA, JS, JE, &
228  POTT, Rtot, &
229  CZ, FZ, F2H, &
230  N2 )
231  use scale_const, only: &
232  grav => const_grav
233  integer, intent(in) :: ka, ks, ke
234  integer, intent(in) :: ia, is, ie
235  integer, intent(in) :: ja, js, je
236 
237  real(rp), intent(in) :: pott(ka,ia,ja)
238  real(rp), intent(in) :: rtot(ka,ia,ja)
239 
240  real(rp), intent(in) :: cz ( ka,ia,ja)
241  real(rp), intent(in) :: fz (0:ka,ia,ja)
242  real(rp), intent(in) :: f2h(ka,2,ia,ja)
243 
244  real(rp), intent(out) :: n2 (ka,ia,ja)
245 
246  real(rp) :: rpt(ka)
247  real(rp) :: rpt_h(ka)
248 
249  integer :: k, i, j
250  !---------------------------------------------------------------------------
251 
252  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
253  !$omp private(i,j,k) &
254  !$omp private(RPT,RPT_h) &
255  !$omp shared(N2,POTT,Rtot,CZ,FZ,F2H,GRAV) &
256  !$omp shared(KS,KE,IS,IE,JS,JE)
257  !$acc kernels
258  !$acc loop private(RPT,RPT_h)
259  do j = js, je
260  !$acc loop private(RPT,RPT_h)
261  do i = is, ie
262  do k = ks, ke
263  rpt(k) = rtot(k,i,j) * pott(k,i,j)
264  end do
265 
266  do k = ks, ke-1
267  rpt_h(k) = f2h(k,1,i,j) * rpt(k+1) + f2h(k,2,i,j) * rpt(k)
268  end do
269 
270  n2(ks,i,j) = grav * ( rpt_h(ks) - rpt(ks) ) / ( ( fz(ks,i,j) - cz(ks,i,j) ) * rpt(ks) )
271  do k = ks+1,ke-1
272  n2(k,i,j) = grav * ( rpt_h(k) - rpt_h(k-1) ) / ( ( fz(k,i,j) - fz(k-1,i,j) ) * rpt(k) )
273  end do
274  n2(ke,i,j) = grav * ( rpt(ke) - rpt_h(ke-1) ) / ( ( cz(ke,i,j) - fz(ke-1,i,j) ) * rpt(ke) )
275  end do
276  end do
277  !$acc end kernels
278 
279  return
280  end subroutine atmos_diagnostic_get_n2
281 
282  !-----------------------------------------------------------------------------
286  subroutine atmos_diagnostic_get_potv( &
287  KA, KS, KE, &
288  IA, IS, IE, &
289  JA, JS, JE, &
290  POTT, &
291  Rtot, &
292  POTV )
293  use scale_const, only: &
294  rdry => const_rdry
295  integer, intent(in) :: ka, ks, ke
296  integer, intent(in) :: ia, is, ie
297  integer, intent(in) :: ja, js, je
298 
299  real(rp), intent(in) :: pott(ka,ia,ja)
300  real(rp), intent(in) :: rtot(ka,ia,ja)
301 
302  real(rp), intent(out) :: potv(ka,ia,ja)
303 
304  integer :: k, i, j
305  !---------------------------------------------------------------------------
306 
307 !OCL XFILL
308  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
309  !$omp private(i,j,k) &
310  !$omp shared(POTV,POTT,Rtot,Rdry) &
311  !$omp shared(KS,KE,IS,IE,JS,JE)
312  !$acc kernels
313  do j = js, je
314  do i = is, ie
315  do k = ks, ke
316  potv(k,i,j) = pott(k,i,j) * rtot(k,i,j) / rdry
317  end do
318  end do
319  end do
320  !$acc end kernels
321 
322  return
323  end subroutine atmos_diagnostic_get_potv
324 
325  !-----------------------------------------------------------------------------
329  subroutine atmos_diagnostic_get_teml( &
330  KA, KS, KE, &
331  IA, IS, IE, &
332  JA, JS, JE, &
333  TEMP, &
334  LHV, LHS, &
335  QC, QI, &
336  CPtot, &
337  TEML )
338  integer, intent(in) :: ka, ks, ke
339  integer, intent(in) :: ia, is, ie
340  integer, intent(in) :: ja, js, je
341 
342  real(rp), intent(in) :: temp (ka,ia,ja)
343  real(rp), intent(in) :: lhv (ka,ia,ja)
344  real(rp), intent(in) :: lhs (ka,ia,ja)
345  real(rp), intent(in) :: qc (ka,ia,ja)
346  real(rp), intent(in) :: qi (ka,ia,ja)
347  real(rp), intent(in) :: cptot(ka,ia,ja)
348 
349  real(rp), intent(out) :: teml(ka,ia,ja)
350 
351  integer :: k, i, j
352  !---------------------------------------------------------------------------
353 
354 !OCL XFILL
355  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
356  !$omp private(i,j,k) &
357  !$omp shared(TEML,TEMP,LHV,LHS,QC,QI,CPtot) &
358  !$omp shared(KS,KE,IS,IE,JS,JE)
359  !$acc kernels
360  do j = js, je
361  do i = is, ie
362  do k = ks, ke
363  teml(k,i,j) = temp(k,i,j) &
364  - ( lhv(k,i,j) * qc(k,i,j) + lhs(k,i,j) * qi(k,i,j) ) / cptot(k,i,j)
365  end do
366  end do
367  end do
368  !$acc end kernels
369 
370  return
371  end subroutine atmos_diagnostic_get_teml
372 
373 end module scale_atmos_diagnostic
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
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:338
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:293
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
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:165
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:231
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:59
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:116
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97