SCALE-RM
Functions/Subroutines
scale_atmos_phy_rd_mstrnx Module Reference

module atmosphere / physics / radiation / mstrnX More...

Functions/Subroutines

subroutine, public atmos_phy_rd_mstrnx_setup (KA, KS, KE, CZ, FZ)
 Setup. More...
 
subroutine, public atmos_phy_rd_mstrnx_finalize
 finalize More...
 
subroutine, public atmos_phy_rd_mstrnx_flux (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, TEMP, PRES, QV, CZ, FZ, fact_ocean, fact_land, fact_urban, temp_sfc, albedo_sfc, solins, cosSZA, CLDFRAC, MP_Re, MP_Qe, AE_Re, AE_Qe, flux_rad, flux_rad_top, flux_rad_sfc_dn, dtau_s, dem_s)
 Radiation main. More...
 

Detailed Description

module atmosphere / physics / radiation / mstrnX

Description
Atmospheric radiation transfer process mstrnX Ref: Nakajima and Tanaka(1986) Nakajima et al.(2000) Sekiguchi and Nakajima(2008)
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_RD_MSTRN
    nametypedefault valuecomment
    ATMOS_PHY_RD_MSTRN_TOA real(RP)
    ATMOS_PHY_RD_MSTRN_KADD integer
    ATMOS_PHY_RD_MSTRN_GASPARA_IN_FILENAME character(len=H_LONG)
    ATMOS_PHY_RD_MSTRN_AEROPARA_IN_FILENAME character(len=H_LONG)
    ATMOS_PHY_RD_MSTRN_HYGROPARA_IN_FILENAME character(len=H_LONG)
    ATMOS_PHY_RD_MSTRN_NBAND integer
    ATMOS_PHY_RD_MSTRN_NPTYPE integer
    ATMOS_PHY_RD_MSTRN_NRADIUS integer
    ATMOS_PHY_RD_MSTRN_NRADIUS_CLOUD integer
    ATMOS_PHY_RD_MSTRN_NRADIUS_AERO integer
    ATMOS_PHY_RD_MSTRN_ONLY_QCI logical .false.
    ATMOS_PHY_RD_MSTRN_ONLY_TROPOCLOUD logical .false.
    ATMOS_PHY_RD_MSTRN_USE_AERO logical .false.

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_rd_mstrnx_setup()

subroutine, public scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_setup ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension(ka), intent(in)  CZ,
real(rp), dimension(0:ka), intent(in)  FZ 
)

Setup.

Definition at line 222 of file scale_atmos_phy_rd_mstrnx.F90.

