SCALE-RM
Functions/Subroutines
scale_atmos_phy_rd_mstrnx Module Reference

module ATMOSPHERE / Physics Radiation More...

Functions/Subroutines

subroutine, public atmos_phy_rd_mstrnx_setup (RD_TYPE)
 Setup. More...
 
subroutine, public atmos_phy_rd_mstrnx (DENS, RHOT, QTRC, CZ, FZ, fact_ocean, fact_land, fact_urban, temp_sfc, albedo_land, solins, cosSZA, flux_rad, flux_rad_top, flux_rad_sfc_dn)
 Radiation main. More...
 

Detailed Description

module ATMOSPHERE / Physics Radiation

Description
Atmospheric radiation transfer process mstrnX Ref: Nakajima and Tanaka(1986) Nakajima et al.(2000) Sekiguchi and Nakajima(2008)
Author
Team SCALE
History
  • 2013-02-06 (H.Yashiro) [new]

Function/Subroutine Documentation

◆ atmos_phy_rd_mstrnx_setup()

subroutine, public scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_setup ( character(len=*), intent(in)  RD_TYPE)

Setup.

Definition at line 221 of file scale_atmos_phy_rd_mstrnx.F90.

References scale_tracer::ae_qa, 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_tracer::i_ae2rd, scale_tracer::i_mp2rd, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::kmax, scale_tracer::mp_qa, scale_process::prc_mpistop(), scale_grid_real::real_basepoint_lat, and scale_time::time_nowdate.

Referenced by scale_atmos_phy_rd::atmos_phy_rd_setup().

