27 public :: atmos_adiabat_cape
28 public :: atmos_adiabat_liftparcel
30 interface atmos_adiabat_cape
33 end interface atmos_adiabat_cape
35 interface atmos_adiabat_liftparcel
38 end interface atmos_adiabat_liftparcel
78 QV, QC, Qdry, Rtot, CPtot, &
80 CAPE, CIN, LCL, LFC, LNB, &
81 DENS_p, TEMP_p, BUOY_p, QV_p, &
86 integer,
intent(in) :: KA, KS, KE
88 integer,
intent(in) :: Kstr
89 real(RP),
intent(in) :: DENS (ka)
90 real(RP),
intent(in) :: TEMP (ka)
91 real(RP),
intent(in) :: PRES (ka)
92 real(RP),
intent(in) :: QV (ka)
93 real(RP),
intent(in) :: QC (ka)
94 real(RP),
intent(in) :: Qdry (ka)
95 real(RP),
intent(in) :: Rtot (ka)
96 real(RP),
intent(in) :: CPtot(ka)
97 real(RP),
intent(in) :: CZ ( ka)
98 real(RP),
intent(in) :: FZ (0:ka)
100 real(RP),
intent(out) :: CAPE
101 real(RP),
intent(out) :: CIN
102 real(RP),
intent(out) :: LCL
103 real(RP),
intent(out) :: LFC
104 real(RP),
intent(out) :: LNB
105 real(RP),
intent(out) :: DENS_p(ka)
106 real(RP),
intent(out) :: TEMP_p(ka)
107 real(RP),
intent(out) :: BUOY_p(ka)
108 real(RP),
intent(out) :: QV_p (ka)
109 logical,
intent(out) :: converged
111 real(RP) :: BUOY_pf(ka)
112 integer :: kLCL, kLFC, kLNB
120 dens(:), temp(:), pres(:), &
122 qdry(:), rtot(:), cptot(:), &
123 dens_p(:), temp_p(:), qv_p(:), &
126 if ( .not. converged )
return 130 buoy_p(k) = grav * ( dens(k) - dens_p(k) ) / dens_p(k)
134 buoy_pf(k) = 0.5_rp * ( buoy_p(k+1) + buoy_p(k) )
136 buoy_p(ke) = buoy_p(ke-1)
148 if ( buoy_p(k) >= 0.0_rp )
then 155 if ( buoy_p(k) >= 0.0_rp )
then 161 if ( klcl >= kstr )
then 165 if ( klfc >= kstr )
then 169 if ( klnb >= kstr )
then 173 if ( klfc >= kstr .AND. klnb > kstr )
then 175 cape = cape + buoy_pf(k-1) * ( fz(k)-fz(k-1) )
179 if ( klfc >= kstr )
then 181 cin = cin + buoy_pf(k-1) * ( fz(k)-fz(k-1) )
189 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
192 QV, QC, Qdry, Rtot, CPtot, &
194 CAPE, CIN, LCL, LFC, LNB, &
195 DENS_p, TEMP_p, BUOY_p, QV_p )
201 integer,
intent(in) :: KA, KS, KE
202 integer,
intent(in) :: IA, IS, IE
203 integer,
intent(in) :: JA, JS, JE
205 integer,
intent(in) :: Kstr
206 real(RP),
intent(in) :: DENS (ka,ia,ja)
207 real(RP),
intent(in) :: TEMP (ka,ia,ja)
208 real(RP),
intent(in) :: PRES (ka,ia,ja)
209 real(RP),
intent(in) :: QV (ka,ia,ja)
210 real(RP),
intent(in) :: QC (ka,ia,ja)
211 real(RP),
intent(in) :: Qdry (ka,ia,ja)
212 real(RP),
intent(in) :: Rtot (ka,ia,ja)
213 real(RP),
intent(in) :: CPtot(ka,ia,ja)
214 real(RP),
intent(in) :: CZ ( ka,ia,ja)
215 real(RP),
intent(in) :: FZ (0:ka,ia,ja)
217 real(RP),
intent(out) :: CAPE(ia,ja)
218 real(RP),
intent(out) :: CIN (ia,ja)
219 real(RP),
intent(out) :: LCL (ia,ja)
220 real(RP),
intent(out) :: LFC (ia,ja)
221 real(RP),
intent(out) :: LNB (ia,ja)
223 real(RP),
intent(out),
optional,
target :: DENS_p(ka,ia,ja)
224 real(RP),
intent(out),
optional,
target :: TEMP_p(ka,ia,ja)
225 real(RP),
intent(out),
optional,
target :: BUOY_p(ka,ia,ja)
226 real(RP),
intent(out),
optional,
target :: QV_p (ka,ia,ja)
228 logical :: converged, error
230 real(RP),
pointer :: P_DENS(:,:,:), P_TEMP(:,:,:), P_BUOY(:,:,:), P_QV(:,:,:)
237 if (
present(dens_p) )
then 240 allocate( p_dens(ka,ia,ja) )
242 if (
present(temp_p) )
then 245 allocate( p_temp(ka,ia,ja) )
247 if (
present(buoy_p) )
then 250 allocate( p_buoy(ka,ia,ja) )
252 if (
present(qv_p) )
then 255 allocate( p_qv(ka,ia,ja) )
266 dens(:,i,j), temp(:,i,j), pres(:,i,j), &
267 qv(:,i,j), qc(:,i,j), qdry(:,i,j), &
268 rtot(:,i,j), cptot(:,i,j), &
269 cz(:,i,j), fz(:,i,j), &
270 cape(i,j), cin(i,j), &
271 lcl(i,j), lfc(i,j), lnb(i,j), &
272 p_dens(:,i,j), p_temp(:,i,j), &
273 p_buoy(:,i,j), p_qv(:,i,j), &
276 if ( .not. converged )
then 277 log_error(
"ATMOS_ADIABAT_cape_3D",*)
'[liftparcel] not converged! ', i, j
286 call file_history_in( p_dens(:,:,:),
'DENS_parcel',
'density profile in lifting parcel',
'kg/m3' )
287 call file_history_in( p_temp(:,:,:),
'TEMP_parcel',
'temperature profile in lifting parcel',
'K' )
288 call file_history_in( p_buoy(:,:,:),
'BUOY_parcel',
'buoyancy profile in lifting parcel',
'm/s2' )
289 call file_history_in( p_qv(:,:,:),
'QV_parcel',
'humidity profile in lifting parcel',
'kg/kg' )
291 if ( .not.
present(dens_p) )
deallocate( p_dens )
292 if ( .not.
present(temp_p) )
deallocate( p_temp )
293 if ( .not.
present(buoy_p) )
deallocate( p_buoy )
294 if ( .not.
present(qv_p ) )
deallocate( p_qv )
307 DENS, TEMP, PRES, QV, QC, &
309 DENS_p1D, TEMP_p1D, QV_p1D, &
312 hydrometeor_entr => atmos_hydrometeor_entr, &
315 atmos_saturation_moist_conversion_pres_liq
317 integer,
intent(in) :: KA, KS, KE
319 integer,
intent(in) :: Kstr
320 real(RP),
intent(in) :: DENS (ka)
321 real(RP),
intent(in) :: TEMP (ka)
322 real(RP),
intent(in) :: PRES (ka)
323 real(RP),
intent(in) :: QV (ka)
324 real(RP),
intent(in) :: QC (ka)
325 real(RP),
intent(in) :: QDRY (ka)
326 real(RP),
intent(in) :: Rtot (ka)
327 real(RP),
intent(in) :: CPtot(ka)
329 real(RP),
intent(out) :: DENS_p1D(ka)
330 real(RP),
intent(out) :: TEMP_p1D(ka)
331 real(RP),
intent(out) :: QV_p1D (ka)
332 integer,
intent(out) :: kLCL
333 logical,
intent(out) :: converged
336 real(RP) :: QV_p, QC_p, Qdry_p
337 real(RP) :: Rtot_p, CPtot_p
346 dens_p1d(k) = dens(k)
347 temp_p1d(k) = temp(k)
356 cptot_p = cptot(kstr)
358 call hydrometeor_entr( temp(kstr), pres(kstr), &
359 qv_p, qc_p, qdry_p, &
366 call atmos_saturation_moist_conversion_pres_liq( &
367 pres(k), entr_p, qdry_p, &
368 qv_p, qc_p, rtot_p, cptot_p, &
372 if ( .NOT. converged )
then 376 if ( qc_p > 0.0_rp .and. klcl == -1 ) klcl = k
379 fact = 1.0_rp / ( 1.0_rp - qc_p )
380 cptot_p = ( cptot_p -
cp_water * qc_p ) * fact
381 rtot_p = rtot_p * fact
382 qdry_p = qdry_p * fact
387 dens_p1d(k) = pres(k) / ( rtot_p * temp_p )
397 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
399 DENS, TEMP, PRES, QV, QC, &
401 DENS_p3D, TEMP_p3D, QV_p3D, &
406 hydrometeor_entr => atmos_hydrometeor_entr
408 atmos_saturation_moist_conversion_pres_liq
411 integer,
intent(in) :: KA, KS, KE
412 integer,
intent(in) :: IA, IS, IE
413 integer,
intent(in) :: JA, JS, JE
415 integer,
intent(in) :: Kstr
416 real(RP),
intent(in) :: DENS (ka,ia,ja)
417 real(RP),
intent(in) :: TEMP (ka,ia,ja)
418 real(RP),
intent(in) :: PRES (ka,ia,ja)
419 real(RP),
intent(in) :: QV (ka,ia,ja)
420 real(RP),
intent(in) :: QC (ka,ia,ja)
421 real(RP),
intent(in) :: QDRY (ka,ia,ja)
422 real(RP),
intent(in) :: Rtot (ka,ia,ja)
423 real(RP),
intent(in) :: CPtot(ka,ia,ja)
425 real(RP),
intent(out) :: DENS_p3D(ka,ia,ja)
426 real(RP),
intent(out) :: TEMP_p3D(ka,ia,ja)
427 real(RP),
intent(out) :: QV_p3D (ka,ia,ja)
428 integer,
intent(out) :: kLCL(ia,ja)
429 logical,
intent(out) :: converged
444 dens(:,i,j), temp(:,i,j), pres(:,i,j), &
445 qv(:,i,j), qc(:,i,j), &
446 qdry(:,i,j), rtot(:,i,j), cptot(:,i,j), &
447 dens_p3d(:,i,j), temp_p3d(:,i,j), qv_p3d(:,i,j), &
448 klcl(i,j), converged )
449 if ( .not. converged )
then 450 log_error(
"ATMOS_ADIABAT_liftparcel_3D",*)
'[liftparcel] not converged! ', i, j
module atmosphere / adiabat
module atmosphere / saturation
subroutine atmos_adiabat_cape_1d(KA, KS, KE, Kstr, DENS, TEMP, PRES, QV, QC, Qdry, Rtot, CPtot, CZ, FZ, CAPE, CIN, LCL, LFC, LNB, DENS_p, TEMP_p, BUOY_p, QV_p, converged)
Calc CAPE and CIN Type of parcel method: Pseudo-adiabatic ascend from lowermost layer of the model Re...
subroutine, public atmos_saturation_setup
Setup.
subroutine atmos_adiabat_cape_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, Kstr, DENS, TEMP, PRES, QV, QC, Qdry, Rtot, CPtot, CZ, FZ, CAPE, CIN, LCL, LFC, LNB, DENS_p, TEMP_p, BUOY_p, QV_p)
subroutine, public atmos_adiabat_setup
Setup.
module atmosphere / hydrometeor
real(rp), public const_grav
standard acceleration of gravity [m/s2]
subroutine, public prc_abort
Abort Process.
subroutine, public atmos_hydrometeor_setup
Setup.
subroutine atmos_adiabat_liftparcel_1d(KA, KS, KE, Kstr, DENS, TEMP, PRES, QV, QC, QDRY, Rtot, CPtot, DENS_p1D, TEMP_p1D, QV_p1D, kLCL, converged)
Calc temperature profile with lifting parcel Method: Pseudo-adiabatic ascend from lowermost layer of ...
subroutine atmos_adiabat_liftparcel_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, Kstr, DENS, TEMP, PRES, QV, QC, QDRY, Rtot, CPtot, DENS_p3D, TEMP_p3D, QV_p3D, kLCL, converged)
real(rp), public cp_water
CP for water [J/kg/K].