SCALE-RM
scale_atmos_phy_rd_offline.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
12 #include "inc_openmp.h"
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_stdio
20  use scale_prof
22  use scale_tracer
23  !-----------------------------------------------------------------------------
24  implicit none
25  private
26  !-----------------------------------------------------------------------------
27  !
28  !++ Public procedure
29  !
31  public :: atmos_phy_rd_offline
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private procedure
40  !
41  !-----------------------------------------------------------------------------
42  !
43  !++ Private parameters & variables
44  !
45  !-----------------------------------------------------------------------------
46  integer, private, parameter :: num_vars_3d = 4
47  integer, private, parameter :: num_vars_2d = 4
48  integer, private, parameter :: num_vars_2d_op = 1 ! optional
49 
50  real, private :: ATMOS_PHY_RD_offline_diffuse_rate = 0.5_rp
51 
52  logical, private :: vars_2d_exist(num_vars_2d_op)
53 
54 contains
55  !-----------------------------------------------------------------------------
57  subroutine atmos_phy_rd_offline_setup( RD_TYPE )
58  use scale_process, only: &
60  use scale_external_input, only: &
62  use scale_const, only: &
63  undef => const_undef
64  implicit none
65 
66  character(len=*), intent(in) :: rd_type
67 
68  character(len=H_SHORT) :: vars_3d (num_vars_3d)
69  character(len=H_SHORT) :: vars_2d (num_vars_2d)
70  character(len=H_SHORT) :: vars_2d_op(num_vars_2d_op)
71 
72  data vars_3d / 'RFLX_LW_up', 'RFLX_LW_dn', 'RFLX_SW_up', 'RFLX_SW_dn' /
73  data vars_2d / 'SFLX_LW_up', 'SFLX_LW_dn', 'SFLX_SW_up', 'SFLX_SW_dn' /
74  data vars_2d_op / 'SFLX_SW_dn_dir' /
75 
76  character(len=H_LONG) :: atmos_phy_rd_offline_basename = ''
77  character(len=H_SHORT) :: atmos_phy_rd_offline_axistype = 'XYZ'
78  logical :: atmos_phy_rd_offline_enable_periodic_year = .false.
79  logical :: atmos_phy_rd_offline_enable_periodic_month = .false.
80  logical :: atmos_phy_rd_offline_enable_periodic_day = .false.
81  integer :: atmos_phy_rd_offline_step_fixed = 0
82  real(RP) :: atmos_phy_rd_offline_offset = 0.0_rp
83  real(RP) :: atmos_phy_rd_offline_defval
84  logical :: atmos_phy_rd_offline_check_coordinates = .true.
85  integer :: atmos_phy_rd_offline_step_limit = 0
86 
87  namelist / param_atmos_phy_rd_offline / &
88  atmos_phy_rd_offline_basename, &
89  atmos_phy_rd_offline_axistype, &
90  atmos_phy_rd_offline_enable_periodic_year, &
91  atmos_phy_rd_offline_enable_periodic_month, &
92  atmos_phy_rd_offline_enable_periodic_day, &
93  atmos_phy_rd_offline_step_fixed, &
94  atmos_phy_rd_offline_offset, &
95  atmos_phy_rd_offline_defval, &
96  atmos_phy_rd_offline_check_coordinates, &
97  atmos_phy_rd_offline_step_limit, &
98  atmos_phy_rd_offline_diffuse_rate
99 
100  integer :: n, ierr
101  !---------------------------------------------------------------------------
102 
103  if( io_l ) write(io_fid_log,*)
104  if( io_l ) write(io_fid_log,*) '++++++ Module[RADIATION] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
105  if( io_l ) write(io_fid_log,*) '*** Offline radiation process'
106 
107  if ( rd_type /= 'OFFLINE' ) then
108  write(*,*) 'xxx RD_TYPE is not OFFLINE. Check!'
109  call prc_mpistop
110  endif
111 
112 
113  atmos_phy_rd_offline_defval = undef
114 
115  !--- read namelist
116  rewind(io_fid_conf)
117  read(io_fid_conf,nml=param_atmos_phy_rd_offline,iostat=ierr)
118  if( ierr < 0 ) then !--- missing
119  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
120  elseif( ierr > 0 ) then !--- fatal error
121  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_RD_OFFLINE. Check!'
122  call prc_mpistop
123  endif
124  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_rd_offline)
125 
126  if ( atmos_phy_rd_offline_basename == '' ) then
127  write(*,*) 'xxx ATMOS_PHY_RD_offline_basename is necessary'
128  call prc_mpistop
129  end if
130 
131  do n = 1, num_vars_3d
132  call extin_regist( atmos_phy_rd_offline_basename, & ! [IN]
133  vars_3d(n), & ! [IN]
134  atmos_phy_rd_offline_axistype, & ! [IN]
135  atmos_phy_rd_offline_enable_periodic_year, & ! [IN]
136  atmos_phy_rd_offline_enable_periodic_month, & ! [IN]
137  atmos_phy_rd_offline_enable_periodic_day, & ! [IN]
138  atmos_phy_rd_offline_step_fixed, & ! [IN]
139  atmos_phy_rd_offline_offset, & ! [IN]
140  atmos_phy_rd_offline_defval, & ! [IN]
141  atmos_phy_rd_offline_check_coordinates, & ! [IN]
142  atmos_phy_rd_offline_step_limit ) ! [IN]
143  end do
144 
145  do n = 1, num_vars_2d
146  call extin_regist( atmos_phy_rd_offline_basename, & ! [IN]
147  vars_2d(n), & ! [IN]
148  'XY', & ! [IN]
149  atmos_phy_rd_offline_enable_periodic_year, & ! [IN]
150  atmos_phy_rd_offline_enable_periodic_month, & ! [IN]
151  atmos_phy_rd_offline_enable_periodic_day, & ! [IN]
152  atmos_phy_rd_offline_step_fixed, & ! [IN]
153  atmos_phy_rd_offline_offset, & ! [IN]
154  atmos_phy_rd_offline_defval, & ! [IN]
155  atmos_phy_rd_offline_check_coordinates, & ! [IN]
156  atmos_phy_rd_offline_step_limit ) ! [IN]
157  end do
158 
159  do n = 1, num_vars_2d_op
160  call extin_regist( atmos_phy_rd_offline_basename, & ! [IN]
161  vars_2d_op(n), & ! [IN]
162  'XY', & ! [IN]
163  atmos_phy_rd_offline_enable_periodic_year, & ! [IN]
164  atmos_phy_rd_offline_enable_periodic_month, & ! [IN]
165  atmos_phy_rd_offline_enable_periodic_day, & ! [IN]
166  atmos_phy_rd_offline_step_fixed, & ! [IN]
167  atmos_phy_rd_offline_offset, & ! [IN]
168  atmos_phy_rd_offline_defval, & ! [IN]
169  atmos_phy_rd_offline_check_coordinates, & ! [IN]
170  atmos_phy_rd_offline_step_limit, & ! [IN]
171  exist = vars_2d_exist(n) ) ! [OUT]
172  if ( vars_2d_exist(n) ) then
173  if( io_l ) write(io_fid_log,*) '*** ', trim(vars_2d_op(n)), ' found.'
174  else
175  if( io_l ) write(io_fid_log,*) '*** ', trim(vars_2d_op(n)), ' not found.'
176  end if
177  end do
178 
179  return
180  end subroutine atmos_phy_rd_offline_setup
181 
182  !-----------------------------------------------------------------------------
184  subroutine atmos_phy_rd_offline( &
185  DENS, RHOT, QTRC, &
186  CZ, FZ, &
187  fact_ocean, &
188  fact_land, &
189  fact_urban, &
190  temp_sfc, albedo_land, &
191  solins, cosSZA, &
192  flux_rad, &
193  flux_rad_top, &
194  SFLX_rad_dn )
195 ! Jval )
196  use scale_grid_index
197  use scale_tracer
198  use scale_process, only: &
200  use scale_external_input, only: &
201  extin_update
202  use scale_time, only: &
204  use scale_atmos_phy_rd_common, only: &
205  i_sw, &
206  i_lw, &
207  i_dn, &
208  i_up, &
209  i_direct, &
210  i_diffuse
211  implicit none
212  real(RP), intent(in) :: dens (ka,ia,ja)
213  real(RP), intent(in) :: rhot (ka,ia,ja)
214  real(RP), intent(in) :: qtrc (ka,ia,ja,qa)
215  real(RP), intent(in) :: cz ( ka,ia,ja) ! UNUSED
216  real(RP), intent(in) :: fz (0:ka,ia,ja)
217  real(RP), intent(in) :: fact_ocean (ia,ja)
218  real(RP), intent(in) :: fact_land (ia,ja)
219  real(RP), intent(in) :: fact_urban (ia,ja)
220  real(RP), intent(in) :: temp_sfc (ia,ja)
221  real(RP), intent(in) :: albedo_land (ia,ja,2)
222  real(RP), intent(in) :: solins (ia,ja)
223  real(RP), intent(in) :: cossza (ia,ja)
224  real(RP), intent(out) :: flux_rad (ka,ia,ja,2,2,2)
225  real(RP), intent(out) :: flux_rad_top(ia,ja,2,2,2)
226  real(RP), intent(out) :: sflx_rad_dn (ia,ja,2,2)
227 ! real(RP), intent(out) :: Jval (KA,IA,JA,CH_QA_photo)
228 
229  real(RP) :: buffer(ia,ja)
230  logical :: error, error_sum
231 
232  integer :: i, j
233  !---------------------------------------------------------------------------
234 
235  if( io_l ) write(io_fid_log,*) '*** Atmos physics step: Radiation(offline)'
236 
237  ! [note] external data input is now support only SCALE-RM
238 
239  error_sum = .false.
240 
241  ! 3D
242  call extin_update( flux_rad(:,:,:,i_lw,i_up,2), 'RFLX_LW_up', time_nowdaysec, error )
243  error_sum = ( error .OR. error_sum )
244 
245  call extin_update( flux_rad(:,:,:,i_lw,i_dn,2), 'RFLX_LW_dn', time_nowdaysec, error )
246  error_sum = ( error .OR. error_sum )
247 
248  call extin_update( flux_rad(:,:,:,i_sw,i_up,2), 'RFLX_SW_up', time_nowdaysec, error )
249  error_sum = ( error .OR. error_sum )
250 
251  call extin_update( flux_rad(:,:,:,i_sw,i_dn,2), 'RFLX_SW_dn', time_nowdaysec, error )
252  error_sum = ( error .OR. error_sum )
253 
254 
255  ! 2D
256  call extin_update( buffer(:,:), 'SFLX_LW_up', time_nowdaysec, error )
257  if ( error ) then
258  error_sum = .true.
259  else
260  !$omp parallel do default(none) OMP_SCHEDULE_ &
261  !$omp private(i,j) &
262  !$omp shared(KS,IS,IE,JS,JE) &
263  !$omp shared(flux_rad,buffer)
264  do j = js, je
265  do i = is, ie
266  flux_rad(ks-1,i,j,i_lw,i_up,2) = buffer(i,j)
267  end do
268  end do
269  end if
270 
271  call extin_update( buffer(:,:), 'SFLX_LW_dn', time_nowdaysec, error )
272  if ( error ) then
273  error_sum = .true.
274  else
275  !$omp parallel do default(none) OMP_SCHEDULE_ &
276  !$omp private(i,j) &
277  !$omp shared(KS,IS,IE,JS,JE) &
278  !$omp shared(flux_rad,buffer)
279  do j = js, je
280  do i = is, ie
281  flux_rad(ks-1,i,j,i_lw,i_dn,2) = buffer(i,j)
282  end do
283  end do
284  end if
285 
286  call extin_update( buffer(:,:), 'SFLX_SW_up', time_nowdaysec, error )
287  if ( error ) then
288  error_sum = .true.
289  else
290  !$omp parallel do default(none) OMP_SCHEDULE_ &
291  !$omp private(i,j) &
292  !$omp shared(KS,IS,IE,JS,JE) &
293  !$omp shared(flux_rad,buffer)
294  do j = js, je
295  do i = is, ie
296  flux_rad(ks-1,i,j,i_sw,i_up,2) = buffer(i,j)
297  end do
298  end do
299  end if
300 
301  call extin_update( buffer(:,:), 'SFLX_SW_dn', time_nowdaysec, error )
302  if ( error ) then
303  error_sum = .true.
304  else
305  !$omp parallel do default(none) OMP_SCHEDULE_ &
306  !$omp private(i,j) &
307  !$omp shared(KS,IS,IE,JS,JE) &
308  !$omp shared(flux_rad,buffer)
309  do j = js, je
310  do i = is, ie
311  flux_rad(ks-1,i,j,i_sw,i_dn,2) = buffer(i,j)
312  end do
313  end do
314  end if
315 
316  !$omp parallel do default(none) OMP_SCHEDULE_ &
317  !$omp private(i,j) &
318  !$omp shared(KS,IS,IE,JS,JE,ATMOS_PHY_RD_offline_diffuse_rate) &
319  !$omp shared(SFLX_rad_dn,flux_rad)
320  do j = js, je
321  do i = is, ie
322  sflx_rad_dn(i,j,i_lw,i_diffuse) = flux_rad(ks-1,i,j,i_lw,i_dn,2)
323  sflx_rad_dn(i,j,i_lw,i_direct ) = 0.0_rp
324  end do
325  end do
326 
327  ! 2D optional
328  if ( vars_2d_exist(1) ) then
329  call extin_update( sflx_rad_dn(:,:,i_sw,i_direct), 'SFLX_SW_dn_dir', time_nowdaysec, error )
330  if ( error ) then
331  error_sum = .true.
332  else
333  !$omp parallel do default(none) OMP_SCHEDULE_ &
334  !$omp private(i,j) &
335  !$omp shared(KS,IS,IE,JS,JE) &
336  !$omp shared(SFLX_rad_dn,flux_rad)
337  do j = js, je
338  do i = is, ie
339  sflx_rad_dn(i,j,i_sw,i_diffuse) = flux_rad(ks-1,i,j,i_sw,i_dn,2) - sflx_rad_dn(i,j,i_sw,i_direct)
340  end do
341  end do
342  end if
343  else
344  !$omp parallel do default(none) OMP_SCHEDULE_ &
345  !$omp private(i,j) &
346  !$omp shared(KS,IS,IE,JS,JE,ATMOS_PHY_RD_offline_diffuse_rate) &
347  !$omp shared(SFLX_rad_dn,flux_rad)
348  do j = js, je
349  do i = is, ie
350  sflx_rad_dn(i,j,i_sw,i_diffuse) = ( atmos_phy_rd_offline_diffuse_rate ) * flux_rad(ks-1,i,j,i_sw,i_dn,2)
351  sflx_rad_dn(i,j,i_sw,i_direct ) = ( 1.0_rp - atmos_phy_rd_offline_diffuse_rate ) * flux_rad(ks-1,i,j,i_sw,i_dn,2)
352  end do
353  end do
354  end if
355 
356  if ( error_sum ) then
357  write(*,*) 'xxx Requested data is not found!'
358  call prc_mpistop
359  endif
360 
361  ! clearsky and TOA value are not defined
362  flux_rad(:,:,:,:,:,1) = 0.0_rp
363  flux_rad_top(:,:,:,:,:) = 0.0_rp
364 
365  return
366  end subroutine atmos_phy_rd_offline
367 
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
integer, parameter, public i_direct
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
subroutine, public atmos_phy_rd_offline_setup(RD_TYPE)
Setup.
integer, parameter, public i_diffuse
integer, parameter, public i_lw
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
integer, parameter, public i_sw
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:62
module TRACER
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_dn
integer, public ka
of whole cells: z, local, with HALO
module ATMOSPHERE / Physics Radiation
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module PROCESS
subroutine, public extin_regist(basename, varname, axistype, enable_periodic_year, enable_periodic_month, enable_periodic_day, step_fixed, offset, defval, check_coordinates, step_limit, exist)
Regist data.
module CONSTANT
Definition: scale_const.F90:14
module EXTERNAL INPUT
integer, public ks
start point of inner domain: z, local
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
subroutine, public atmos_phy_rd_offline(DENS, RHOT, QTRC, CZ, FZ, fact_ocean, fact_land, fact_urban, temp_sfc, albedo_land, solins, cosSZA, flux_rad, flux_rad_top, SFLX_rad_dn)
Radiation main.
module ATMOSPHERE / Physics Radiation
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
integer, parameter, public i_up
integer, public io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
integer, public ja
of whole cells: y, local, with HALO