221  use scale_process, only: &
223  use scale_time, only: &
225  use scale_grid_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  implicit none
232 
233  character(len=*), intent(in) :: rd_type
234 
235  real(RP) :: atmos_phy_rd_mstrn_toa
236  integer :: atmos_phy_rd_mstrn_kadd
237  character(len=H_LONG) :: atmos_phy_rd_mstrn_gaspara_in_filename
238  character(len=H_LONG) :: atmos_phy_rd_mstrn_aeropara_in_filename
239  character(len=H_LONG) :: atmos_phy_rd_mstrn_hygropara_in_filename
240  integer :: atmos_phy_rd_mstrn_nband
241  integer :: atmos_phy_rd_mstrn_nptype
242  integer :: atmos_phy_rd_mstrn_nradius
243 
244  namelist / param_atmos_phy_rd_mstrn / &
245  atmos_phy_rd_mstrn_toa, &
246  atmos_phy_rd_mstrn_kadd, &
247  atmos_phy_rd_mstrn_gaspara_in_filename, &
248  atmos_phy_rd_mstrn_aeropara_in_filename, &
249  atmos_phy_rd_mstrn_hygropara_in_filename, &
250  atmos_phy_rd_mstrn_nband, &
251  atmos_phy_rd_mstrn_nptype, &
252  atmos_phy_rd_mstrn_nradius, &
253  atmos_phy_rd_mstrn_only_qci
254 
255  integer :: ngas, ncfc
256  integer :: ihydro, iaero
257  integer :: ierr
258  !---------------------------------------------------------------------------
259 
260  if( io_l ) write(io_fid_log,*)
261  if( io_l ) write(io_fid_log,*) '++++++ Module[RADIATION] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
262  if( io_l ) write(io_fid_log,*) '+++ MstrnX radiation process'
263 
264  if ( rd_type /= 'MSTRNX' ) then
265  write(*,*) 'xxx RD_TYPE is not MSTRNX. Check!'
266  call prc_mpistop
267  endif
268 
269  !--- read namelist
270  atmos_phy_rd_mstrn_toa = rd_toa
271  atmos_phy_rd_mstrn_kadd = rd_kadd
272  atmos_phy_rd_mstrn_gaspara_in_filename = mstrn_gaspara_inputfile
273  atmos_phy_rd_mstrn_aeropara_in_filename = mstrn_aeropara_inputfile
274  atmos_phy_rd_mstrn_hygropara_in_filename = mstrn_hygropara_inputfile
275  atmos_phy_rd_mstrn_nband = mstrn_nband
276  atmos_phy_rd_mstrn_nptype = mstrn_nptype
277  atmos_phy_rd_mstrn_nradius = mstrn_nradius
278 
279  rewind(io_fid_conf)
280  read(io_fid_conf,nml=param_atmos_phy_rd_mstrn,iostat=ierr)
281  if( ierr < 0 ) then !--- missing
282  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
283  elseif( ierr > 0 ) then !--- fatal error
284  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_RD_MSTRN. Check!'
285  call prc_mpistop
286  endif
287  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_rd_mstrn)
288 
289  rd_toa = atmos_phy_rd_mstrn_toa
290  rd_kadd = atmos_phy_rd_mstrn_kadd
291  mstrn_gaspara_inputfile = atmos_phy_rd_mstrn_gaspara_in_filename
292  mstrn_aeropara_inputfile = atmos_phy_rd_mstrn_aeropara_in_filename
293  mstrn_hygropara_inputfile = atmos_phy_rd_mstrn_hygropara_in_filename
294  mstrn_nband = atmos_phy_rd_mstrn_nband
295  mstrn_nptype = atmos_phy_rd_mstrn_nptype
296  mstrn_nradius = atmos_phy_rd_mstrn_nradius
297 
298  !--- setup MSTRN parameter
299  call rd_mstrn_setup( ngas, & ! [OUT]
300  ncfc ) ! [OUT]
301 
302  !--- setup climatological profile
303  call rd_profile_setup
304 
305  rd_kmax = kmax + rd_kadd
306  rd_naero = mp_qa + ae_qa
307  rd_hydro_str = 1
308  rd_hydro_end = mp_qa
309  rd_aero_str = mp_qa + 1
310  rd_aero_end = mp_qa + ae_qa
311 
312  !--- allocate arrays
313  ! input
314  allocate( rd_zh(rd_kmax+1) )
315  allocate( rd_z(rd_kmax ) )
316 
317  allocate( rd_rhodz(rd_kmax ) )
318  allocate( rd_pres(rd_kmax ) )
319  allocate( rd_presh(rd_kmax+1) )
320  allocate( rd_temp(rd_kmax ) )
321  allocate( rd_temph(rd_kmax+1) )
322 
323  allocate( rd_gas(rd_kmax,ngas ) )
324  allocate( rd_cfc(rd_kmax,ncfc ) )
325  allocate( rd_aerosol_conc(rd_kmax,rd_naero) )
326  allocate( rd_aerosol_radi(rd_kmax,rd_naero) )
327  allocate( rd_cldfrac(rd_kmax ) )
328 
329  allocate( i_mpae2rd(rd_naero) )
330 
331  ! make look-up table between hydrometeor tracer and mstrn particle type
332  do ihydro = 1, mp_qa
333  i_mpae2rd(ihydro) = i_mp2rd(ihydro)
334  enddo
335  ! make look-up table between aerosol tracer and mstrn particle type
336  do iaero = 1, ae_qa
337  i_mpae2rd(mp_qa+iaero) = i_ae2rd(iaero)
338  enddo
339 
340  !--- setup vartical grid for radiation (larger TOA than Model domain)
341  call rd_profile_setup_zgrid( rd_toa, rd_kmax, rd_kadd, & ! [IN]
342  rd_zh(:), rd_z(:) ) ! [INOUT]
343 
344  !--- read climatological profile
345  call rd_profile_read( rd_kmax, & ! [IN]
346  ngas, & ! [IN]
347  ncfc, & ! [IN]
348  rd_naero, & ! [IN]
349  real_basepoint_lat, & ! [IN]
350  time_nowdate(:), & ! [IN]
351  rd_zh(:), & ! [IN]
352  rd_z(:), & ! [IN]
353  rd_rhodz(:), & ! [OUT]
354  rd_pres(:), & ! [OUT]
355  rd_presh(:), & ! [OUT]
356  rd_temp(:), & ! [OUT]
357  rd_temph(:), & ! [OUT]
358  rd_gas(:,:), & ! [OUT]
359  rd_cfc(:,:), & ! [OUT]
360  rd_aerosol_conc(:,:), & ! [OUT]
361  rd_aerosol_radi(:,:), & ! [OUT]
362  rd_cldfrac(:) ) ! [OUT]
363 
364  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), public real_basepoint_lat
position of base point in real world [rad,-pi,pi]
module ATMOSPHERE / Physics Radiation / Vertical profile
subroutine, public atmos_phy_rd_profile_setup
Setup.
module GRID (real space)
integer, public mp_qa
integer, public kmax
of computational cells: z
module TIME
Definition: scale_time.F90:15
subroutine, public atmos_phy_rd_profile_setup_zgrid(toa, kmax, kadd, zh, z)
Setup vertical grid for radiation.
module PROCESS
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.
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_rd_mstrnx()

