SCALE-RM
Functions/Subroutines
scale_atmos_sfc_ch_rn222 Module Reference

module atmosphere / surface / chemistry / RN222 More...

Functions/Subroutines

subroutine, public atmos_sfc_ch_rn222_setup (IA, JA, real_lon, real_lat)
 Setup. More...
 
subroutine, public atmos_sfc_ch_rn222_ocean_flux (IA, IS, IE, JA, JS, JE, QA_CH, SFLX_QTRC)
 Emission from the ocean surface. More...
 
subroutine, public atmos_sfc_ch_rn222_land_flux (IA, IS, IE, JA, JS, JE, QA_CH, NOWDATE, SFLX_QTRC)
 Emission from the land surface. More...
 

Detailed Description

module atmosphere / surface / chemistry / RN222

Description
Surface emission component for rn222 tracer
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_SFC_CH_RN222
    nametypedefault valuecomment
    ATMOS_SFC_CH_RN222_EMISSION_TYPE character(len=H_SHORT) 'CONST' Emission type
    ATMOS_SFC_CH_RN222_CONST_EMISSION_LAND real(RP) 20.8E-3_RP Surface flux from land [Bq/m2/s]
    ATMOS_SFC_CH_RN222_CONST_EMISSION_OCEAN real(RP) 0.14E-3_RP Surface flux from ocean [Bq/m2/s]
    ATMOS_SFC_CH_RN222_SCHERY1998_DIRPATH character(len=H_LONG) '.'
    ATMOS_SFC_CH_RN222_HIRAO2010_DIRPATH character(len=H_LONG) '.'
    ATMOS_SFC_CH_RN222_HIRAO2010_YSTART integer 1979
    ATMOS_SFC_CH_RN222_HIRAO2010_YEND integer 2012
    ATMOS_SFC_CH_RN222_NINTRP integer 5

History Output
No history output

Function/Subroutine Documentation

◆ atmos_sfc_ch_rn222_setup()

subroutine, public scale_atmos_sfc_ch_rn222::atmos_sfc_ch_rn222_setup ( integer, intent(in)  IA,
integer, intent(in)  JA,
real(rp), dimension(ia,ja), intent(in)  real_lon,
real(rp), dimension(ia,ja), intent(in)  real_lat 
)

Setup.

Definition at line 81 of file scale_atmos_sfc_ch_rn222.F90.

References scale_const::const_d2r, scale_interp::interp_factor2d(), scale_io::io_fid_conf, scale_io::io_get_available_fid(), scale_prc::prc_abort(), and scale_prc::prc_ismaster.

Referenced by mod_atmos_phy_ch_driver::atmos_phy_ch_driver_setup().

