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 220 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().

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

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