222  use scale_prc, only: &
223  prc_abort
224  use scale_time, only: &
226  use scale_atmos_grid_cartesc_real, only: &
228  use scale_atmos_phy_rd_profile, only: &
229  rd_profile_setup => atmos_phy_rd_profile_setup, &
230  rd_profile_setup_zgrid => atmos_phy_rd_profile_setup_zgrid, &
231  rd_profile_read => atmos_phy_rd_profile_read
232  use scale_atmos_hydrometeor, only: &
233  n_hyd
234  use scale_atmos_aerosol, only: &
235  n_ae
236  implicit none
237  integer, intent(in) :: KA, KS, KE
238 
239  real(RP), intent(in) :: CZ(KA)
240  real(RP), intent(in) :: FZ(0:KA)
241 
242  real(RP) :: ATMOS_PHY_RD_MSTRN_TOA
243  integer :: ATMOS_PHY_RD_MSTRN_KADD
244  character(len=H_LONG) :: ATMOS_PHY_RD_MSTRN_GASPARA_IN_FILENAME
245  character(len=H_LONG) :: ATMOS_PHY_RD_MSTRN_AEROPARA_IN_FILENAME
246  character(len=H_LONG) :: ATMOS_PHY_RD_MSTRN_HYGROPARA_IN_FILENAME
247  integer :: ATMOS_PHY_RD_MSTRN_nband
248  integer :: ATMOS_PHY_RD_MSTRN_nptype
249  integer :: ATMOS_PHY_RD_MSTRN_nradius
250  integer :: ATMOS_PHY_RD_MSTRN_nradius_cloud
251  integer :: ATMOS_PHY_RD_MSTRN_nradius_aero
252 
253  namelist / param_atmos_phy_rd_mstrn / &
254  atmos_phy_rd_mstrn_toa, &
255  atmos_phy_rd_mstrn_kadd, &
256  atmos_phy_rd_mstrn_gaspara_in_filename, &
257  atmos_phy_rd_mstrn_aeropara_in_filename, &
258  atmos_phy_rd_mstrn_hygropara_in_filename, &
259  atmos_phy_rd_mstrn_nband, &
260  atmos_phy_rd_mstrn_nptype, &
261  atmos_phy_rd_mstrn_nradius, &
262  atmos_phy_rd_mstrn_nradius_cloud, &
263  atmos_phy_rd_mstrn_nradius_aero, &
264  atmos_phy_rd_mstrn_only_qci, &
265  atmos_phy_rd_mstrn_only_tropocloud, &
266  atmos_phy_rd_mstrn_use_aero
267 
268  integer :: KMAX
269  integer :: ngas, ncfc
270  integer :: ierr
271  !---------------------------------------------------------------------------
272 
273  log_newline
274  log_info("ATMOS_PHY_RD_mstrnx_setup",*) 'Setup'
275  log_info("ATMOS_PHY_RD_mstrnx_setup",*) 'Sekiguchi and Nakajima (2008) mstrnX radiation process'
276 
277  !--- read namelist
278  atmos_phy_rd_mstrn_toa = rd_toa
279  atmos_phy_rd_mstrn_kadd = rd_kadd
280  atmos_phy_rd_mstrn_gaspara_in_filename = mstrn_gaspara_inputfile
281  atmos_phy_rd_mstrn_aeropara_in_filename = mstrn_aeropara_inputfile
282  atmos_phy_rd_mstrn_hygropara_in_filename = mstrn_hygropara_inputfile
283  atmos_phy_rd_mstrn_nband = mstrn_nband
284  atmos_phy_rd_mstrn_nptype = mstrn_nptype
285  atmos_phy_rd_mstrn_nradius = mstrn_nradius
286  atmos_phy_rd_mstrn_nradius_cloud = -1
287  atmos_phy_rd_mstrn_nradius_aero = -1
288 
289  rewind(io_fid_conf)
290  read(io_fid_conf,nml=param_atmos_phy_rd_mstrn,iostat=ierr)
291  if( ierr < 0 ) then !--- missing
292  log_info("ATMOS_PHY_RD_mstrnx_setup",*) 'Not found namelist. Default used.'
293  elseif( ierr > 0 ) then !--- fatal error
294  log_error("ATMOS_PHY_RD_mstrnx_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_RD_MSTRN. Check!'
295  call prc_abort
296  endif
297  log_nml(param_atmos_phy_rd_mstrn)
298 
299  rd_toa = atmos_phy_rd_mstrn_toa
300  rd_kadd = atmos_phy_rd_mstrn_kadd
301  mstrn_gaspara_inputfile = atmos_phy_rd_mstrn_gaspara_in_filename
302  mstrn_aeropara_inputfile = atmos_phy_rd_mstrn_aeropara_in_filename
303  mstrn_hygropara_inputfile = atmos_phy_rd_mstrn_hygropara_in_filename
304  mstrn_nband = atmos_phy_rd_mstrn_nband
305  mstrn_nptype = atmos_phy_rd_mstrn_nptype
306 
307  if ( atmos_phy_rd_mstrn_nradius_cloud < 0 ) then ! set default
308  atmos_phy_rd_mstrn_nradius_cloud = atmos_phy_rd_mstrn_nradius
309  endif
310  if ( atmos_phy_rd_mstrn_nradius_aero < 0 ) then ! set default
311  atmos_phy_rd_mstrn_nradius_aero = atmos_phy_rd_mstrn_nradius
312  endif
313 
314  mstrn_nradius = max( atmos_phy_rd_mstrn_nradius, &
315  atmos_phy_rd_mstrn_nradius_cloud, &
316  atmos_phy_rd_mstrn_nradius_aero )
317 
318  !--- setup MP/AE - RD(particle) coupling parameter
319  allocate( ptype_nradius(mstrn_nptype) )
320 
321  if ( mstrn_nptype == 9 ) then
322 
323  i_mpae2rd( 1) = 1 ! cloud water Qc -> Water
324  i_mpae2rd( 2) = 1 ! rain Qr -> Water
325  i_mpae2rd( 3) = 2 ! cloud ice Qi -> Ice
326  i_mpae2rd( 4) = 2 ! snow Qs -> Ice
327  i_mpae2rd( 5) = 2 ! graupel Qg -> Ice
328  i_mpae2rd( 6) = 2 ! hail Qh -> Ice
329  i_mpae2rd( 7) = 3 ! aerosol type1 -> Soil dust
330  i_mpae2rd( 8) = 4 ! aerosol type2 -> Carbonacerous (BC/OC=0.3)
331  i_mpae2rd( 9) = 5 ! aerosol type3 -> Carbonacerous (BC/OC=0.15)
332  i_mpae2rd(10) = 6 ! aerosol type4 -> Carbonacerous (BC/OC=0.)
333  i_mpae2rd(11) = 7 ! aerosol type5 -> Black carbon
334  i_mpae2rd(12) = 8 ! aerosol type6 -> Sulfate
335  i_mpae2rd(13) = 9 ! aerosol type7 -> Sea salt
336 
337  ptype_nradius( 1) = atmos_phy_rd_mstrn_nradius_cloud
338  ptype_nradius( 2) = atmos_phy_rd_mstrn_nradius_cloud
339  ptype_nradius( 3) = atmos_phy_rd_mstrn_nradius_aero
340  ptype_nradius( 4) = atmos_phy_rd_mstrn_nradius_aero
341  ptype_nradius( 5) = atmos_phy_rd_mstrn_nradius_aero
342  ptype_nradius( 6) = atmos_phy_rd_mstrn_nradius_aero
343  ptype_nradius( 7) = atmos_phy_rd_mstrn_nradius_aero
344  ptype_nradius( 8) = atmos_phy_rd_mstrn_nradius_aero
345  ptype_nradius( 9) = atmos_phy_rd_mstrn_nradius_aero
346 
347  elseif( mstrn_nptype == 12 ) then
348 
349  i_mpae2rd( 1) = 1 ! cloud water Qc -> CLOUD
350  i_mpae2rd( 2) = 2 ! rain Qr -> RAIN
351  i_mpae2rd( 3) = 3 ! cloud ice Qi -> ICE
352  i_mpae2rd( 4) = 4 ! snow Qs -> SNOW
353  i_mpae2rd( 5) = 5 ! graupel Qg -> GRAUPEL
354  i_mpae2rd( 6) = 5 ! hail Qh -> GRAUPEL
355  i_mpae2rd( 7) = 6 ! aerosol type1 -> Soil dust
356  i_mpae2rd( 8) = 7 ! aerosol type2 -> Carbonacerous (BC/OC=0.3)
357  i_mpae2rd( 9) = 8 ! aerosol type3 -> Carbonacerous (BC/OC=0.15)
358  i_mpae2rd(10) = 9 ! aerosol type4 -> Carbonacerous (BC/OC=0.)
359  i_mpae2rd(11) = 10 ! aerosol type5 -> Black carbon
360  i_mpae2rd(12) = 11 ! aerosol type6 -> Sulfate
361  i_mpae2rd(13) = 12 ! aerosol type7 -> Sea salt
362 
363  ptype_nradius( 1) = atmos_phy_rd_mstrn_nradius_cloud
364  ptype_nradius( 2) = atmos_phy_rd_mstrn_nradius_cloud
365  ptype_nradius( 3) = atmos_phy_rd_mstrn_nradius_cloud
366  ptype_nradius( 4) = atmos_phy_rd_mstrn_nradius_cloud
367  ptype_nradius( 5) = atmos_phy_rd_mstrn_nradius_cloud
368  ptype_nradius( 6) = atmos_phy_rd_mstrn_nradius_aero
369  ptype_nradius( 7) = atmos_phy_rd_mstrn_nradius_aero
370  ptype_nradius( 8) = atmos_phy_rd_mstrn_nradius_aero
371  ptype_nradius( 9) = atmos_phy_rd_mstrn_nradius_aero
372  ptype_nradius(10) = atmos_phy_rd_mstrn_nradius_aero
373  ptype_nradius(11) = atmos_phy_rd_mstrn_nradius_aero
374  ptype_nradius(12) = atmos_phy_rd_mstrn_nradius_aero
375 
376  endif
377  !$acc enter data copyin(I_MPAE2RD, ptype_nradius)
378 
379  !--- setup MSTRN parameter
380  call rd_mstrn_setup( ngas, & ! [OUT]
381  ncfc ) ! [OUT]
382 
383  !--- setup climatological profile
384  call rd_profile_setup
385 
386  kmax = ke - ks + 1
387  rd_kmax = kmax + rd_kadd
388 
389  !--- allocate arrays
390  ! input
391  allocate( rd_zh(rd_kmax+1) )
392  allocate( rd_z(rd_kmax ) )
393 
394  allocate( rd_rhodz(rd_kmax ) )
395  allocate( rd_pres(rd_kmax ) )
396  allocate( rd_presh(rd_kmax+1) )
397  allocate( rd_temp(rd_kmax ) )
398  allocate( rd_temph(rd_kmax+1) )
399 
400  allocate( rd_gas(rd_kmax,ngas ) )
401  allocate( rd_cfc(rd_kmax,ncfc ) )
402  allocate( rd_aerosol_conc(rd_kmax,rd_naero) )
403  allocate( rd_aerosol_radi(rd_kmax,rd_naero) )
404  allocate( rd_cldfrac(rd_kmax ) )
405 
406  !--- setup vartical grid for radiation (larger TOA than Model domain)
407  call rd_profile_setup_zgrid( &
408  ka, ks, ke, &
409  rd_kmax, rd_kadd, & ! [IN]
410  rd_toa, cz(:), fz(:), & ! [IN]
411  rd_zh(:), rd_z(:) ) ! [OUT]
412 
413  !--- read climatological profile
414  call rd_profile_read( rd_kmax, & ! [IN]
415  ngas, & ! [IN]
416  ncfc, & ! [IN]
417  rd_naero, & ! [IN]
419  time_nowdate(:), & ! [IN]
420  rd_zh(:), & ! [IN]
421  rd_z(:), & ! [IN]
422  rd_rhodz(:), & ! [OUT]
423  rd_pres(:), & ! [OUT]
424  rd_presh(:), & ! [OUT]
425  rd_temp(:), & ! [OUT]
426  rd_temph(:), & ! [OUT]
427  rd_gas(:,:), & ! [OUT]
428  rd_cfc(:,:), & ! [OUT]
429  rd_aerosol_conc(:,:), & ! [OUT]
430  rd_aerosol_radi(:,:), & ! [OUT]
431  rd_cldfrac(:) ) ! [OUT]
432 
433  !$acc enter data copyin(RD_rhodz,RD_pres,RD_presh,RD_temp,RD_temph,RD_gas,RD_cfc,RD_aerosol_conc,RD_aerosol_radi,RD_cldfrac)
434 
435  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_basepoint_lat, scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup_zgrid(), scale_io::io_fid_conf, scale_atmos_aerosol::n_ae, scale_atmos_hydrometeor::n_hyd, scale_prc::prc_abort(), and scale_time::time_nowdate.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_rd_mstrnx_finalize()

subroutine, public scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_finalize

finalize

Definition at line 441 of file scale_atmos_phy_rd_mstrnx.F90.

441 
442  !$acc exit data delete(I_MPAE2RD, ptype_nradius)
443 
444  !$acc exit data &
445  !$acc delete(waveh, &
446  !$acc wgtch, fitPLK, logfitP, logfitT, fitT, &
447  !$acc radmode, ngasabs, igasabs, &
448  !$acc fsol, q, qmol, rayleigh, acfc_pow, nch, AKD, SKD )
449 
450  !$acc exit data delete(RD_rhodz,RD_pres,RD_presh,RD_temp,RD_temph,RD_gas,RD_cfc,RD_aerosol_conc,RD_aerosol_radi,RD_cldfrac)
451 
452  deallocate( ptype_nradius )
453  deallocate( rd_zh )
454  deallocate( rd_z )
455 
456  deallocate( rd_rhodz )
457  deallocate( rd_pres )
458  deallocate( rd_presh )
459  deallocate( rd_temp )
460  deallocate( rd_temph )
461 
462  deallocate( rd_gas )
463  deallocate( rd_cfc )
464  deallocate( rd_aerosol_conc )
465  deallocate( rd_aerosol_radi )
466  deallocate( rd_cldfrac )
467 
468  deallocate( waveh )
469  deallocate( logfitp )
470  deallocate( fitt )
471  deallocate( logfitt )
472 
473  deallocate( iflgb )
474  deallocate( nch )
475  deallocate( wgtch )
476  deallocate( ngasabs )
477  deallocate( igasabs )
478 
479  deallocate( akd )
480  deallocate( skd )
481  deallocate( acfc_pow )
482 
483  deallocate( fitplk )
484  deallocate( fsol )
485  deallocate( sfc )
486  deallocate( rayleigh )
487 
488  deallocate( qmol )
489  deallocate( q )
490 
491  deallocate( hygro_flag )
492  deallocate( radmode )
493 
494  return

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver_finalize().

Here is the caller graph for this function:

◆ atmos_phy_rd_mstrnx_flux()

subroutine, public scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux ( 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)  TEMP,
real(rp), dimension (ka,ia,ja), intent(in)  PRES,
real(rp), dimension (ka,ia,ja), intent(in)  QV,
real(rp), dimension ( ka,ia,ja), intent(in)  CZ,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
real(rp), dimension (ia,ja), intent(in)  fact_ocean,
real(rp), dimension (ia,ja), intent(in)  fact_land,
real(rp), dimension (ia,ja), intent(in)  fact_urban,
real(rp), dimension (ia,ja), intent(in)  temp_sfc,
real(rp), dimension (ia,ja,n_rad_dir,n_rad_rgn), intent(in)  albedo_sfc,
real(rp), dimension (ia,ja), intent(in)  solins,
real(rp), dimension (ia,ja), intent(in)  cosSZA,
real(rp), dimension (ka,ia,ja), intent(in)  CLDFRAC,
real(rp), dimension (ka,ia,ja,n_hyd), intent(in)  MP_Re,
real(rp), dimension (ka,ia,ja,n_hyd), intent(in)  MP_Qe,
real(rp), dimension (ka,ia,ja,n_ae), intent(in)  AE_Re,
real(rp), dimension (ka,ia,ja,n_ae), intent(in)  AE_Qe,
real(rp), dimension (ka,ia,ja,2,2,2), intent(out)  flux_rad,
real(rp), dimension (ia,ja,2,2,2), intent(out)  flux_rad_top,
real(rp), dimension(ia,ja,n_rad_dir,n_rad_rgn), intent(out)  flux_rad_sfc_dn,
real(rp), dimension(ka,ia,ja), intent(out), optional  dtau_s,
real(rp), dimension (ka,ia,ja), intent(out), optional  dem_s 
)