81  use scale_prc, only: &
82  prc_ismaster, &
83  prc_abort
84  use scale_const, only: &
85  const_d2r
86  use scale_comm_cartesc, only: &
87  comm_bcast
88  use scale_interp, only: &
90  implicit none
91 
92  integer, intent(in) :: ia
93  integer, intent(in) :: ja
94  real(RP), intent(in) :: real_lon(ia,ja) ! longitude [rad]
95  real(RP), intent(in) :: real_lat(ia,ja) ! latitude [rad]
96 
97  namelist / param_atmos_sfc_ch_rn222 / &
98  atmos_sfc_ch_rn222_emission_type, &
99  atmos_sfc_ch_rn222_const_emission_land, &
100  atmos_sfc_ch_rn222_const_emission_ocean, &
101  atmos_sfc_ch_rn222_schery1998_dirpath, &
102  atmos_sfc_ch_rn222_hirao2010_dirpath, &
103  atmos_sfc_ch_rn222_hirao2010_ystart, &
104  atmos_sfc_ch_rn222_hirao2010_yend, &
105  atmos_sfc_ch_rn222_nintrp
106 
107  character(len=H_LONG) :: fname
108  real(RP) :: lon, lat, value
109 
110  integer :: ierr, fid
111  integer :: i, j, m, y, yy
112  !---------------------------------------------------------------------------
113 
114  log_newline
115  log_info("ATMOS_SFC_CH_rn222_setup",*) 'Setup'
116  log_info("ATMOS_SFC_CH_rn222_setup",*) 'rn222 surface flux'
117 
118  !--- read namelist
119  rewind(io_fid_conf)
120  read(io_fid_conf,nml=param_atmos_sfc_ch_rn222,iostat=ierr)
121  if( ierr < 0 ) then !--- missing
122  log_info("ATMOS_SFC_CH_rn222_setup",*) 'Not found namelist. Default used.'
123  elseif( ierr > 0 ) then !--- fatal error
124  log_error("ATMOS_SFC_CH_rn222_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_SFC_CH_RN222. Check!'
125  call prc_abort
126  endif
127  log_nml(param_atmos_sfc_ch_rn222)
128 
129 
130 
131  log_newline
132  log_info("ATMOS_SFC_CH_rn222_setup",*) 'Type of emission of Rn222: ', trim(atmos_sfc_ch_rn222_emission_type)
133 
134  select case( atmos_sfc_ch_rn222_emission_type )
135  case( 'CONST' )
136 
137  log_info_cont('(A,ES16.6)') 'From land [Bq/m2/s] : ', atmos_sfc_ch_rn222_const_emission_land
138  log_info_cont('(A,ES16.6)') 'From ocean [Bq/m2/s] : ', atmos_sfc_ch_rn222_const_emission_ocean
139 
140  case( 'SCHERY1998' )
141 
142  log_info_cont(*) 'Flux map by Schery and Wasiolek (1998) is used'
143 
144  nlon = 360
145  nlat = 180
146  nmonth = 12
147  nyear = 1
148 
149  allocate( emission_lon(nlon,nlat) )
150  allocate( emission_lat(nlon,nlat) )
151  allocate( emission_value(nlon,nlat,nmonth,nyear) )
152 
153  if ( prc_ismaster ) then
154  y = 1
155  do m = 1, nmonth
156  write(fname,'(A,A,I2.2,A)') trim(atmos_sfc_ch_rn222_schery1998_dirpath), "/fdh3a.", m
157  log_info_cont(*) 'Read from the ASCII file: ', trim(fname)
158 
159  fid = io_get_available_fid()
160  open( unit = fid, &
161  file = trim(fname), &
162  status = "old", &
163  form = "formatted" )
164 
165  do j = 1, nlat
166  do i = 1, nlon
167  lon = real(i-1,kind=RP) - 180.0_rp ! [180W-179E]
168  lat = 90.0_rp - real(j-1,kind=RP) ! [90N-89S]
169 
170  read(fid,*) value
171 
172  emission_lon(i,j) = ( lon + 0.5_rp ) * const_d2r ! shift +0.5deg, [deg->rad]
173  emission_lat(i,j) = ( lat - 0.5_rp ) * const_d2r ! shift -0.5deg, [deg->rad]
174  emission_value(i,j,m,y) = value * 1.e-3_rp ! [mBq/m2->Bq/m2]
175  enddo
176  enddo
177 
178  close(fid)
179  enddo ! month loop
180  endif
181 
182  case( 'HIRAO2010' )
183 
184  log_info_cont(*) 'Flux map by Hirao et al. (2010) is used'
185  log_info_cont(*) 'Start year: ', atmos_sfc_ch_rn222_hirao2010_ystart
186  log_info_cont(*) 'End year: ', atmos_sfc_ch_rn222_hirao2010_yend
187 
188  if ( atmos_sfc_ch_rn222_hirao2010_ystart < 1979 &
189  .OR. atmos_sfc_ch_rn222_hirao2010_yend > 2012 ) then
190  log_warn('ATMOS_SFC_CH_rn222_setup',*) 'Available period of the data is between 1979 and 2012.'
191  log_warn_cont(*) 'Please check the range of ystart and yend.'
192  atmos_sfc_ch_rn222_hirao2010_ystart = max( atmos_sfc_ch_rn222_hirao2010_ystart, 1979 )
193  atmos_sfc_ch_rn222_hirao2010_yend = min( atmos_sfc_ch_rn222_hirao2010_yend, 2012 )
194  endif
195 
196  nlon = 360
197  nlat = 180
198  nmonth = 12
199  nyear = atmos_sfc_ch_rn222_hirao2010_yend - atmos_sfc_ch_rn222_hirao2010_ystart + 1
200 
201  allocate( emission_lon(nlon,nlat) )
202  allocate( emission_lat(nlon,nlat) )
203  allocate( emission_value(nlon,nlat,nmonth,nyear) )
204 
205  if ( prc_ismaster ) then
206  do y = 1, nyear
207  do m = 1, nmonth
208  yy = y+atmos_sfc_ch_rn222_hirao2010_ystart-1
209  write(fname,'(A,A,I4.4,I2.2)') trim(atmos_sfc_ch_rn222_hirao2010_dirpath), "/flux-hra-revi", yy, m
210  log_info_cont(*) 'Read from the ASCII file: ', trim(fname)
211 
212  fid = io_get_available_fid()
213  open( unit = fid, &
214  file = trim(fname), &
215  status = "old", &
216  form = "formatted" )
217 
218  do j = 1, nlat
219  do i = 1, nlon
220  read(fid,*) lon, lat, value
221 
222  emission_lon(i,j) = ( lon + 0.5_rp ) * const_d2r ! shift +0.5deg, [deg->rad]
223  emission_lat(i,j) = ( lat - 0.5_rp ) * const_d2r ! shift -0.5deg, [deg->rad]
224  emission_value(i,j,m,y) = value * 1.e-3_rp ! [mBq/m2->Bq/m2]
225  enddo
226  enddo
227 
228  close(fid)
229  enddo ! month loop
230  enddo ! year loop
231  endif
232 
233  case default
234  log_error("ATMOS_SFC_CH_rn222_setup",*) 'Not supported type of Rn222 emission! Stop.'
235  call prc_abort
236  end select
237 
238  select case( atmos_sfc_ch_rn222_emission_type )
239  case( 'SCHERY1998', 'HIRAO2010' )
240 
241  call comm_bcast( emission_lon(:,:), nlon, nlat )
242  call comm_bcast( emission_lat(:,:), nlon, nlat )
243  call comm_bcast( emission_value(:,:,:,:), nlon, nlat, nmonth, nyear )
244 
245  allocate( idx_i(ia,ja,atmos_sfc_ch_rn222_nintrp) )
246  allocate( idx_j(ia,ja,atmos_sfc_ch_rn222_nintrp) )
247  allocate( hfact(ia,ja,atmos_sfc_ch_rn222_nintrp) )
248 
249  call interp_factor2d( atmos_sfc_ch_rn222_nintrp, & ! [IN]
250  nlon, nlat, & ! [IN]
251  emission_lon(:,:), & ! [IN]
252  emission_lat(:,:), & ! [IN]
253  ia, ja, & ! [IN]
254  real_lon(:,:), & ! [IN]
255  real_lat(:,:), & ! [IN]
256  idx_i(:,:,:), & ! [OUT]
257  idx_j(:,:,:), & ! [OUT]
258  hfact(:,:,:) ) ! [OUT]
259  end select
260 
261  return
integer, public ia
of whole cells: x, local, with HALO
module INTERPOLATION
integer, public ja
of whole cells: y, local, with HALO
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
module COMMUNICATION
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:313
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:90
subroutine, public interp_factor2d(npoints, IA_ref, JA_ref, lon_ref, lat_ref, IA, JA, lon, lat, idx_i, idx_j, hfact, search_limit, latlon_structure, weight_order)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_sfc_ch_rn222_ocean_flux()

subroutine, public scale_atmos_sfc_ch_rn222::atmos_sfc_ch_rn222_ocean_flux ( integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
integer, intent(in)  QA_CH,
real(rp), dimension(ia,ja,qa_ch), intent(inout)  SFLX_QTRC 
)

Emission from the ocean surface.

Definition at line 271 of file scale_atmos_sfc_ch_rn222.F90.

Referenced by mod_atmos_phy_ch_driver::atmos_phy_ch_driver_ocean_flux().

271  implicit none
272 
273  integer, intent(in) :: ia, is, ie
274  integer, intent(in) :: ja, js, je
275  integer, intent(in) :: qa_ch
276  real(RP), intent(inout) :: sflx_qtrc(ia,ja,qa_ch)
277 
278  integer :: i, j
279  !---------------------------------------------------------------------------
280 
281  select case( atmos_sfc_ch_rn222_emission_type )
282  case( 'CONST', 'SCHERY1998', 'HIRAO2010' )
283 
284  do j = js, je
285  do i = is, ie
286  sflx_qtrc(i,j,i_ch_rn222) = atmos_sfc_ch_rn222_const_emission_ocean
287  enddo
288  enddo
289 
290  end select
291 
292  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public js
start point of inner domain: y, local
Here is the caller graph for this function:

◆ atmos_sfc_ch_rn222_land_flux()

subroutine, public scale_atmos_sfc_ch_rn222::atmos_sfc_ch_rn222_land_flux ( integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
integer, intent(in)  QA_CH,
integer, dimension(6), intent(in)  NOWDATE,
real(rp), dimension(ia,ja,qa_ch), intent(inout)  SFLX_QTRC 
)

Emission from the land surface.

Definition at line 303 of file scale_atmos_sfc_ch_rn222.F90.

References scale_interp::interp_interp2d(), and scale_prc::prc_abort().

Referenced by mod_atmos_phy_ch_driver::atmos_phy_ch_driver_land_flux(), and mod_atmos_phy_ch_driver::atmos_phy_ch_driver_urban_flux().

303  use scale_prc, only: &
304  prc_abort
305  use scale_interp, only: &
307  implicit none
308 
309  integer, intent(in) :: ia, is, ie
310  integer, intent(in) :: ja, js, je
311  integer, intent(in) :: qa_ch
312  integer, intent(in) :: nowdate(6)
313  real(RP), intent(inout) :: sflx_qtrc(ia,ja,qa_ch)
314 
315  integer :: i, j, m, y, yy
316  !---------------------------------------------------------------------------
317 
318  select case( atmos_sfc_ch_rn222_emission_type )
319  case( 'CONST' )
320 
321  do j = js, je
322  do i = is, ie
323  sflx_qtrc(i,j,i_ch_rn222) = atmos_sfc_ch_rn222_const_emission_land
324  enddo
325  enddo
326 
327  case( 'SCHERY1998' )
328 
329  y = 1
330  m = nowdate(2)
331 
332  call interp_interp2d( atmos_sfc_ch_rn222_nintrp, & ! [IN]
333  nlon, nlat, & ! [IN]
334  ia, ja, & ! [IN]
335  idx_i(:,:,:), & ! [IN]
336  idx_j(:,:,:), & ! [IN]
337  hfact(:,:,:), & ! [IN]
338  emission_value(:,:,m,y), & ! [IN]
339  sflx_qtrc(:,:,i_ch_rn222) ) ! [OUT]
340 
341  case( 'HIRAO2010' )
342 
343  yy = nowdate(1)
344  yy = max( yy, 1979 ) ! Use flux of 1979 before 1977
345  yy = min( yy, 2012 ) ! Use flux of 2012 after 2013
346 
347  y = yy-atmos_sfc_ch_rn222_hirao2010_ystart+1
348 
349  if ( y < 1 .OR. y > nyear ) then
350  log_error("ATMOS_SFC_CH_rn222_setup",*) 'emission file does not exist for year=', yy
351  call prc_abort
352  endif
353 
354  m = nowdate(2)
355 
356  call interp_interp2d( atmos_sfc_ch_rn222_nintrp, & ! [IN]
357  nlon, nlat, & ! [IN]
358  ia, ja, & ! [IN]
359  idx_i(:,:,:), & ! [IN]
360  idx_j(:,:,:), & ! [IN]
361  hfact(:,:,:), & ! [IN]
362  emission_value(:,:,m,y), & ! [IN]
363  sflx_qtrc(:,:,i_ch_rn222) ) ! [OUT]
364 
365  end select
366 
367  return
integer, public ia
of whole cells: x, local, with HALO
module INTERPOLATION
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
subroutine, public interp_interp2d(npoints, IA_ref, JA_ref, IA, JA, idx_i, idx_j, hfact, val_ref, val)
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public js
start point of inner domain: y, local
Here is the call graph for this function:
Here is the caller graph for this function: