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_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 -1
    ATMOS_PHY_RD_MSTRN_NRADIUS_AERO integer -1
    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 221 of file scale_atmos_phy_rd_mstrnx.F90.

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

449  !Jval )
450  use scale_const, only: &
451  eps => const_eps, &
452  mdry => const_mdry, &
453  mvap => const_mvap, &
454  ppm => const_ppm
455  use scale_time, only: &
457  use scale_atmos_grid_cartesc_real, only: &
459  use scale_atmos_hydrometeor, only: &
460  n_hyd, &
461  i_hc, &
462  i_hi, &
463  hyd_dens
464  use scale_atmos_aerosol, only: &
465  n_ae, &
466  ae_dens
467  use scale_atmos_phy_rd_profile, only: &
468  rd_profile_read => atmos_phy_rd_profile_read, &
469  rd_profile_use_climatology => atmos_phy_rd_profile_use_climatology
470  implicit none
471  integer, intent(in) :: KA, KS, KE
472  integer, intent(in) :: IA, IS, IE
473  integer, intent(in) :: JA, JS, JE
474 
475  real(RP), intent(in) :: DENS (KA,IA,JA)
476  real(RP), intent(in) :: TEMP (KA,IA,JA)
477  real(RP), intent(in) :: PRES (KA,IA,JA)
478  real(RP), intent(in) :: QV (KA,IA,JA)
479  real(RP), intent(in) :: CZ ( KA,IA,JA) ! UNUSED
480  real(RP), intent(in) :: FZ (0:KA,IA,JA)
481  real(RP), intent(in) :: fact_ocean (IA,JA)
482  real(RP), intent(in) :: fact_land (IA,JA)
483  real(RP), intent(in) :: fact_urban (IA,JA)
484  real(RP), intent(in) :: temp_sfc (IA,JA)
485  real(RP), intent(in) :: albedo_sfc (IA,JA,N_RAD_DIR,N_RAD_RGN)
486  real(RP), intent(in) :: solins (IA,JA)
487  real(RP), intent(in) :: cosSZA (IA,JA)
488  real(RP), intent(in) :: cldfrac (KA,IA,JA)
489  real(RP), intent(in) :: MP_Re (KA,IA,JA,N_HYD)
490  real(RP), intent(in) :: MP_Qe (KA,IA,JA,N_HYD)
491  real(RP), intent(in) :: AE_Re (KA,IA,JA,N_AE)
492  real(RP), intent(in) :: AE_Qe (KA,IA,JA,N_AE)
493  real(RP), intent(out) :: flux_rad (KA,IA,JA,2,2,2)
494  real(RP), intent(out) :: flux_rad_top (IA,JA,2,2,2)
495  real(RP), intent(out) :: flux_rad_sfc_dn(IA,JA,N_RAD_DIR,N_RAD_RGN)
496 
497  real(RP), intent(out), optional :: dtau_s(KA,IA,JA) ! 0.67 micron cloud optical depth
498  real(RP), intent(out), optional :: dem_s (KA,IA,JA) ! 10.5 micron cloud emissivity
499 ! real(RP), intent(out), optional :: Jval (KA,IA,JA,CH_QA_photo)
500 
501  integer :: tropopause(IA,JA)
502  real(RP) :: gamma
503 
504  real(RP), parameter :: min_cldfrac = 1.e-8_rp
505 
506  real(RP) :: rhodz_merge (RD_KMAX,IA,JA)
507  real(RP) :: pres_merge (RD_KMAX,IA,JA)
508  real(RP) :: temp_merge (RD_KMAX,IA,JA)
509  real(RP) :: temph_merge (RD_KMAX+1,IA,JA)
510 
511  real(RP) :: gas_merge (RD_KMAX,IA,JA,MSTRN_ngas)
512  real(RP) :: cfc_merge (RD_KMAX,IA,JA,MSTRN_ncfc)
513  real(RP) :: aerosol_conc_merge(RD_KMAX,IA,JA,RD_naero )
514  real(RP) :: aerosol_radi_merge(RD_KMAX,IA,JA,RD_naero )
515  real(RP) :: cldfrac_merge (RD_KMAX,IA,JA)
516 
517  ! output
518  real(RP) :: flux_rad_merge(RD_KMAX+1,IA,JA,2,2,MSTRN_ncloud)
519  real(RP) :: tauCLD_067u (RD_KMAX,IA,JA) ! 0.67 micron cloud optical depth
520  real(RP) :: emisCLD_105u (RD_KMAX,IA,JA) ! 10.5 micron cloud emissivity
521 
522  real(RP) :: zerosw
523 
524  integer :: ihydro, iaero
525  integer :: RD_k, k, i, j, v, ic
526  !---------------------------------------------------------------------------
527 
528  log_progress(*) 'atmosphere / physics / radiation / mstrnX'
529 
530  call prof_rapstart('RD_Profile', 3)
531 
532  if ( atmos_phy_rd_mstrn_only_tropocloud ) then
533  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
534  !$omp private(k,i,j, &
535  !$omp gamma) &
536  !$omp shared (tropopause,pres,temp,CZ, &
537  !$omp KS,KE,IS,IE,JS,JE)
538  do j = js, je
539  do i = is, ie
540  tropopause(i,j) = ke+1
541  do k = ke, ks, -1
542  if ( pres(k,i,j) >= 300.e+2_rp ) then
543  exit
544  elseif( pres(k,i,j) < 50.e+2_rp ) then
545  tropopause(i,j) = k
546  else
547  gamma = ( temp(k+1,i,j) - temp(k,i,j) ) / ( cz(k+1,i,j) - cz(k,i,j) )
548  if( gamma > 1.e-4_rp ) tropopause(i,j) = k
549  endif
550  enddo
551  enddo
552  enddo
553  else
554  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
555  !$omp private(i,j) &
556  !$omp shared (tropopause, &
557  !$omp KE,IS,IE,JS,JE)
558 !OCL XFILL
559  do j = js, je
560  do i = is, ie
561  tropopause(i,j) = ke+1
562  end do
563  end do
564  endif
565 
566  ! marge basic profile and value in model domain
567 
568  if ( rd_profile_use_climatology ) then
569  call rd_profile_read( rd_kmax, & ! [IN]
570  mstrn_ngas, & ! [IN]
571  mstrn_ncfc, & ! [IN]
572  rd_naero, & ! [IN]
574  time_nowdate(:), & ! [IN]
575  rd_zh(:), & ! [IN]
576  rd_z(:), & ! [IN]
577  rd_rhodz(:), & ! [OUT]
578  rd_pres(:), & ! [OUT]
579  rd_presh(:), & ! [OUT]
580  rd_temp(:), & ! [OUT]
581  rd_temph(:), & ! [OUT]
582  rd_gas(:,:), & ! [OUT]
583  rd_cfc(:,:), & ! [OUT]
584  rd_aerosol_conc(:,:), & ! [OUT]
585  rd_aerosol_radi(:,:), & ! [OUT]
586  rd_cldfrac(:) ) ! [OUT]
587  endif
588 
589 !OCL XFILL
590  !$omp parallel do default(none) &
591  !$omp shared(JS,JE,IS,IE,RD_KADD,temph_merge,RD_temph,KE,RD_KMAX,KS,temp,CZ,FZ) &
592  !$omp private(i,j,k,RD_k) OMP_SCHEDULE_ collapse(2)
593  do j = js, je
594  do i = is, ie
595  do rd_k = 1, rd_kadd
596  temph_merge(rd_k,i,j) = rd_temph(rd_k)
597  enddo
598 
599  temph_merge(rd_kadd+1,i,j) = temp(ke,i,j)
600  do rd_k = rd_kadd+2, rd_kmax
601  k = ks + rd_kmax - rd_k ! reverse axis
602 
603  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)) &
604  + temp(k+1,i,j) * (fz(k ,i,j)-cz(k,i,j)) / (cz(k+1,i,j)-cz(k,i,j)) )
605  enddo
606  temph_merge(rd_kmax+1,i,j) = temp(ks,i,j)
607  enddo
608  enddo
609 
610 !OCL XFILL
611  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
612  !$omp private(k,i,j,RD_k) &
613  !$omp shared (RD_temp,dens,FZ,pres,temp, &
614  !$omp rhodz_merge,RD_rhodz,pres_merge,RD_pres,temp_merge, &
615  !$omp KS,JS,JE,IS,IE,RD_KMAX,RD_KADD)
616  do j = js, je
617  do i = is, ie
618  do rd_k = 1, rd_kadd
619  rhodz_merge(rd_k,i,j) = rd_rhodz(rd_k)
620  pres_merge(rd_k,i,j) = rd_pres(rd_k)
621  temp_merge(rd_k,i,j) = rd_temp(rd_k)
622  enddo
623 
624  do rd_k = rd_kadd+1, rd_kmax
625  k = ks + rd_kmax - rd_k ! reverse axis
626 
627  rhodz_merge(rd_k,i,j) = dens(k,i,j) * ( fz(k,i,j)-fz(k-1,i,j) ) ! [kg/m2]
628  pres_merge(rd_k,i,j) = pres(k,i,j) * 1.e-2_rp ! [hPa]
629  temp_merge(rd_k,i,j) = temp(k,i,j)
630  enddo
631  enddo
632  enddo
633 
634 !OCL XFILL
635 !OCL SERIAL
636  do v = 1, mstrn_ngas
637 !OCL PARALLEL
638  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
639  !$omp private(i,j,RD_k) &
640  !$omp shared (gas_merge,RD_gas, &
641  !$omp IS,IE,JS,JE,RD_KMAX,v)
642  do j = js, je
643  do i = is, ie
644  do rd_k = 1, rd_kmax
645  gas_merge(rd_k,i,j,v) = rd_gas(rd_k,v)
646  enddo
647  enddo
648  enddo
649  enddo
650 
651  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
652  !$omp private(k,i,j,RD_k, &
653  !$omp zerosw ) &
654  !$omp shared (gas_merge,QV,EPS,Mvap,Mdry, &
655  !$omp KS,IS,IE,JS,JE,RD_KADD,RD_KMAX)
656  do j = js, je
657  do i = is, ie
658  do rd_k = rd_kadd+1, rd_kmax
659  k = ks + rd_kmax - rd_k ! reverse axis
660  zerosw = sign(0.5_rp, qv(k,i,j)-eps) + 0.5_rp
661  gas_merge(rd_k,i,j,1) = qv(k,i,j) / mvap * mdry / ppm * zerosw ! [PPM]
662  enddo
663  enddo
664  enddo
665 
666 !OCL XFILL
667 !OCL SERIAL
668  do v = 1, mstrn_ncfc
669  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
670  !$omp private(i,j,RD_k) &
671  !$omp shared (cfc_merge,RD_cfc, &
672  !$omp IS,IE,JS,JE,RD_KMAX,v)
673 !OCL PARALLEL
674  do j = js, je
675  do i = is, ie
676  do rd_k = 1, rd_kmax
677  cfc_merge(rd_k,i,j,v) = rd_cfc(rd_k,v)
678  enddo
679  enddo
680  enddo
681  enddo
682 
683  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
684  !$omp private(k,i,j,RD_k) &
685  !$omp shared (cldfrac_merge,RD_cldfrac,cldfrac, &
686  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
687  do j = js, je
688  do i = is, ie
689  do rd_k = 1, rd_kadd
690  cldfrac_merge(rd_k,i,j) = rd_cldfrac(rd_k)
691  enddo
692 
693  do rd_k = rd_kadd+1, rd_kmax
694  k = ks + rd_kmax - rd_k ! reverse axis
695 
696  cldfrac_merge(rd_k,i,j) = min( max( cldfrac(k,i,j), 0.0_rp ), 1.0_rp )
697 ! cldfrac_merge(RD_k,i,j) = 0.5_RP + sign( 0.5_RP, cldfrac(k,i,j)-min_cldfrac )
698  enddo
699  enddo
700  enddo
701 
702 !OCL XFILL
703 !OCL SERIAL
704  do v = 1, rd_naero
705  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
706  !$omp private(i,j,RD_k) &
707  !$omp shared (aerosol_conc_merge,RD_aerosol_conc,aerosol_radi_merge,RD_aerosol_radi, &
708  !$omp IS,IE,JS,JE,RD_KADD,v)
709 !OCL PARALLEL
710  do j = js, je
711  do i = is, ie
712  do rd_k = 1, rd_kadd
713  aerosol_conc_merge(rd_k,i,j,v) = rd_aerosol_conc(rd_k,v)
714  aerosol_radi_merge(rd_k,i,j,v) = rd_aerosol_radi(rd_k,v)
715  enddo
716  enddo
717  enddo
718  enddo
719 
720 !OCL SERIAL
721  do ihydro = 1, n_hyd
722  if ( atmos_phy_rd_mstrn_only_qci .and. &
723  ( ihydro /= i_hc .and. ihydro /= i_hi ) ) then
724 !OCL XFILL
725  aerosol_conc_merge(:,:,:,ihydro) = 0.0_rp
726  cycle
727  end if
728  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
729  !$omp private(k,i,j,RD_k) &
730  !$omp shared (aerosol_conc_merge,tropopause,MP_Qe,HYD_DENS,RHO_std, &
731  !$omp KS,KE,IS,IE,JS,JE,RD_KMAX,RD_KADD,ihydro)
732 !OCL PARALLEL
733  do j = js, je
734  do i = is, ie
735  do rd_k = rd_kadd+1, rd_kadd+1 + ke - tropopause(i,j)
736  aerosol_conc_merge(rd_k,i,j,ihydro) = 0.0_rp
737  end do
738  do rd_k = rd_kadd+1 + ke - tropopause(i,j) + 1, rd_kmax
739  k = ks + rd_kmax - rd_k ! reverse axis
740  aerosol_conc_merge(rd_k,i,j,ihydro) = max( mp_qe(k,i,j,ihydro), 0.0_rp ) &
741  / hyd_dens(ihydro) * rho_std / ppm ! [PPM to standard air]
742  enddo
743  enddo
744  enddo
745  enddo
746 
747 !OCL SERIAL
748  do ihydro = 1, n_hyd
749  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
750  !$omp private(k,i,j,RD_k) &
751  !$omp shared (aerosol_radi_merge,MP_Re, &
752  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD,ihydro)
753 !OCL PARALLEL
754  do j = js, je
755  do i = is, ie
756  do rd_k = rd_kadd+1, rd_kmax
757  k = ks + rd_kmax - rd_k ! reverse axis
758  aerosol_radi_merge(rd_k,i,j,ihydro) = mp_re(k,i,j,ihydro)
759  enddo
760  enddo
761  enddo
762  enddo
763 
764 !OCL SERIAL
765  do iaero = 1, n_ae
766 
767  if ( atmos_phy_rd_mstrn_use_aero ) then
768  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
769  !$omp private(k,i,j,RD_k) &
770  !$omp shared (aerosol_conc_merge,aerosol_radi_merge,AE_Qe,RHO_std,AE_Re, &
771  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD,iaero)
772 !OCL PARALLEL
773  do j = js, je
774  do i = is, ie
775  do rd_k = rd_kadd+1, rd_kmax
776  k = ks + rd_kmax - rd_k ! reverse axis
777  aerosol_conc_merge(rd_k,i,j,n_hyd+iaero) = max( ae_qe(k,i,j,iaero), 0.0_rp ) &
778  / ae_dens(iaero) * rho_std / ppm ! [PPM to standard air]
779  aerosol_radi_merge(rd_k,i,j,n_hyd+iaero) = ae_re(k,i,j,iaero)
780  enddo
781  enddo
782  enddo
783  else
784  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
785  !$omp private(i,j,RD_k) &
786  !$omp shared (aerosol_conc_merge,RD_aerosol_conc,aerosol_radi_merge,RD_aerosol_radi, &
787  !$omp IS,IE,JS,JE,RD_KMAX,RD_KADD,iaero)
788 !OCL PARALLEL
789  do j = js, je
790  do i = is, ie
791  do rd_k = rd_kadd+1, rd_kmax
792  aerosol_conc_merge(rd_k,i,j,n_hyd+iaero) = rd_aerosol_conc(rd_k,n_hyd+iaero)
793  aerosol_radi_merge(rd_k,i,j,n_hyd+iaero) = rd_aerosol_radi(rd_k,n_hyd+iaero)
794  enddo
795  enddo
796  enddo
797  endif
798 
799  enddo
800 
801  call prof_rapend ('RD_Profile', 3)
802  call prof_rapstart('RD_MSTRN_DTRN3', 3)
803 
804  ! calc radiative transfer
805  call rd_mstrn_dtrn3( rd_kmax, ia, is, ie, ja, js, je, & ! [IN]
806  mstrn_ngas, mstrn_ncfc, rd_naero, & ! [IN]
807  rd_hydro_str, rd_hydro_end, rd_aero_str, rd_aero_end, & ! [IN]
808  solins(:,:), cossza(:,:), & ! [IN]
809  rhodz_merge(:,:,:), pres_merge(:,:,:), & ! [IN]
810  temp_merge(:,:,:), temph_merge(:,:,:), temp_sfc(:,:), & ! [IN]
811  gas_merge(:,:,:,:), cfc_merge(:,:,:,:), & ! [IN]
812  aerosol_conc_merge(:,:,:,:), aerosol_radi_merge(:,:,:,:), & ! [IN]
813  i_mpae2rd(:), & ! [IN]
814  cldfrac_merge(:,:,:), & ! [IN]
815  albedo_sfc(:,:,:,:), & ! [IN]
816  fact_ocean(:,:), fact_land(:,:), fact_urban(:,:), & ! [IN]
817  flux_rad_merge(:,:,:,:,:,:), flux_rad_sfc_dn(:,:,:,:), & ! [OUT]
818  taucld_067u(:,:,:), emiscld_105u(:,:,:) ) ! [OUT]
819 
820  call prof_rapend ('RD_MSTRN_DTRN3', 3)
821 
822  ! return to grid coordinate of model domain
823 !OCL SERIAL
824  do ic = 1, 2
825  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
826  !$omp private(k,i,j,RD_k) &
827  !$omp shared (flux_rad,flux_rad_merge, &
828  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD,ic)
829 !OCL PARALLEL
830  do j = js, je
831  do i = is, ie
832  do rd_k = rd_kadd+1, rd_kmax+1
833  k = ks + rd_kmax - rd_k ! reverse axis
834 
835  flux_rad(k,i,j,i_lw,i_up,ic) = flux_rad_merge(rd_k,i,j,i_lw,i_up,ic)
836  flux_rad(k,i,j,i_lw,i_dn,ic) = flux_rad_merge(rd_k,i,j,i_lw,i_dn,ic)
837  flux_rad(k,i,j,i_sw,i_up,ic) = flux_rad_merge(rd_k,i,j,i_sw,i_up,ic)
838  flux_rad(k,i,j,i_sw,i_dn,ic) = flux_rad_merge(rd_k,i,j,i_sw,i_dn,ic)
839  enddo
840  enddo
841  enddo
842  enddo
843 
844 !OCL XFILL
845 !OCL SERIAL
846  do ic = 1, 2
847  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
848  !$omp private(i,j) &
849  !$omp shared (flux_rad_top,flux_rad_merge, &
850  !$omp IS,IE,JS,JE,ic)
851 !OCL PARALLEL
852  do j = js, je
853  do i = is, ie
854  flux_rad_top(i,j,i_lw,i_up,ic) = flux_rad_merge(1,i,j,i_lw,i_up,ic)
855  flux_rad_top(i,j,i_lw,i_dn,ic) = flux_rad_merge(1,i,j,i_lw,i_dn,ic)
856  flux_rad_top(i,j,i_sw,i_up,ic) = flux_rad_merge(1,i,j,i_sw,i_up,ic)
857  flux_rad_top(i,j,i_sw,i_dn,ic) = flux_rad_merge(1,i,j,i_sw,i_dn,ic)
858  enddo
859  enddo
860  enddo
861 
862  if ( present( dtau_s ) ) then
863  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
864  !$omp private(k,i,j,RD_k) &
865  !$omp shared (dtau_s,tauCLD_067u, &
866  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
867  do j = js, je
868  do i = is, ie
869  do rd_k = rd_kadd+1, rd_kmax
870  k = ks + rd_kmax - rd_k ! reverse axis
871  dtau_s(k,i,j) = taucld_067u(rd_k,i,j)
872  enddo
873  enddo
874  enddo
875  end if
876 
877  if ( present( dem_s ) ) then
878  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
879  !$omp private(k,i,j,RD_k) &
880  !$omp shared (dem_s,emisCLD_105u, &
881  !$omp KS,IS,IE,JS,JE,RD_KMAX,RD_KADD)
882  do j = js, je
883  do i = is, ie
884  do rd_k = rd_kadd+1, rd_kmax
885  k = ks + rd_kmax - rd_k ! reverse axis
886  dem_s(k,i,j) = emiscld_105u(rd_k,i,j)
887  enddo
888  enddo
889  enddo
890  end if
891 
892  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_atmos_aerosol::n_ae, scale_atmos_hydrometeor::n_hyd, 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:342
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:38
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:33
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:36
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:83
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:54
scale_const::const_ppm
real(rp), parameter, public const_ppm
parts par million
Definition: scale_const.F90:91
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:1086
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:81
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:66
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:536
scale_const::const_mvap
real(rp), public const_mvap
mass weight (water vapor) [g/mol]
Definition: scale_const.F90:62
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:134
scale_atmos_hydrometeor::hyd_dens
real(rp), dimension(n_hyd), public hyd_dens
Definition: scale_atmos_hydrometeor.F90:95
scale_atmos_phy_rd_profile
module atmosphere / physics/ radiation / profile
Definition: scale_atmos_phy_rd_profile.F90:15