Radiation main.

Definition at line 514 of file scale_atmos_phy_rd_mstrnx.F90.

514  !Jval )
515 #ifdef _OPENACC
516  use openacc
517 #endif
518  use scale_const, only: &
519  eps => const_eps, &
520  mdry => const_mdry, &
521  mvap => const_mvap, &
522  ppm => const_ppm
523  use scale_time, only: &
525  use scale_atmos_grid_cartesc_real, only: &
527  use scale_atmos_hydrometeor, only: &
528  n_hyd, &
529  i_hc, &
530  i_hi, &
531  hyd_dens
532  use scale_atmos_aerosol, only: &
533  n_ae, &
534  ae_dens
535  use scale_atmos_phy_rd_profile, only: &
536  rd_profile_read => atmos_phy_rd_profile_read, &
537  rd_profile_use_climatology => atmos_phy_rd_profile_use_climatology
538  implicit none
539  integer, intent(in) :: KA, KS, KE
540  integer, intent(in) :: IA, IS, IE
541  integer, intent(in) :: JA, JS, JE
542 
543  real(RP), intent(in) :: DENS (KA,IA,JA)
544  real(RP), intent(in) :: TEMP (KA,IA,JA)
545  real(RP), intent(in) :: PRES (KA,IA,JA)
546  real(RP), intent(in) :: QV (KA,IA,JA)
547  real(RP), intent(in) :: CZ ( KA,IA,JA)
548  real(RP), intent(in) :: FZ (0:KA,IA,JA)
549  real(RP), intent(in) :: fact_ocean (IA,JA)
550  real(RP), intent(in) :: fact_land (IA,JA)
551  real(RP), intent(in) :: fact_urban (IA,JA)
552  real(RP), intent(in) :: temp_sfc (IA,JA)
553  real(RP), intent(in) :: albedo_sfc (IA,JA,N_RAD_DIR,N_RAD_RGN)
554  real(RP), intent(in) :: solins (IA,JA)
555  real(RP), intent(in) :: cosSZA (IA,JA)
556  real(RP), intent(in) :: cldfrac (KA,IA,JA)
557  real(RP), intent(in) :: MP_Re (KA,IA,JA,N_HYD)
558  real(RP), intent(in) :: MP_Qe (KA,IA,JA,N_HYD)
559  real(RP), intent(in) :: AE_Re (KA,IA,JA,N_AE)
560  real(RP), intent(in) :: AE_Qe (KA,IA,JA,N_AE)
561  real(RP), intent(out) :: flux_rad (KA,IA,JA,2,2,2)
562  real(RP), intent(out) :: flux_rad_top (IA,JA,2,2,2)
563  real(RP), intent(out) :: flux_rad_sfc_dn(IA,JA,N_RAD_DIR,N_RAD_RGN)
564 
565  real(RP), intent(out), optional :: dtau_s(KA,IA,JA) ! 0.67 micron cloud optical depth
566  real(RP), intent(out), optional :: dem_s (KA,IA,JA) ! 10.5 micron cloud emissivity
567 ! real(RP), intent(out), optional :: Jval (KA,IA,JA,CH_QA_photo)
568 
569  integer :: tropopause(IA,JA)
570  real(RP) :: gamma
571 
572  real(RP), parameter :: min_cldfrac = 1.e-8_rp
573 
574  real(RP) :: rhodz_merge (RD_KMAX,IA,JA)
575  real(RP) :: pres_merge (RD_KMAX,IA,JA)
576  real(RP) :: temp_merge (RD_KMAX,IA,JA)
577  real(RP) :: temph_merge (RD_KMAX+1,IA,JA)
578 
579  real(RP) :: gas_merge (RD_KMAX,IA,JA,MSTRN_ngas)
580  real(RP) :: cfc_merge (RD_KMAX,IA,JA,MSTRN_ncfc)
581  real(RP) :: aerosol_conc_merge(RD_KMAX,IA,JA,RD_naero )
582  real(RP) :: aerosol_radi_merge(RD_KMAX,IA,JA,RD_naero )
583  real(RP) :: cldfrac_merge (RD_KMAX,IA,JA)
584 
585  ! output
586  real(RP) :: flux_rad_merge(RD_KMAX+1,IA,JA,2,2,MSTRN_ncloud)
587  real(RP) :: tauCLD_067u (RD_KMAX,IA,JA) ! 0.67 micron cloud optical depth
588  real(RP) :: emisCLD_105u (RD_KMAX,IA,JA) ! 10.5 micron cloud emissivity
589 
590  real(RP) :: zerosw
591 
592  integer :: ihydro, iaero
593  integer :: RD_k, k, i, j, v, ic
594  !---------------------------------------------------------------------------
595 
596  log_progress(*) 'atmosphere / physics / radiation / mstrnX'
597 
598  call prof_rapstart('RD_Profile', 3)
599 
600  !$acc data copyin(DENS,TEMP,PRES,QV,CZ,FZ,fact_ocean,fact_land,fact_urban,temp_sfc,albedo_sfc,solins,cosSZA,cldfrac,MP_Re,AE_Re,AE_Qe) &
601  !$acc copyout(flux_rad,flux_rad_top,flux_rad_sfc_dn) &
602  !$acc create(tropopause,rhodz_merge,pres_merge,temp_merge,temph_merge,gas_merge,cfc_merge,aerosol_conc_merge,aerosol_radi_merge,cldfrac_merge,flux_rad_merge,tauCLD_067u,emisCLD_105u)
603 
604  !$acc data copyout(dtau_s) if(acc_is_present(dtau_s))
605  !$acc data copyout(dem_s) if(acc_is_present(dem_s))
606 
607  if ( atmos_phy_rd_mstrn_only_tropocloud ) then
608  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
609  !$omp private(k,i,j, &
610  !$omp gamma) &
611  !$omp shared (tropopause,pres,temp,CZ, &
612  !$omp KS,KE,IS,IE,JS,JE)
613  !$acc kernels
614  do j = js, je
615  do i = is, ie
616  tropopause(i,j) = ke+1
617  do k = ke, ks, -1
618  if ( pres(k,i,j) >= 300.e+2_rp ) then
619  exit
620  elseif( pres(k,i,j) < 50.e+2_rp ) then
621  tropopause(i,j) = k
622  else
623  gamma = ( temp(k+1,i,j) - temp(k,i,j) ) / ( cz(k+1,i,j) - cz(k,i,j) )
624  if( gamma > 1.e-4_rp ) tropopause(i,j) = k
625  endif
626  enddo
627  enddo
628  enddo
629  !$acc end kernels
630  else
631  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
632  !$omp private(i,j) &
633  !$omp shared (tropopause, &
634  !$omp KE,IS,IE,JS,JE)
635 !OCL XFILL
636  !$acc kernels
637  do j = js, je
638  do i = is, ie
639  tropopause(i,j) = ke+1
640  end do
641  end do
642  !$acc end kernels
643  endif
644 
645  ! marge basic profile and value in model domain
646 
647  if ( rd_profile_use_climatology ) then
648  call rd_profile_read( rd_kmax, & ! [IN]
649  mstrn_ngas, & ! [IN]
650  mstrn_ncfc, & ! [IN]
651  rd_naero, & ! [IN]
653  time_nowdate(:), & ! [IN]
654  rd_zh(:), & ! [IN]
655  rd_z(:), & ! [IN]
656  rd_rhodz(:), & ! [OUT]
657  rd_pres(:), & ! [OUT]
658  rd_presh(:), & ! [OUT]
659  rd_temp(:), & ! [OUT]
660  rd_temph(:), & ! [OUT]
661  rd_gas(:,:), & ! [OUT]
662  rd_cfc(:,:), & ! [OUT]
663  rd_aerosol_conc(:,:), & ! [OUT]
664  rd_aerosol_radi(:,:), & ! [OUT]
665  rd_cldfrac(:) ) ! [OUT]
666  endif
667 
668 !OCL XFILL
669  !$omp parallel do default(none) &
670  !$omp shared(JS,JE,IS,IE,RD_KADD,temph_merge,RD_temph,KE,RD_KMAX,KS,temp,CZ,FZ) &
671  !$omp private(i,j,k,RD_k) OMP_SCHEDULE_ collapse(2)
672  !$acc kernels
673  do j = js, je
674  do i = is, ie
675  do rd_k = 1, rd_kadd
676  temph_merge(rd_k,i,j) = rd_temph(rd_k)
677  enddo
678 
679  temph_merge(rd_kadd+1,i,j) = temp(ke,i,j)
680  do rd_k = rd_kadd+2, rd_kmax
681  k = ks + rd_kmax - rd_k ! reverse axis
682 
683  temph_merge(rd_k,i,j) = ( temp(k ,i,j) * (cz(k+1,i,j)-fz(k,i,j)) / (cz(k+1,i,j)-cz(k,i,j)) &
684  + temp(k+1,i,j) * (fz(k ,i,j)-cz(k,i,j)) / (cz(k+1,i,j)-cz(k,i,j)) )
685  enddo
686  temph_merge(rd_kmax+1,i,j) = temp(ks,i,j)
687  enddo
688  enddo
689  !$acc end kernels
690 
691 !OCL XFILL
692  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
693  !$omp private(k,i,j,RD_k) &
694  !$omp shared (RD_temp,dens,FZ,pres,temp, &
695  !$omp rhodz_merge,RD_rhodz,pres_merge,RD_pres,temp_merge, &
696  !$omp KS,JS,JE,IS,IE,RD_KMAX,RD_KADD)
697  !$acc kernels
698  do j = js, je
699  do i = is, ie
700  do rd_k = 1, rd_kadd
701  rhodz_merge(rd_k,i,j) = rd_rhodz(rd_k)
702  pres_merge(rd_k,i,j) = rd_pres(rd_k)
703  temp_merge(rd_k,i,j) = rd_temp(rd_k)
704  enddo
705 
706  do rd_k = rd_kadd+1, rd_kmax
707  k = ks + rd_kmax - rd_k ! reverse axis
708 
709  rhodz_merge(rd_k,i,j) = dens(k,i,j) * ( fz(k,i,j)-fz(k-1,i,j) ) ! [kg/m2]
710  pres_merge(rd_k,i,j) = pres(k,i,j) * 1.e-2_rp ! [hPa]
711  temp_merge(rd_k,i,j) = temp(k,i,j)
712  enddo
713  enddo
714  enddo
715  !$acc end kernels
716 
717 !OCL XFILL
718  !$omp parallel default(none) &
719  !$omp private(v,i,j,RD_k) &
720  !$omp shared (gas_merge,RD_gas, &
721  !$omp IS,IE,JS,JE,RD_KMAX)
722  !$acc kernels
723  do v = 1, mstrn_ngas
724  !$omp do OMP_SCHEDULE_ collapse(2)
725  do j = js, je
726  do i = is, ie
727  do rd_k = 1, rd_kmax
728  gas_merge(rd_k,i,j,v) = rd_gas(rd_k,v)
729  enddo
730  enddo
731  enddo
732  !$omp end do nowait
733  enddo
734  !$acc end kernels
735  !$omp end parallel
736 
737  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
738  !$omp private(k,i,j,RD_k, &
739  !$omp zerosw ) &
740  !$omp shared (gas_merge,QV,EPS,Mvap,Mdry, &
741  !$omp KS,IS,IE,JS,JE,RD_KADD,RD_KMAX)
742  !$acc kernels
743  do j = js, je
744  do i = is, ie
745  do rd_k = rd_kadd+1, rd_kmax
746  k = ks + rd_kmax - rd_k ! reverse axis
747  zerosw = sign(0.5_rp, qv(k,i,j)-eps) + 0.5_rp
748  gas_merge(rd_k,i,j,1) = qv(k,i,j) / mvap * mdry / ppm * zerosw ! [PPM]
749  enddo
750  enddo
751  enddo
752  !$acc end kernels
753 
754 !OCL XFILL
755  !$omp parallel default(none) &
756  !$omp private(v,i,j,RD_k) &
757  !$omp shared (cfc_merge,RD_cfc, &
758  !$omp IS,IE,JS,JE,RD_KMAX)
759  !$acc kernels
760  do v = 1, mstrn_ncfc
761  !$omp do OMP_SCHEDULE_ collapse(2)
762  do j = js, je
763  do i = is, ie
764  do rd_k = 1, rd_kmax
765  cfc_merge(rd_k,i,j,v) = rd_cfc(rd_k,v)
766  enddo
767  enddo
768  enddo
769  !$omp end do nowait
770  enddo
771  !$acc end kernels
772  !$omp end parallel
773 
774  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
775  !$omp private(k,i,j,RD_k) &
776  !$omp shared (cldfrac_merge,RD_cldfrac,cldfrac, &
777  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
778  !$acc kernels
779  do j = js, je
780  do i = is, ie
781  do rd_k = 1, rd_kadd
782  cldfrac_merge(rd_k,i,j) = rd_cldfrac(rd_k)
783  enddo
784 
785  do rd_k = rd_kadd+1, rd_kmax
786  k = ks + rd_kmax - rd_k ! reverse axis
787 
788  cldfrac_merge(rd_k,i,j) = min( max( cldfrac(k,i,j), 0.0_rp ), 1.0_rp )
789 ! cldfrac_merge(RD_k,i,j) = 0.5_RP + sign( 0.5_RP, cldfrac(k,i,j)-min_cldfrac )
790  enddo
791  enddo
792  enddo
793  !$acc end kernels
794 
795 !OCL XFILL
796  !$omp parallel default(none) &
797  !$omp private(v,i,j,RD_k) &
798  !$omp shared (aerosol_conc_merge,RD_aerosol_conc,aerosol_radi_merge,RD_aerosol_radi, &
799  !$omp IS,IE,JS,JE,RD_KADD)
800  !$acc kernels
801  do v = 1, rd_naero
802  !$omp do OMP_SCHEDULE_ collapse(2)
803  do j = js, je
804  do i = is, ie
805  do rd_k = 1, rd_kadd
806  aerosol_conc_merge(rd_k,i,j,v) = rd_aerosol_conc(rd_k,v)
807  aerosol_radi_merge(rd_k,i,j,v) = rd_aerosol_radi(rd_k,v)
808  enddo
809  enddo
810  enddo
811  !$omp end do nowait
812  enddo
813  !$acc end kernels
814  !$omp end parallel
815 
816  !$omp parallel default(none) &
817  !$omp private(ihydro,k,i,j,RD_k) &
818  !$omp shared (aerosol_conc_merge,tropopause,MP_Qe,HYD_DENS,RHO_std, &
819  !$omp ATMOS_PHY_RD_MSTRN_ONLY_QCI, &
820  !$omp KS,KE,IS,IE,JS,JE,RD_KMAX,RD_KADD)
821  do ihydro = 1, n_hyd
822  if ( atmos_phy_rd_mstrn_only_qci .and. &
823  ( ihydro /= i_hc .and. ihydro /= i_hi ) ) then
824 !OCL XFILL
825  !$omp do OMP_SCHEDULE_ collapse(2)
826  !$acc kernels
827  do j = js, je
828  do i = is, ie
829  do rd_k = rd_kadd+1, rd_kmax
830  aerosol_conc_merge(rd_k,i,j,ihydro) = 0.0_rp
831  end do
832  end do
833  end do
834  !$acc end kernels
835  !$omp end do nowait
836  else
837  !$omp do OMP_SCHEDULE_ collapse(2)
838  !$acc kernels
839  do j = js, je
840  do i = is, ie
841  do rd_k = rd_kadd+1, rd_kadd+1 + ke - tropopause(i,j)
842  aerosol_conc_merge(rd_k,i,j,ihydro) = 0.0_rp
843  end do
844  !$acc loop independent
845  do rd_k = rd_kadd+1 + ke - tropopause(i,j) + 1, rd_kmax
846  k = ks + rd_kmax - rd_k ! reverse axis
847  aerosol_conc_merge(rd_k,i,j,ihydro) = max( mp_qe(k,i,j,ihydro), 0.0_rp ) &
848  / hyd_dens(ihydro) * rho_std / ppm ! [PPM to standard air]
849  enddo
850  enddo
851  enddo
852  !$acc end kernels
853  !$omp end do nowait
854  end if
855  enddo
856  !$omp end parallel
857 
858  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(3) &
859  !$omp private(ihydro,k,i,j,RD_k) &
860  !$omp shared (aerosol_radi_merge,MP_Re, &
861  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
862  !$acc kernels
863  do ihydro = 1, n_hyd
864  do j = js, je
865  do i = is, ie
866  do rd_k = rd_kadd+1, rd_kmax
867  k = ks + rd_kmax - rd_k ! reverse axis
868  aerosol_radi_merge(rd_k,i,j,ihydro) = mp_re(k,i,j,ihydro)
869  enddo
870  enddo
871  enddo
872  enddo
873  !$acc end kernels
874 
875  !$omp parallel default(none) &
876  !$omp private(iaero,k,i,j,RD_k) &
877  !$omp shared (aerosol_conc_merge,RD_aerosol_conc,aerosol_radi_merge,RD_aerosol_radi, &
878  !$omp ATMOS_PHY_RD_MSTRN_USE_AERO, &
879  !$omp AE_Qe,RHO_std,AE_Re, &
880  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
881  do iaero = 1, n_ae
882 
883  if ( atmos_phy_rd_mstrn_use_aero ) then
884  !$omp do OMP_SCHEDULE_ collapse(2)
885  !$acc kernels
886  do j = js, je
887  do i = is, ie
888  do rd_k = rd_kadd+1, rd_kmax
889  k = ks + rd_kmax - rd_k ! reverse axis
890  aerosol_conc_merge(rd_k,i,j,n_hyd+iaero) = max( ae_qe(k,i,j,iaero), 0.0_rp ) &
891  / ae_dens(iaero) * rho_std / ppm ! [PPM to standard air]
892  aerosol_radi_merge(rd_k,i,j,n_hyd+iaero) = ae_re(k,i,j,iaero)
893  enddo
894  enddo
895  enddo
896  !$acc end kernels
897  !$omp end do nowait
898  else
899  !$omp do OMP_SCHEDULE_ collapse(2)
900  !$acc kernels
901  do j = js, je
902  do i = is, ie
903  do rd_k = rd_kadd+1, rd_kmax
904  aerosol_conc_merge(rd_k,i,j,n_hyd+iaero) = rd_aerosol_conc(rd_k,n_hyd+iaero)
905  aerosol_radi_merge(rd_k,i,j,n_hyd+iaero) = rd_aerosol_radi(rd_k,n_hyd+iaero)
906  enddo
907  enddo
908  enddo
909  !$acc end kernels
910  !$omp end do nowait
911  endif
912 
913  enddo
914  !$omp end parallel
915 
916  call prof_rapend ('RD_Profile', 3)
917  call prof_rapstart('RD_MSTRN_DTRN3', 3)
918 
919  ! calc radiative transfer
920 #ifndef _OPENACC
921  !$omp parallel do
922  do j = js, je
923  do i = is, ie
924 #endif
925  call rd_mstrn_dtrn3( rd_kmax, &
926 #ifdef _OPENACC
927  ia, is, ie, ja, js, je, &
928 #else
929  ia, i, i, ja, j, j, &
930  i, j, &
931 #endif
932  mstrn_ngas, mstrn_ncfc, rd_naero, & ! [IN]
933  rd_hydro_str, rd_hydro_end, rd_aero_str, rd_aero_end, & ! [IN]
934  solins(:,:), cossza(:,:), & ! [IN]
935  rhodz_merge(:,:,:), pres_merge(:,:,:), & ! [IN]
936  temp_merge(:,:,:), temph_merge(:,:,:), temp_sfc(:,:), & ! [IN]
937  gas_merge(:,:,:,:), cfc_merge(:,:,:,:), & ! [IN]
938  aerosol_conc_merge(:,:,:,:), aerosol_radi_merge(:,:,:,:), & ! [IN]
939  i_mpae2rd(:), & ! [IN]
940  cldfrac_merge(:,:,:), & ! [IN]
941  albedo_sfc(:,:,:,:), & ! [IN]
942  fact_ocean(:,:), fact_land(:,:), fact_urban(:,:), & ! [IN]
943  flux_rad_merge(:,:,:,:,:,:), flux_rad_sfc_dn(:,:,:,:), & ! [OUT]
944  taucld_067u(:,:,:), emiscld_105u(:,:,:) ) ! [OUT]
945 #ifndef _OPENACC
946  end do
947  end do
948 #endif
949 
950 
951  call prof_rapend ('RD_MSTRN_DTRN3', 3)
952 
953  ! return to grid coordinate of model domain
954  !$omp parallel default(none) &
955  !$omp private(ic,k,i,j,RD_k) &
956  !$omp shared (flux_rad,flux_rad_merge, &
957  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
958  !$acc kernels
959  !$acc loop independent collapse(4)
960  do ic = 1, 2
961  !$omp do OMP_SCHEDULE_ collapse(2)
962  do j = js, je
963  do i = is, ie
964  do rd_k = rd_kadd+1, rd_kmax+1
965  k = ks + rd_kmax - rd_k ! reverse axis
966 
967  flux_rad(k,i,j,i_lw,i_up,ic) = flux_rad_merge(rd_k,i,j,i_lw,i_up,ic)
968  flux_rad(k,i,j,i_lw,i_dn,ic) = flux_rad_merge(rd_k,i,j,i_lw,i_dn,ic)
969  flux_rad(k,i,j,i_sw,i_up,ic) = flux_rad_merge(rd_k,i,j,i_sw,i_up,ic)
970  flux_rad(k,i,j,i_sw,i_dn,ic) = flux_rad_merge(rd_k,i,j,i_sw,i_dn,ic)
971  enddo
972  enddo
973  enddo
974  !$omp end do nowait
975  enddo
976  !$acc end kernels
977  !$omp end parallel
978 
979 !OCL XFILL
980  !$omp parallel default(none) &
981  !$omp private(ic,i,j) &
982  !$omp shared (flux_rad_top,flux_rad_merge, &
983  !$omp IS,IE,JS,JE)
984  !$acc kernels
985  do ic = 1, 2
986  !$omp do OMP_SCHEDULE_
987  do j = js, je
988  do i = is, ie
989  flux_rad_top(i,j,i_lw,i_up,ic) = flux_rad_merge(1,i,j,i_lw,i_up,ic)
990  flux_rad_top(i,j,i_lw,i_dn,ic) = flux_rad_merge(1,i,j,i_lw,i_dn,ic)
991  flux_rad_top(i,j,i_sw,i_up,ic) = flux_rad_merge(1,i,j,i_sw,i_up,ic)
992  flux_rad_top(i,j,i_sw,i_dn,ic) = flux_rad_merge(1,i,j,i_sw,i_dn,ic)
993  enddo
994  enddo
995  !$omp end do nowait
996  enddo
997  !$acc end kernels
998  !$omp end parallel
999 
1000  if ( present( dtau_s ) ) then
1001  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1002  !$omp private(k,i,j,RD_k) &
1003  !$omp shared (dtau_s,tauCLD_067u, &
1004  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
1005  !$acc kernels
1006  !$acc loop independent collapse(3)
1007  do j = js, je
1008  do i = is, ie
1009  do rd_k = rd_kadd+1, rd_kmax
1010  k = ks + rd_kmax - rd_k ! reverse axis
1011  dtau_s(k,i,j) = taucld_067u(rd_k,i,j)
1012  enddo
1013  enddo
1014  enddo
1015  !$acc end kernels
1016  end if
1017 
1018  if ( present( dem_s ) ) then
1019  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1020  !$omp private(k,i,j,RD_k) &
1021  !$omp shared (dem_s,emisCLD_105u, &
1022  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
1023  !$acc kernels
1024  !$acc loop independent collapse(3)
1025  do j = js, je
1026  do i = is, ie
1027  do rd_k = rd_kadd+1, rd_kmax
1028  k = ks + rd_kmax - rd_k ! reverse axis
1029  dem_s(k,i,j) = emiscld_105u(rd_k,i,j)
1030  enddo
1031  enddo
1032  enddo
1033  !$acc end kernels
1034  end if
1035 
1036  !$acc end data
1037  !$acc end data
1038  !$acc end data
1039 
1040  return

References scale_atmos_aerosol::ae_dens, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_basepoint_lat, scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_use_climatology, scale_const::const_eps, scale_const::const_eps1, scale_const::const_grav, scale_const::const_mdry, scale_const::const_mvap, scale_const::const_pi, scale_const::const_ppm, scale_const::const_pstd, scale_const::const_rdry, scale_const::const_tem00, scale_atmos_hydrometeor::hyd_dens, scale_atmos_phy_rd_common::i_dn, scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hi, scale_atmos_phy_rd_common::i_lw, scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_atmos_phy_rd_common::i_sw, scale_atmos_phy_rd_common::i_up, scale_io::io_get_available_fid(), scale_io::io_get_fname(), scale_atmos_aerosol::n_ae, scale_atmos_hydrometeor::n_hyd, scale_cpl_sfc_index::n_rad_dir, scale_cpl_sfc_index::n_rad_rgn, scale_prc::prc_abort(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_time::time_nowdate.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver_calc_tendency().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_phy_rd_profile::atmos_phy_rd_profile_use_climatology
logical, public atmos_phy_rd_profile_use_climatology
use climatology?
Definition: scale_atmos_phy_rd_profile.F90:39
scale_atmos_aerosol::n_ae
integer, parameter, public n_ae
Definition: scale_atmos_aerosol.F90:33
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_basepoint_lat
real(rp), public atmos_grid_cartesc_real_basepoint_lat
position of base point in real world [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:37
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:99
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_mdry
real(rp), public const_mdry
mass weight (dry air) [g/mol]
Definition: scale_const.F90:58
scale_const::const_ppm
real(rp), parameter, public const_ppm
parts par million
Definition: scale_const.F90:100
scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup_zgrid
subroutine, public atmos_phy_rd_profile_setup_zgrid(KA, KS, KE, KMAX, KADD, toa, CZ, FZ, zh, z)
Setup vertical grid for radiation.
Definition: scale_atmos_phy_rd_profile.F90:1099
scale_time
module TIME
Definition: scale_time.F90:11
scale_atmos_hydrometeor::i_hc
integer, parameter, public i_hc
liquid water cloud
Definition: scale_atmos_hydrometeor.F90:97
scale_atmos_aerosol::ae_dens
real(rp), dimension(n_ae), parameter, public ae_dens
Definition: scale_atmos_aerosol.F90:53
scale_atmos_aerosol
module atmosphere / aerosol
Definition: scale_atmos_aerosol.F90:12
scale_time::time_nowdate
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:68
scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read
subroutine, public atmos_phy_rd_profile_read(kmax, ngas, ncfc, naero, real_lat, now_date, zh, z, rhodz, pres, presh, temp, temph, gas, cfc, aerosol_conc, aerosol_radi, cldfrac)
Read profile for radiation.
Definition: scale_atmos_phy_rd_profile.F90:546
scale_const::const_mvap
real(rp), public const_mvap
mass weight (water vapor) [g/mol]
Definition: scale_const.F90:67
scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup
subroutine, public atmos_phy_rd_profile_setup
Setup.
Definition: scale_atmos_phy_rd_profile.F90:135
scale_atmos_hydrometeor::hyd_dens
real(rp), dimension(n_hyd), public hyd_dens
Definition: scale_atmos_hydrometeor.F90:111
scale_atmos_phy_rd_profile
module atmosphere / physics/ radiation / profile
Definition: scale_atmos_phy_rd_profile.F90:15