subroutine, public scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx ( real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  RHOT,
real(rp), dimension (ka,ia,ja,qa), intent(in)  QTRC,
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,2), intent(in)  albedo_land,
real(rp), dimension (ia,ja), intent(in)  solins,
real(rp), dimension (ia,ja), intent(in)  cosSZA,
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,2,2), intent(out)  flux_rad_sfc_dn 
)

Radiation main.

Definition at line 380 of file scale_atmos_phy_rd_mstrnx.F90.

References scale_atmos_phy_ae::ae_dens, scale_tracer::ae_qa, scale_atmos_phy_ae::atmos_phy_ae_effectiveradius, scale_atmos_phy_mp::atmos_phy_mp_cloudfraction, scale_atmos_phy_mp::atmos_phy_mp_dens, scale_atmos_phy_mp::atmos_phy_mp_effectiveradius, scale_atmos_phy_mp::atmos_phy_mp_mixingratio, 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_huge, 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_tracer::i_ae2all, scale_atmos_phy_rd_common::i_dn, scale_atmos_phy_rd_common::i_lw, scale_tracer::i_mp2all, scale_tracer::i_qc, scale_tracer::i_qi, scale_tracer::i_qv, scale_atmos_phy_rd_common::i_sw, scale_atmos_phy_rd_common::i_up, scale_grid_index::ie, scale_grid_index::imax, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::jmax, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, scale_tracer::mp_qa, scale_process::prc_mpistop(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_grid_real::real_basepoint_lat, and scale_time::time_nowdate.

Referenced by scale_atmos_phy_rd::atmos_phy_rd_setup().

380 ! Jval )
381  use scale_grid_index
382  use scale_tracer
383  use scale_const, only: &
384  eps => const_eps, &
385  mdry => const_mdry, &
386  mvap => const_mvap, &
387  ppm => const_ppm
388  use scale_time, only: &
390  use scale_grid_real, only: &
392  use scale_atmos_thermodyn, only: &
393  thermodyn_temp_pres => atmos_thermodyn_temp_pres
394  use scale_atmos_saturation, only: &
395  saturation_dens2qsat_liq => atmos_saturation_dens2qsat_liq
396  use scale_atmos_phy_mp, only: &
397  mp_effectiveradius => atmos_phy_mp_effectiveradius, &
398  mp_cloudfraction => atmos_phy_mp_cloudfraction, &
399  mp_mixingratio => atmos_phy_mp_mixingratio, &
400  mp_dens => atmos_phy_mp_dens
401  use scale_atmos_phy_ae, only: &
402  ae_effectiveradius => atmos_phy_ae_effectiveradius, &
403  ae_dens
404  use scale_atmos_phy_rd_profile, only: &
405  rd_profile_read => atmos_phy_rd_profile_read, &
406  rd_profile_use_climatology => atmos_phy_rd_profile_use_climatology
407  implicit none
408 
409  real(RP), intent(in) :: dens (ka,ia,ja)
410  real(RP), intent(in) :: rhot (ka,ia,ja)
411  real(RP), intent(in) :: qtrc (ka,ia,ja,qa)
412  real(RP), intent(in) :: cz ( ka,ia,ja) ! UNUSED
413  real(RP), intent(in) :: fz (0:ka,ia,ja)
414  real(RP), intent(in) :: fact_ocean (ia,ja)
415  real(RP), intent(in) :: fact_land (ia,ja)
416  real(RP), intent(in) :: fact_urban (ia,ja)
417  real(RP), intent(in) :: temp_sfc (ia,ja)
418  real(RP), intent(in) :: albedo_land (ia,ja,2)
419  real(RP), intent(in) :: solins (ia,ja)
420  real(RP), intent(in) :: cossza (ia,ja)
421  real(RP), intent(out) :: flux_rad (ka,ia,ja,2,2,2)
422  real(RP), intent(out) :: flux_rad_top (ia,ja,2,2,2)
423  real(RP), intent(out) :: flux_rad_sfc_dn(ia,ja,2,2)
424 ! real(RP), intent(out) :: Jval (KA,IA,JA,CH_QA_photo)
425 
426  real(RP) :: temp (ka,ia,ja)
427  real(RP) :: pres (ka,ia,ja)
428  real(RP) :: qsat (ka,ia,ja)
429  real(RP) :: rh (ka,ia,ja)
430  real(RP) :: cldfrac(ka,ia,ja)
431  real(RP) :: mp_re (ka,ia,ja,mp_qa)
432  real(RP) :: mp_qe (ka,ia,ja,mp_qa)
433  real(RP) :: ae_re (ka,ia,ja,ae_qa)
434 
435  real(RP), parameter :: min_cldfrac = 1.e-8_rp
436 
437  real(RP) :: rhodz_merge (rd_kmax,ia,ja)
438  real(RP) :: pres_merge (rd_kmax,ia,ja)
439  real(RP) :: temp_merge (rd_kmax,ia,ja)
440  real(RP) :: temph_merge (rd_kmax+1,ia,ja)
441 
442  real(RP) :: gas_merge (rd_kmax,ia,ja,mstrn_ngas)
443  real(RP) :: cfc_merge (rd_kmax,ia,ja,mstrn_ncfc)
444  real(RP) :: aerosol_conc_merge(rd_kmax,ia,ja,rd_naero )
445  real(RP) :: aerosol_radi_merge(rd_kmax,ia,ja,rd_naero )
446  real(RP) :: cldfrac_merge (rd_kmax,ia,ja)
447 
448  ! output
449  real(RP) :: flux_rad_merge(rd_kmax+1,ia,ja,2,2,mstrn_ncloud)
450 
451  real(RP) :: zerosw
452 
453  integer :: ihydro, iaero, iq
454  integer :: rd_k, k, i, j, v, ic
455  !---------------------------------------------------------------------------
456 
457  if( io_l ) write(io_fid_log,*) '*** Physics step: Radiation(mstrnX)'
458 
459  call prof_rapstart('RD_Profile', 3)
460 
461  call thermodyn_temp_pres( temp(:,:,:), & ! [OUT]
462  pres(:,:,:), & ! [OUT]
463  dens(:,:,:), & ! [IN]
464  rhot(:,:,:), & ! [IN]
465  qtrc(:,:,:,:) ) ! [IN]
466 
467  call saturation_dens2qsat_liq( qsat(:,:,:), & ! [OUT]
468  temp(:,:,:), & ! [IN]
469  dens(:,:,:) ) ! [IN]
470 
471 !OCL XFILL
472  do j = js, je
473  do i = is, ie
474  do k = ks, ke
475  rh(k,i,j) = qtrc(k,i,j,i_qv) / qsat(k,i,j)
476  enddo
477  enddo
478  enddo
479 
480  call mp_cloudfraction( cldfrac(:,:,:), & ! [OUT]
481  qtrc(:,:,:,:), & ! [IN]
482  eps ) ! [IN]
483 
484  call mp_effectiveradius( mp_re(:,:,:,:), & ! [OUT]
485  qtrc(:,:,:,:), & ! [IN]
486  dens(:,:,:) , & ! [IN]
487  temp(:,:,:) ) ! [IN]
488 
489  call ae_effectiveradius( ae_re(:,:,:,:), & ! [OUT]
490  qtrc(:,:,:,:), & ! [IN]
491  rh(:,:,:) ) ! [IN]
492 
493  call mp_mixingratio( mp_qe(:,:,:,:), & ! [OUT]
494  qtrc(:,:,:,:) ) ! [IN]
495 
496  if ( atmos_phy_rd_mstrn_only_qci ) then
497  do ihydro = 1, mp_qa
498  iq = i_mp2all(ihydro)
499  if ( iq /= i_qc .AND. iq /= i_qi ) then
500 !OCL XFILL
501  mp_qe(:,:,:,ihydro) = 0.0_rp
502  end if
503  end do
504  end if
505 
506  ! marge basic profile and value in model domain
507 
508  if ( rd_profile_use_climatology ) then
509  call rd_profile_read( rd_kmax, & ! [IN]
510  mstrn_ngas, & ! [IN]
511  mstrn_ncfc, & ! [IN]
512  rd_naero, & ! [IN]
513  real_basepoint_lat, & ! [IN]
514  time_nowdate(:), & ! [IN]
515  rd_zh(:), & ! [IN]
516  rd_z(:), & ! [IN]
517  rd_rhodz(:), & ! [OUT]
518  rd_pres(:), & ! [OUT]
519  rd_presh(:), & ! [OUT]
520  rd_temp(:), & ! [OUT]
521  rd_temph(:), & ! [OUT]
522  rd_gas(:,:), & ! [OUT]
523  rd_cfc(:,:), & ! [OUT]
524  rd_aerosol_conc(:,:), & ! [OUT]
525  rd_aerosol_radi(:,:), & ! [OUT]
526  rd_cldfrac(:) ) ! [OUT]
527  endif
528 
529 !OCL XFILL
530  do j = js, je
531  do i = is, ie
532  do rd_k = 1, rd_kadd
533  temph_merge(rd_k,i,j) = rd_temph(rd_k)
534  enddo
535 
536  temp(ke+1,i,j) = temp(ke,i,j)
537  do rd_k = rd_kadd+1, rd_kmax
538  k = ks + rd_kmax - rd_k ! reverse axis
539 
540  temph_merge(rd_k,i,j) = 0.5_rp * ( temp(k,i,j) + temp(k+1,i,j) )
541  enddo
542  temph_merge(rd_kmax+1,i,j) = temp(ks,i,j)
543  enddo
544  enddo
545 
546 !OCL XFILL
547  do j = js, je
548  do i = is, ie
549  do rd_k = 1, rd_kadd
550  rhodz_merge(rd_k,i,j) = rd_rhodz(rd_k)
551  pres_merge(rd_k,i,j) = rd_pres(rd_k)
552  temp_merge(rd_k,i,j) = rd_temp(rd_k)
553  enddo
554 
555  do rd_k = rd_kadd+1, rd_kmax
556  k = ks + rd_kmax - rd_k ! reverse axis
557 
558  rhodz_merge(rd_k,i,j) = dens(k,i,j) * ( fz(k,i,j)-fz(k-1,i,j) ) ! [kg/m2]
559  pres_merge(rd_k,i,j) = pres(k,i,j) * 1.e-2_rp ! [hPa]
560  temp_merge(rd_k,i,j) = temp(k,i,j)
561  enddo
562  enddo
563  enddo
564 
565 !OCL XFILL
566  do v = 1, mstrn_ngas
567  do j = js, je
568  do i = is, ie
569  do rd_k = 1, rd_kmax
570  gas_merge(rd_k,i,j,v) = rd_gas(rd_k,v)
571  enddo
572  enddo
573  enddo
574  enddo
575 
576  do j = js, je
577  do i = is, ie
578  do rd_k = rd_kadd+1, rd_kmax
579  k = ks + rd_kmax - rd_k ! reverse axis
580  zerosw = sign(0.5_rp, qtrc(k,i,j,i_qv)-eps) + 0.5_rp
581  gas_merge(rd_k,i,j,1) = qtrc(k,i,j,i_qv) / mvap * mdry / ppm * zerosw ! [PPM]
582  enddo
583  enddo
584  enddo
585 
586 !OCL XFILL
587  do v = 1, mstrn_ncfc
588  do j = js, je
589  do i = is, ie
590  do rd_k = 1, rd_kmax
591  cfc_merge(rd_k,i,j,v) = rd_cfc(rd_k,v)
592  enddo
593  enddo
594  enddo
595  enddo
596 
597  do j = js, je
598  do i = is, ie
599  do rd_k = 1, rd_kadd
600  cldfrac_merge(rd_k,i,j) = rd_cldfrac(rd_k)
601  enddo
602 
603  do rd_k = rd_kadd+1, rd_kmax
604  k = ks + rd_kmax - rd_k ! reverse axis
605 
606  cldfrac_merge(rd_k,i,j) = 0.5_rp + sign( 0.5_rp, cldfrac(k,i,j)-min_cldfrac )
607  enddo
608  enddo
609  enddo
610 
611 !OCL XFILL
612  do v = 1, rd_naero
613  do j = js, je
614  do i = is, ie
615  do rd_k = 1, rd_kadd
616  aerosol_conc_merge(rd_k,i,j,v) = rd_aerosol_conc(rd_k,v)
617  aerosol_radi_merge(rd_k,i,j,v) = rd_aerosol_radi(rd_k,v)
618  enddo
619  enddo
620  enddo
621  enddo
622 
623  do ihydro = 1, mp_qa
624  do j = js, je
625  do i = is, ie
626  do rd_k = rd_kadd+1, rd_kmax
627  k = ks + rd_kmax - rd_k ! reverse axis
628 
629  aerosol_conc_merge(rd_k,i,j,ihydro) = max( mp_qe(k,i,j,ihydro), 0.0_rp ) &
630  / mp_dens(ihydro) * rho_std / ppm ! [PPM to standard air]
631  aerosol_radi_merge(rd_k,i,j,ihydro) = mp_re(k,i,j,ihydro)
632  enddo
633  enddo
634  enddo
635  enddo
636 
637  do iaero = 1, ae_qa
638  if ( i_ae2all(iaero) > 0 ) then
639  do j = js, je
640  do i = is, ie
641  do rd_k = rd_kadd+1, rd_kmax
642  k = ks + rd_kmax - rd_k ! reverse axis
643 
644  aerosol_conc_merge(rd_k,i,j,mp_qa+iaero) = max( qtrc(k,i,j,i_ae2all(iaero)), 0.0_rp ) &
645  / ae_dens(iaero) * rho_std / ppm ! [PPM to standard air]
646  aerosol_radi_merge(rd_k,i,j,mp_qa+iaero) = ae_re(k,i,j,iaero)
647  enddo
648  enddo
649  enddo
650  else
651 !OCL XFILL
652  do j = js, je
653  do i = is, ie
654  do rd_k = rd_kadd+1, rd_kmax
655  aerosol_conc_merge(rd_k,i,j,mp_qa+iaero) = rd_aerosol_conc(rd_k,mp_qa+iaero)
656  aerosol_radi_merge(rd_k,i,j,mp_qa+iaero) = rd_aerosol_radi(rd_k,mp_qa+iaero)
657  enddo
658  enddo
659  enddo
660  endif
661  enddo
662 
663  call prof_rapend ('RD_Profile', 3)
664  call prof_rapstart('RD_MSTRN_DTRN3', 3)
665 
666  ! calc radiative transfer
667  call rd_mstrn_dtrn3( rd_kmax, & ! [IN]
668  mstrn_ngas, & ! [IN]
669  mstrn_ncfc, & ! [IN]
670  rd_naero, & ! [IN]
671  rd_hydro_str, & ! [IN]
672  rd_hydro_end, & ! [IN]
673  rd_aero_str, & ! [IN]
674  rd_aero_end, & ! [IN]
675  solins(:,:), & ! [IN]
676  cossza(:,:), & ! [IN]
677  rhodz_merge(:,:,:), & ! [IN]
678  pres_merge(:,:,:), & ! [IN]
679  temp_merge(:,:,:), & ! [IN]
680  temph_merge(:,:,:), & ! [IN]
681  temp_sfc(:,:), & ! [IN]
682  gas_merge(:,:,:,:), & ! [IN]
683  cfc_merge(:,:,:,:), & ! [IN]
684  aerosol_conc_merge(:,:,:,:), & ! [IN]
685  aerosol_radi_merge(:,:,:,:), & ! [IN]
686  i_mpae2rd(:), & ! [IN]
687  cldfrac_merge(:,:,:), & ! [IN]
688  albedo_land(:,:,:), & ! [IN]
689  fact_ocean(:,:), & ! [IN]
690  fact_land(:,:), & ! [IN]
691  fact_urban(:,:), & ! [IN]
692  flux_rad_merge(:,:,:,:,:,:), & ! [OUT]
693  flux_rad_sfc_dn(:,:,:,:) ) ! [OUT]
694 
695  call prof_rapend ('RD_MSTRN_DTRN3', 3)
696 
697  ! return to grid coordinate of model domain
698  do ic = 1, 2
699  do j = js, je
700  do i = is, ie
701  do rd_k = rd_kadd+1, rd_kmax+1
702  k = ks + rd_kmax - rd_k ! reverse axis
703 
704  flux_rad(k,i,j,i_lw,i_up,ic) = flux_rad_merge(rd_k,i,j,i_lw,i_up,ic)
705  flux_rad(k,i,j,i_lw,i_dn,ic) = flux_rad_merge(rd_k,i,j,i_lw,i_dn,ic)
706  flux_rad(k,i,j,i_sw,i_up,ic) = flux_rad_merge(rd_k,i,j,i_sw,i_up,ic)
707  flux_rad(k,i,j,i_sw,i_dn,ic) = flux_rad_merge(rd_k,i,j,i_sw,i_dn,ic)
708  enddo
709  enddo
710  enddo
711  enddo
712 
713 !OCL XFILL
714  do ic = 1, 2
715  do j = js, je
716  do i = is, ie
717  flux_rad_top(i,j,i_lw,i_up,ic) = flux_rad_merge(1,i,j,i_lw,i_up,ic)
718  flux_rad_top(i,j,i_lw,i_dn,ic) = flux_rad_merge(1,i,j,i_lw,i_dn,ic)
719  flux_rad_top(i,j,i_sw,i_up,ic) = flux_rad_merge(1,i,j,i_sw,i_up,ic)
720  flux_rad_top(i,j,i_sw,i_dn,ic) = flux_rad_merge(1,i,j,i_sw,i_dn,ic)
721  enddo
722  enddo
723  enddo
724 
725  return
integer, public is
start point of inner domain: x, local
procedure(cf), pointer, public atmos_phy_mp_cloudfraction
integer, public je
end point of inner domain: y, local
real(rp), parameter, public const_ppm
parts par million
Definition: scale_const.F90:96
module ATMOSPHERE / Saturation adjustment
module ATMOSPHERE / Physics Aerosol Microphysics
real(rp), public real_basepoint_lat
position of base point in real world [rad,-pi,pi]
module ATMOSPHERE / Physics Cloud Microphysics
integer, public ke
end point of inner domain: z, local
integer, public qa
module ATMOSPHERE / Physics Radiation / Vertical profile
integer, parameter, public i_lw
module grid index
integer, parameter, public i_sw
real(rp), public const_mvap
mass weight (water vapor) [g/mol]
Definition: scale_const.F90:64
module TRACER
integer, public ia
of x whole cells (local, with HALO)
module GRID (real space)
integer, public mp_qa
integer, parameter, public i_dn
integer, public ka
of z whole cells (local, with HALO)
integer, public i_qv
procedure(mr), pointer, public atmos_phy_mp_mixingratio
integer, public ae_qa
integer, dimension(:), allocatable, public i_mp2all
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
procedure(er), pointer, public atmos_phy_mp_effectiveradius
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public i_qi
integer, public ie
end point of inner domain: x, local
real(rp), public const_eps
small number
Definition: scale_const.F90:36
module ATMOSPHERE / Thermodynamics
real(rp), dimension(:), pointer, public ae_dens
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.
procedure(er), pointer, public atmos_phy_ae_effectiveradius
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
integer, dimension(:), allocatable, public i_ae2all
real(rp), public const_mdry
mass weight (dry air) [g/mol]
Definition: scale_const.F90:56
real(rp), dimension(:), pointer, public atmos_phy_mp_dens
logical, public atmos_phy_rd_profile_use_climatology
use climatology?
integer, parameter, public i_up
integer, public i_qc
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function: