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  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
190  end subroutine atmos_diagnostic_get_phyd
191 
192  !-----------------------------------------------------------------------------
196  subroutine atmos_diagnostic_get_n2( &
197  KA, KS, KE, &
198  IA, IS, IE, &
199  JA, JS, JE, &
200  POTT, &
201  Rtot, &
202  CZ, &
203  N2 )
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
242  end subroutine atmos_diagnostic_get_n2
243 
244  !-----------------------------------------------------------------------------
248  subroutine atmos_diagnostic_get_potv( &
249  KA, KS, KE, &
250  IA, IS, IE, &
251  JA, JS, JE, &
252  POTT, &
253  Rtot, &
254  POTV )
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
283  end subroutine atmos_diagnostic_get_potv
284 
285  !-----------------------------------------------------------------------------
289  subroutine atmos_diagnostic_get_teml( &
290  KA, KS, KE, &
291  IA, IS, IE, &
292  JA, JS, JE, &
293  TEMP, &
294  LHV, LHS, &
295  QC, QI, &
296  CPtot, &
297  TEML )
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
329  end subroutine atmos_diagnostic_get_teml
330 
331 end module scale_atmos_diagnostic
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.
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.
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.
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
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.
module atmosphere / diagnostic
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
module CONSTANT
Definition: scale_const.F90:11
module profiler
Definition: scale_prof.F90:11
module atmosphere / thermodyn
module PRECISION
module STDIO
Definition: scale_io.F90:10
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.
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.