SCALE-RM
mod_atmos_bnd_driver.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
12 #include "scalelib.h"
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_io
20  use scale_prof
22  use scale_index
23  use scale_tracer
24 
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public procedure
31  !
37 
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public parameters & variables
41  !
42  integer, public :: bnd_qa
43  integer, public, allocatable :: bnd_iq(:)
44 
45  real(rp), public, allocatable :: atmos_boundary_dens(:,:,:)
46  real(rp), public, allocatable :: atmos_boundary_velz(:,:,:)
47  real(rp), public, allocatable :: atmos_boundary_velx(:,:,:)
48  real(rp), public, allocatable :: atmos_boundary_vely(:,:,:)
49  real(rp), public, allocatable :: atmos_boundary_pott(:,:,:)
50  real(rp), public, allocatable :: atmos_boundary_qtrc(:,:,:,:)
51 
52  real(rp), public, allocatable :: atmos_boundary_alpha_dens(:,:,:)
53  real(rp), public, allocatable :: atmos_boundary_alpha_velz(:,:,:)
54  real(rp), public, allocatable :: atmos_boundary_alpha_velx(:,:,:)
55  real(rp), public, allocatable :: atmos_boundary_alpha_vely(:,:,:)
56  real(rp), public, allocatable :: atmos_boundary_alpha_pott(:,:,:)
57  real(rp), public, allocatable :: atmos_boundary_alpha_qtrc(:,:,:,:)
58 
59  real(rp), public, allocatable :: atmos_boundary_mflux_offset_x(:,:,:)
60  real(rp), public, allocatable :: atmos_boundary_mflux_offset_y(:,:,:)
61 
62  real(rp), public :: atmos_boundary_smoother_fact = 0.2_rp
63 
64  logical, public :: atmos_boundary_update_flag = .false.
65 
66  !-----------------------------------------------------------------------------
67  !
68  !++ Private procedure
69  !
70  private :: atmos_boundary_var_fillhalo
71  private :: atmos_boundary_alpha_fillhalo
72  private :: atmos_boundary_ref_fillhalo
73  private :: atmos_boundary_setalpha
74  private :: atmos_boundary_setinitval
75  private :: atmos_boundary_read
76  private :: atmos_boundary_write
77  private :: atmos_boundary_generate
78  private :: atmos_boundary_initialize_file
79  private :: atmos_boundary_initialize_online
80  private :: atmos_boundary_update_file
81  private :: atmos_boundary_update_online_parent
82  private :: atmos_boundary_update_online_daughter
83  private :: atmos_boundary_firstsend
84  private :: atmos_boundary_send
85  private :: atmos_boundary_recv
86 
87  abstract interface
88  subroutine getbnd( &
89  bnd_DENS, &
90  bnd_VELZ, &
91  bnd_VELX, &
92  bnd_VELY, &
93  bnd_POTT, &
94  bnd_QTRC, &
95  now_step, &
96  update_step )
97  use scale_precision
98  implicit none
99 
100  real(rp), intent(out) :: bnd_dens(:,:,:)
101  real(rp), intent(out) :: bnd_velz(:,:,:)
102  real(rp), intent(out) :: bnd_velx(:,:,:)
103  real(rp), intent(out) :: bnd_vely(:,:,:)
104  real(rp), intent(out) :: bnd_pott(:,:,:)
105  real(rp), intent(out) :: bnd_qtrc(:,:,:,:)
106  integer, intent(in) :: now_step
107  integer, intent(in) :: update_step
108  end subroutine getbnd
109  end interface
110 
111  procedure(getbnd), pointer :: get_boundary => null()
112  private :: get_boundary
113  private :: get_boundary_same_parent
114  private :: get_boundary_nearest_neighbor
115  private :: get_boundary_lerp_initpoint
116  private :: get_boundary_lerp_midpoint
117 
118  !
119  !-----------------------------------------------------------------------------
120  !
121  !++ Private parameters & variables
122  !
123  character(len=H_SHORT), private :: atmos_boundary_type = 'NONE'
124  character(len=H_LONG), private :: atmos_boundary_in_basename = ''
125  logical, private :: atmos_boundary_in_check_coordinates = .true.
126  logical, private :: atmos_boundary_in_aggregate
127  character(len=H_LONG), private :: atmos_boundary_out_basename = ''
128  character(len=H_MID), private :: atmos_boundary_out_title = 'SCALE-RM BOUNDARY CONDITION'
129  character(len=H_SHORT), private :: atmos_boundary_out_dtype = 'DEFAULT'
130  logical, private :: atmos_boundary_out_aggregate
131 
132  logical, private :: atmos_boundary_use_dens = .false. ! read from file?
133  logical, private :: atmos_boundary_use_velz = .false. ! read from file?
134  logical, private :: atmos_boundary_use_velx = .false. ! read from file?
135  logical, private :: atmos_boundary_use_vely = .false. ! read from file?
136  logical, private :: atmos_boundary_use_pt = .false. ! read from file?
137  logical, private :: atmos_boundary_use_qv = .false. ! read from file?
138  logical, private :: atmos_boundary_use_qhyd = .false. ! read from file?
139  logical, private :: atmos_boundary_use_chem = .false. ! read from file?
140 
141  real(rp), private :: atmos_boundary_value_velz = 0.0_rp ! velocity w at boundary, 0 [m/s]
142  real(rp), private :: atmos_boundary_value_velx = 0.0_rp ! velocity u at boundary, 0 [m/s]
143  real(rp), private :: atmos_boundary_value_vely = 0.0_rp ! velocity v at boundary, 0 [m/s]
144  real(rp), private :: atmos_boundary_value_pt = 300.0_rp ! potential temp. at boundary, 300 [K]
145  real(rp), private :: atmos_boundary_value_qtrc = 0.0_rp ! tracer at boundary, 0 [kg/kg]
146 
147  real(rp), private :: atmos_boundary_alphafact_dens = 1.0_rp ! alpha factor again default
148  real(rp), private :: atmos_boundary_alphafact_velz = 1.0_rp ! alpha factor again default
149  real(rp), private :: atmos_boundary_alphafact_velx = 1.0_rp ! alpha factor again default
150  real(rp), private :: atmos_boundary_alphafact_vely = 1.0_rp ! alpha factor again default
151  real(rp), private :: atmos_boundary_alphafact_pt = 1.0_rp ! alpha factor again default
152  real(rp), private :: atmos_boundary_alphafact_qtrc = 1.0_rp ! alpha factor again default
153 
154  real(rp), private :: atmos_boundary_fracz = 1.0_rp ! fraction of boundary region for dumping (z) (0-1)
155  real(rp), private :: atmos_boundary_fracx = 1.0_rp ! fraction of boundary region for dumping (x) (0-1)
156  real(rp), private :: atmos_boundary_fracy = 1.0_rp ! fraction of boundary region for dumping (y) (0-1)
157  real(rp), private :: atmos_boundary_tauz ! maximum value for damping tau (z) [s]
158  real(rp), private :: atmos_boundary_taux ! maximum value for damping tau (x) [s]
159  real(rp), private :: atmos_boundary_tauy ! maximum value for damping tau (y) [s]
160 
161  real(dp), private :: atmos_boundary_update_dt = 0.0_dp ! inteval time of boudary data update [s]
162  integer, private :: update_nstep
163 
164  logical, private :: atmos_grid_nudging_uv = .false. ! grid nudging
165  logical, private :: atmos_grid_nudging_pt = .false. ! grid nudging
166  logical, private :: atmos_grid_nudging_qv = .false. ! grid nudging
167  real(rp), private :: atmos_grid_nudging_tau ! Damping tau for grid nudging [s]
168 
169  real(rp), private, allocatable :: atmos_boundary_ref_dens(:,:,:,:) ! reference DENS (with HALO)
170  real(rp), private, allocatable :: atmos_boundary_ref_velz(:,:,:,:) ! reference VELZ (with HALO)
171  real(rp), private, allocatable :: atmos_boundary_ref_velx(:,:,:,:) ! reference VELX (with HALO)
172  real(rp), private, allocatable :: atmos_boundary_ref_vely(:,:,:,:) ! reference VELY (with HALO)
173  real(rp), private, allocatable :: atmos_boundary_ref_pott(:,:,:,:) ! reference POTT (with HALO)
174  real(rp), private, allocatable, target :: atmos_boundary_ref_qtrc(:,:,:,:,:) ! reference QTRC (with HALO)
175 
176  ! work
177  real(rp), private, allocatable, target :: q_work(:,:,:,:) ! QV + Qe
178 
179 
180  character(len=H_SHORT), private :: atmos_boundary_interp_type = 'lerp_initpoint' ! type of boundary interporation
181 
182  integer, private :: atmos_boundary_start_date(6) = (/ -9999, 0, 0, 0, 0, 0 /) ! boundary initial date
183 
184  integer, private :: atmos_boundary_fid = -1
185 
186  integer, private :: now_step
187  integer, private :: boundary_timestep = 0
188  logical, private :: atmos_boundary_linear_v = .false. ! linear or non-linear profile of relax region
189  logical, private :: atmos_boundary_linear_h = .true. ! linear or non-linear profile of relax region
190  real(rp), private :: atmos_boundary_exp_h = 2.0_rp ! factor of non-linear profile of relax region
191  logical, private :: atmos_boundary_online = .false. ! boundary online update by communicate inter-domain
192  logical, private :: atmos_boundary_online_master = .false. ! master domain in communicate inter-domain
193 
194  logical, private :: atmos_boundary_dens_adjust = .false.
195  real(rp), private :: atmos_boundary_dens_adjust_tau = -1.0_rp
196 
197  logical, private :: do_parent_process = .false.
198  logical, private :: do_daughter_process = .false.
199  logical, private :: l_bnd = .false.
200 
201  real(dp), private :: boundary_time_initdaysec
202 
203  integer, private :: ref_size = 3
204  integer, private :: ref_old = 1
205  integer, private :: ref_now = 2
206  integer, private :: ref_new = 3
207 
208  ! for mass flux offset
209  real(dp), private :: masstot_now = 0.0_dp
210  real(dp), private :: massflx_now = 0.0_dp
211  real(rp), allocatable, private :: areazuy_w(:,:), areazuy_e(:,:)
212  real(rp), allocatable, private :: offset_time_fact(:)
213  real(rp), allocatable, private :: mflux_offset_x(:,:,:,:)
214  real(rp), allocatable, private :: mflux_offset_y(:,:,:,:)
215  real(rp), allocatable, private, target :: zero_x(:,:), zero_y(:,:)
216 
217 
218  !-----------------------------------------------------------------------------
219 contains
220  !-----------------------------------------------------------------------------
222  subroutine atmos_boundary_driver_setup
223  use scale_prc, only: &
224  prc_abort
225  use scale_const, only: &
226  undef => const_undef
227  use scale_time, only: &
228  dt => time_dtsec
229  use scale_file, only: &
231  use scale_comm_cartesc_nest, only: &
234  use_nesting, &
235  offline, &
238  nestqa => comm_cartesc_nest_bnd_qa
239  use scale_atmos_hydrometeor, only: &
241  i_qv
242  use mod_atmos_phy_mp_vars, only: &
243  qs_mp, &
244  qe_mp
245  use mod_atmos_phy_ch_vars, only: &
246  qs_ch, &
247  qe_ch
248  use scale_atmos_grid_cartesc_real, only: &
250  implicit none
251 
252  namelist / param_atmos_boundary / &
253  atmos_boundary_type, &
254  atmos_boundary_in_basename, &
255  atmos_boundary_in_check_coordinates, &
256  atmos_boundary_in_aggregate, &
257  atmos_boundary_out_basename, &
258  atmos_boundary_out_title, &
259  atmos_boundary_out_dtype, &
260  atmos_boundary_out_aggregate, &
261  atmos_boundary_use_velz, &
262  atmos_boundary_use_velx, &
263  atmos_boundary_use_vely, &
264  atmos_boundary_use_pt, &
265  atmos_boundary_use_dens, &
266  atmos_boundary_use_qv, &
267  atmos_boundary_use_qhyd, &
268  atmos_boundary_use_chem, &
269  atmos_boundary_dens_adjust, &
270  atmos_boundary_dens_adjust_tau, &
271  atmos_boundary_value_velz, &
272  atmos_boundary_value_velx, &
273  atmos_boundary_value_vely, &
274  atmos_boundary_value_pt, &
275  atmos_boundary_value_qtrc, &
276  atmos_boundary_alphafact_dens, &
277  atmos_boundary_alphafact_velz, &
278  atmos_boundary_alphafact_velx, &
279  atmos_boundary_alphafact_vely, &
280  atmos_boundary_alphafact_pt, &
281  atmos_boundary_alphafact_qtrc, &
283  atmos_boundary_fracz, &
284  atmos_boundary_fracx, &
285  atmos_boundary_fracy, &
286  atmos_boundary_tauz, &
287  atmos_boundary_taux, &
288  atmos_boundary_tauy, &
289  atmos_boundary_update_dt, &
290  atmos_boundary_start_date, &
291  atmos_boundary_linear_v, &
292  atmos_boundary_linear_h, &
293  atmos_boundary_exp_h, &
294  atmos_boundary_interp_type, &
295  atmos_grid_nudging_uv, &
296  atmos_grid_nudging_pt, &
297  atmos_grid_nudging_qv, &
298  atmos_grid_nudging_tau
299 
300  integer :: k, i, j, iq
301  integer :: ierr
302  !---------------------------------------------------------------------------
303 
304  log_newline
305  log_info("ATMOS_BOUNDARY_setup",*) 'Setup'
306 
307 
308  atmos_boundary_tauz = dt * 10.0_rp
309  atmos_boundary_taux = dt * 10.0_rp
310  atmos_boundary_tauy = dt * 10.0_rp
311 
312  atmos_boundary_in_aggregate = file_aggregate
313  atmos_boundary_out_aggregate = file_aggregate
314 
315  atmos_grid_nudging_tau = 10.0_rp * 24.0_rp * 3600.0_rp ! 10days [s]
316 
317  !--- read namelist
318  rewind(io_fid_conf)
319  read(io_fid_conf,nml=param_atmos_boundary,iostat=ierr)
320  if( ierr < 0 ) then !--- missing
321  log_info("ATMOS_BOUNDARY_setup",*) 'Not found namelist. Default used.'
322  elseif( ierr > 0 ) then !--- fatal error
323  log_error("ATMOS_BOUNDARY_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_BOUNDARY. Check!'
324  call prc_abort
325  endif
326  log_nml(param_atmos_boundary)
327 
328  ! setting switches
329  if( .NOT. use_nesting ) then
330  atmos_boundary_online = .false.
331  else
332  if( offline ) then
333  atmos_boundary_online = .false.
334  else
335  atmos_boundary_online = .true.
336  endif
337  endif
338  do_parent_process = .false.
339  do_daughter_process = .false.
340  atmos_boundary_online_master = .false.
341  if ( atmos_boundary_online ) then
342  if ( online_iam_parent ) then
343  do_parent_process = .true.
344  if ( .NOT. online_iam_daughter ) then
345  atmos_boundary_online_master = .true.
346  endif
347  endif
348  if ( online_iam_daughter ) then
349  do_daughter_process = .true.
350  atmos_boundary_use_qhyd = online_boundary_use_qhyd
351  endif
352  endif
353 
354  allocate( bnd_iq(qa) )
355  bnd_iq(:) = -1
356  bnd_qa = 0
357  if ( .not. atmos_hydrometeor_dry ) then
358  bnd_qa = bnd_qa + 1
359  bnd_iq(i_qv) = bnd_qa
360  if( atmos_boundary_use_qhyd ) then
361  do iq = qs_mp+1, qe_mp
362  bnd_qa = bnd_qa + 1
363  bnd_iq(iq) = bnd_qa
364  end do
365  end if
366  end if
367  if( atmos_boundary_use_chem ) then
368  do iq = qs_ch, qe_ch
369  bnd_qa = bnd_qa + 1
370  bnd_iq(iq) = bnd_qa
371  end do
372  endif
373 
374  if ( atmos_boundary_dens_adjust_tau <= 0.0_rp ) then
375  atmos_boundary_dens_adjust_tau = max( real(atmos_boundary_update_dt,kind=rp) / 6.0_rp, &
376  atmos_boundary_taux, atmos_boundary_tauy )
377  end if
378 
379 
380 
381  allocate( atmos_boundary_dens(ka,ia,ja) )
382  allocate( atmos_boundary_velz(ka,ia,ja) )
383  allocate( atmos_boundary_velx(ka,ia,ja) )
384  allocate( atmos_boundary_vely(ka,ia,ja) )
385  allocate( atmos_boundary_pott(ka,ia,ja) )
386  allocate( atmos_boundary_qtrc(ka,ia,ja,bnd_qa) )
387  atmos_boundary_dens(:,:,:) = undef
388  atmos_boundary_velz(:,:,:) = undef
389  atmos_boundary_velx(:,:,:) = undef
390  atmos_boundary_vely(:,:,:) = undef
391  atmos_boundary_pott(:,:,:) = undef
392  atmos_boundary_qtrc(:,:,:,:) = undef
393 
394  allocate( atmos_boundary_alpha_dens(ka,ia,ja) )
395  allocate( atmos_boundary_alpha_velz(ka,ia,ja) )
396  allocate( atmos_boundary_alpha_velx(ka,ia,ja) )
397  allocate( atmos_boundary_alpha_vely(ka,ia,ja) )
398  allocate( atmos_boundary_alpha_pott(ka,ia,ja) )
400  atmos_boundary_alpha_dens(:,:,:) = 0.0_rp
401  atmos_boundary_alpha_velz(:,:,:) = 0.0_rp
402  atmos_boundary_alpha_velx(:,:,:) = 0.0_rp
403  atmos_boundary_alpha_vely(:,:,:) = 0.0_rp
404  atmos_boundary_alpha_pott(:,:,:) = 0.0_rp
405  atmos_boundary_alpha_qtrc(:,:,:,:) = 0.0_rp
406 
407  allocate( atmos_boundary_mflux_offset_x(ka,ja,2) )
408  allocate( atmos_boundary_mflux_offset_y(ka,ia,2) )
409  atmos_boundary_mflux_offset_x(:,:,:) = 0.0_rp
410  atmos_boundary_mflux_offset_y(:,:,:) = 0.0_rp
411 
412  if ( atmos_boundary_type == 'REAL' .OR. do_daughter_process ) then
413  l_bnd = .true.
414  else
415  l_bnd = .false.
416  end if
417 
418  if ( l_bnd ) then
419 
420  select case(atmos_boundary_interp_type)
421  case('same_parent')
422  get_boundary => get_boundary_same_parent
423  case('nearest_neighbor')
424  get_boundary => get_boundary_nearest_neighbor
425  case('lerp_initpoint')
426  get_boundary => get_boundary_lerp_initpoint
427  case('lerp_midpoint')
428  get_boundary => get_boundary_lerp_midpoint
429  case default
430  log_error("ATMOS_BOUNDARY_setup",*) 'Wrong parameter in ATMOS_BOUNDARY_interp_TYPE. Check!'
431  call prc_abort
432  end select
433 
434  allocate( atmos_boundary_ref_dens(ka,ia,ja,ref_size) )
435  allocate( atmos_boundary_ref_velz(ka,ia,ja,ref_size) )
436  allocate( atmos_boundary_ref_velx(ka,ia,ja,ref_size) )
437  allocate( atmos_boundary_ref_vely(ka,ia,ja,ref_size) )
438  allocate( atmos_boundary_ref_pott(ka,ia,ja,ref_size) )
439  allocate( atmos_boundary_ref_qtrc(ka,ia,ja,bnd_qa,ref_size) )
440  atmos_boundary_ref_dens(:,:,:,:) = undef
441  atmos_boundary_ref_velz(:,:,:,:) = undef
442  atmos_boundary_ref_velx(:,:,:,:) = undef
443  atmos_boundary_ref_vely(:,:,:,:) = undef
444  atmos_boundary_ref_pott(:,:,:,:) = undef
445  atmos_boundary_ref_qtrc(:,:,:,:,:) = undef
446 
447  ! initialize boundary value (reading file or waiting parent domain)
448  if ( do_daughter_process ) then
449  call atmos_boundary_initialize_online
450  else
451  if ( atmos_boundary_in_basename /= '' ) then
452  call atmos_boundary_initialize_file
453  else
454  log_error("ATMOS_BOUNDARY_setup",*) 'You need specify ATMOS_BOUNDARY_IN_BASENAME'
455  call prc_abort
456  endif
457  endif
458 
459  call atmos_boundary_setalpha
460 
462 
463  ! for mass flux offset
464  allocate( areazuy_w(ka,ja), areazuy_e(ka,ja) )
465  allocate( mflux_offset_x(ka,ja,2,2), mflux_offset_y(ka,ia,2,2) )
466 
467  !$omp parallel do
468  do j = js, je
469  do k = ks, ke
470  areazuy_w(k,j) = atmos_grid_cartesc_real_areazuy_x(k,is-1,j)
471  areazuy_e(k,j) = atmos_grid_cartesc_real_areazuy_x(k,ie ,j)
472  end do
473  end do
474  !$omp parallel do
475  do j = js, je
476  do k = ks, ke
477  mflux_offset_x(k,j,:,:) = 0.0_rp
478  end do
479  end do
480  !$omp parallel do
481  do i = is, ie
482  do k = ks, ke
483  mflux_offset_y(k,i,:,:) = 0.0_rp
484  end do
485  end do
486 
487  elseif ( atmos_boundary_type == 'NONE' ) then
488 
490 
491  elseif ( atmos_boundary_type == 'CONST' ) then
492 
493  call atmos_boundary_setalpha
494 
496 
497  elseif ( atmos_boundary_type == 'INIT' ) then
498 
499  call atmos_boundary_setalpha
500 
502 
503  elseif ( atmos_boundary_type == 'OFFLINE' ) then
504 
505  if ( atmos_boundary_in_basename /= '' ) then
506  call atmos_boundary_read
507  else
508  log_error("ATMOS_BOUNDARY_setup",*) 'You need specify ATMOS_BOUNDARY_IN_BASENAME'
509  call prc_abort
510  endif
511 
513 
514  else
515  log_error("ATMOS_BOUNDARY_setup",*) 'unsupported ATMOS_BOUNDARY_TYPE. Check!', trim(atmos_boundary_type)
516  call prc_abort
517  endif
518 
520 
521 
522  !----- report data -----
523  log_newline
524  log_info("ATMOS_BOUNDARY_setup",*) 'Atmospheric boundary parameters '
525  log_info_cont(*) 'Atmospheric boundary type : ', atmos_boundary_type
526  log_newline
527  log_info_cont(*) 'Is VELZ used in atmospheric boundary? : ', atmos_boundary_use_velz
528  log_info_cont(*) 'Is VELX used in atmospheric boundary? : ', atmos_boundary_use_velx
529  log_info_cont(*) 'Is VELY used in atmospheric boundary? : ', atmos_boundary_use_vely
530  log_info_cont(*) 'Is PT used in atmospheric boundary? : ', atmos_boundary_use_pt
531  log_info_cont(*) 'Is DENS used in atmospheric boundary? : ', atmos_boundary_use_dens
532  log_info_cont(*) 'Is QV used in atmospheric boundary? : ', atmos_boundary_use_qv
533  log_info_cont(*) 'Is QHYD used in atmospheric boundary? : ', atmos_boundary_use_qhyd
534  log_info_cont(*) 'Is CHEM used in atmospheric boundary? : ', atmos_boundary_use_chem
535  log_newline
536  log_info_cont(*) 'Atmospheric boundary VELZ values : ', atmos_boundary_value_velz
537  log_info_cont(*) 'Atmospheric boundary VELX values : ', atmos_boundary_value_velx
538  log_info_cont(*) 'Atmospheric boundary VELY values : ', atmos_boundary_value_vely
539  log_info_cont(*) 'Atmospheric boundary PT values : ', atmos_boundary_value_pt
540  log_info_cont(*) 'Atmospheric boundary QTRC values : ', atmos_boundary_value_qtrc
541  log_newline
542  log_info_cont(*) 'Atmospheric boundary smoother factor : ', atmos_boundary_smoother_fact
543  log_info_cont(*) 'Atmospheric boundary z-fraction : ', atmos_boundary_fracz
544  log_info_cont(*) 'Atmospheric boundary x-fraction : ', atmos_boundary_fracx
545  log_info_cont(*) 'Atmospheric boundary y-fraction : ', atmos_boundary_fracy
546  log_info_cont(*) 'Atmospheric boundary z-relaxation time : ', atmos_boundary_tauz
547  log_info_cont(*) 'Atmospheric boundary x-relaxation time : ', atmos_boundary_taux
548  log_info_cont(*) 'Atmospheric boundary y-relaxation time : ', atmos_boundary_tauy
549  log_newline
550  log_info_cont(*) 'Atmospheric boundary update dt : ', atmos_boundary_update_dt
551  log_info_cont(*) 'Atmospheric boundary start date : ', atmos_boundary_start_date(:)
552  log_newline
553  log_info_cont(*) 'Linear profile in vertically relax region : ', atmos_boundary_linear_v
554  log_info_cont(*) 'Linear profile in horizontally relax region : ', atmos_boundary_linear_h
555  log_info_cont(*) 'Non-linear factor in horizontally relax region : ', atmos_boundary_exp_h
556  log_newline
557  log_info_cont(*) 'Online nesting for lateral boundary : ', atmos_boundary_online
558 
559  log_info_cont(*) 'Does lateral boundary exist in this domain? : ', l_bnd
560  if ( l_bnd ) then
561  log_info_cont(*) 'Lateral boundary interporation type : ', trim(atmos_boundary_interp_type)
562  endif
563  log_newline
564  log_info_cont(*) 'Is grid nudging used for VELX & VELY? : ', atmos_grid_nudging_uv
565  log_info_cont(*) 'Is grid nudging used for POTT? : ', atmos_grid_nudging_pt
566  log_info_cont(*) 'Is grid nudging used for QV? : ', atmos_grid_nudging_qv
567  log_info_cont(*) 'Relaxation time for grid nudging : ', atmos_grid_nudging_tau
568 
569  log_info_cont(*) 'Density adjustment : ', atmos_boundary_dens_adjust
570  if ( atmos_boundary_dens_adjust ) then
571  log_info_cont(*) 'Density relaxation time : ', atmos_boundary_dens_adjust_tau
572  end if
573 
574  if ( online_boundary_diagqhyd ) then
575  allocate( q_work(ka,ia,ja,nestqa) )
576  end if
577 
578  allocate( zero_x(ka,ja), zero_y(ka,ia) )
579  !$omp parallel do
580  do j = js, je
581  do k = ks, ke
582  zero_x(k,j) = 0.0_rp
583  end do
584  end do
585  !$omp parallel do
586  do i = is, ie
587  do k = ks, ke
588  zero_y(k,i) = 0.0_rp
589  end do
590  end do
591 
592 
593 
594  return
595  end subroutine atmos_boundary_driver_setup
596 
597  !-----------------------------------------------------------------------------
599  subroutine atmos_boundary_driver_set
600  use mod_atmos_vars, only: &
601  dens, &
602  momz, &
603  momx, &
604  momy, &
605  rhot, &
606  qtrc, &
607  qv, &
608  qe
609  use mod_atmos_phy_mp_vars, only: &
610  qs_mp, &
611  qe_mp
612  implicit none
613 
614  if ( do_parent_process ) then !online [parent]
615  call atmos_boundary_firstsend( &
616  dens, momz, momx, momy, rhot, qtrc(:,:,:,qs_mp:qe_mp), qv, qe )
617  end if
618 
619  if ( l_bnd ) then
620 
621  ! initialize boundary value (reading file or waiting parent domain)
622  if ( do_daughter_process ) then
624  else
625  if ( atmos_boundary_in_basename /= '' ) then
627  endif
628  endif
629 
630  elseif ( atmos_boundary_type == 'CONST' ) then
631 
632  call atmos_boundary_generate
633 
634  elseif ( atmos_boundary_type == 'INIT' ) then
635 
636  call atmos_boundary_setinitval( dens, & ! [IN]
637  momz, & ! [IN]
638  momx, & ! [IN]
639  momy, & ! [IN]
640  rhot, & ! [IN]
641  qtrc ) ! [IN]
642  endif
643 
644  if( atmos_boundary_out_basename /= '' ) then
645  call atmos_boundary_write
646  endif
647 
648  if ( atmos_boundary_update_flag ) then
649 
650  call history_bnd( &
657  end if
658 
659  return
660  end subroutine atmos_boundary_driver_set
661 
662  !-----------------------------------------------------------------------------
664  subroutine atmos_boundary_var_fillhalo
665  use scale_comm_cartesc, only: &
666  comm_vars8, &
667  comm_wait
668  implicit none
669 
670  integer :: i, j, iq
671  !---------------------------------------------------------------------------
672 
673  do j = 1, ja
674  do i = 1, ia
680 
686 
687  do iq = 1, bnd_qa
688  atmos_boundary_qtrc( 1:ks-1,i,j,iq) = atmos_boundary_qtrc(ks,i,j,iq)
689  atmos_boundary_qtrc(ke+1:ka, i,j,iq) = atmos_boundary_qtrc(ke,i,j,iq)
690  end do
691  end do
692  end do
693 
694  call comm_vars8( atmos_boundary_dens(:,:,:), 1 )
695  call comm_vars8( atmos_boundary_velz(:,:,:), 2 )
696  call comm_vars8( atmos_boundary_velx(:,:,:), 3 )
697  call comm_vars8( atmos_boundary_vely(:,:,:), 4 )
698  call comm_vars8( atmos_boundary_pott(:,:,:), 5 )
699  do iq = 1, bnd_qa
700  call comm_vars8( atmos_boundary_qtrc(:,:,:,iq), 5+iq )
701  end do
702 
703  call comm_wait ( atmos_boundary_dens(:,:,:), 1, .false. )
704  call comm_wait ( atmos_boundary_velz(:,:,:), 2, .false. )
705  call comm_wait ( atmos_boundary_velx(:,:,:), 3, .false. )
706  call comm_wait ( atmos_boundary_vely(:,:,:), 4, .false. )
707  call comm_wait ( atmos_boundary_pott(:,:,:), 5, .false. )
708  do iq = 1, bnd_qa
709  call comm_wait ( atmos_boundary_qtrc(:,:,:,iq), 5+iq, .false. )
710  end do
711 
712  return
713  end subroutine atmos_boundary_var_fillhalo
714 
715  !-----------------------------------------------------------------------------
717  subroutine atmos_boundary_alpha_fillhalo
718  use scale_comm_cartesc, only: &
719  comm_vars8, &
720  comm_wait
721  implicit none
722 
723  integer :: i, j, iq
724  !---------------------------------------------------------------------------
725 
726  do j = 1, ja
727  do i = 1, ia
733 
739 
740  do iq = 1, bnd_qa
743  end do
744  enddo
745  enddo
746 
747  call comm_vars8( atmos_boundary_alpha_dens(:,:,:), 1 )
748  call comm_vars8( atmos_boundary_alpha_velz(:,:,:), 2 )
749  call comm_vars8( atmos_boundary_alpha_velx(:,:,:), 3 )
750  call comm_vars8( atmos_boundary_alpha_vely(:,:,:), 4 )
751  call comm_vars8( atmos_boundary_alpha_pott(:,:,:), 5 )
752  do iq = 1, bnd_qa
753  call comm_vars8( atmos_boundary_alpha_qtrc(:,:,:,iq), 5+iq )
754  end do
755 
756  call comm_wait ( atmos_boundary_alpha_dens(:,:,:), 1, .false. )
757  call comm_wait ( atmos_boundary_alpha_velz(:,:,:), 2, .false. )
758  call comm_wait ( atmos_boundary_alpha_velx(:,:,:), 3, .false. )
759  call comm_wait ( atmos_boundary_alpha_vely(:,:,:), 4, .false. )
760  call comm_wait ( atmos_boundary_alpha_pott(:,:,:), 5, .false. )
761  do iq = 1, bnd_qa
762  call comm_wait ( atmos_boundary_alpha_qtrc(:,:,:,iq), 5+iq, .false. )
763  end do
764 
765  return
766  end subroutine atmos_boundary_alpha_fillhalo
767 
768  !-----------------------------------------------------------------------------
770  subroutine atmos_boundary_ref_fillhalo( &
771  ref_idx )
772  use scale_comm_cartesc, only: &
773  comm_vars8, &
774  comm_wait
775  implicit none
776 
777  ! arguments
778  integer, intent(in) :: ref_idx
779 
780  ! works
781  integer :: i, j, iq
782  !---------------------------------------------------------------------------
783 
784  !$omp parallel do default(none) private(i,j,iq) OMP_SCHEDULE_ collapse(2) &
785  !$omp shared(JSB,JEB,ISB,IEB,ATMOS_BOUNDARY_ref_DENS,ref_idx,ATMOS_BOUNDARY_ref_VELZ) &
786  !$omp shared(ATMOS_BOUNDARY_ref_VELX,ATMOS_BOUNDARY_ref_VELY,ATMOS_BOUNDARY_ref_POTT) &
787  !$omp shared(KA,KS,BND_QA,ATMOS_BOUNDARY_ref_QTRC,KE)
788  do j = jsb, jeb
789  do i = isb, ieb
790  atmos_boundary_ref_dens( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_dens(ks,i,j,ref_idx)
791  atmos_boundary_ref_velz( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_velz(ks,i,j,ref_idx)
792  atmos_boundary_ref_velx( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_velx(ks,i,j,ref_idx)
793  atmos_boundary_ref_vely( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_vely(ks,i,j,ref_idx)
794  atmos_boundary_ref_pott( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_pott(ks,i,j,ref_idx)
795 
796  atmos_boundary_ref_dens(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_dens(ke,i,j,ref_idx)
797  atmos_boundary_ref_velz(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_velz(ke,i,j,ref_idx)
798  atmos_boundary_ref_velx(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_velx(ke,i,j,ref_idx)
799  atmos_boundary_ref_vely(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_vely(ke,i,j,ref_idx)
800  atmos_boundary_ref_pott(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_pott(ke,i,j,ref_idx)
801 
802  do iq = 1, bnd_qa
803  atmos_boundary_ref_qtrc( 1:ks-1,i,j,iq,ref_idx) = atmos_boundary_ref_qtrc(ks,i,j,iq,ref_idx)
804  atmos_boundary_ref_qtrc(ke+1:ka, i,j,iq,ref_idx) = atmos_boundary_ref_qtrc(ke,i,j,iq,ref_idx)
805  end do
806  end do
807  end do
808 
809  call comm_vars8( atmos_boundary_ref_dens(:,:,:,ref_idx), 1 )
810  call comm_vars8( atmos_boundary_ref_velz(:,:,:,ref_idx), 2 )
811  call comm_vars8( atmos_boundary_ref_velx(:,:,:,ref_idx), 3 )
812  call comm_vars8( atmos_boundary_ref_vely(:,:,:,ref_idx), 4 )
813  call comm_vars8( atmos_boundary_ref_pott(:,:,:,ref_idx), 5 )
814 
815  do iq = 1, bnd_qa
816  call comm_vars8( atmos_boundary_ref_qtrc(:,:,:,iq,ref_idx), 5+iq )
817  end do
818 
819  call comm_wait ( atmos_boundary_ref_dens(:,:,:,ref_idx), 1, .false. )
820  call comm_wait ( atmos_boundary_ref_velz(:,:,:,ref_idx), 2, .false. )
821  call comm_wait ( atmos_boundary_ref_velx(:,:,:,ref_idx), 3, .false. )
822  call comm_wait ( atmos_boundary_ref_vely(:,:,:,ref_idx), 4, .false. )
823  call comm_wait ( atmos_boundary_ref_pott(:,:,:,ref_idx), 5, .false. )
824 
825  do iq = 1, bnd_qa
826  call comm_wait ( atmos_boundary_ref_qtrc(:,:,:,iq,ref_idx), 5+iq, .false. )
827  end do
828 
829  return
830  end subroutine atmos_boundary_ref_fillhalo
831 
832  !-----------------------------------------------------------------------------
834  subroutine atmos_boundary_setalpha
835  use scale_const, only: &
836  eps => const_eps, &
837  pi => const_pi
838  use scale_atmos_grid_cartesc, only: &
839  cbfz => atmos_grid_cartesc_cbfz, &
840  cbfx => atmos_grid_cartesc_cbfx, &
841  cbfy => atmos_grid_cartesc_cbfy, &
842  fbfz => atmos_grid_cartesc_fbfz, &
843  fbfx => atmos_grid_cartesc_fbfx, &
845  use scale_comm_cartesc_nest, only: &
847  use scale_atmos_hydrometeor, only: &
848  i_qv
849  real(rp) :: coef_z, alpha_z1, alpha_z2
850  real(rp) :: coef_x, alpha_x1, alpha_x2
851  real(rp) :: coef_y, alpha_y1, alpha_y2
852  real(rp) :: alpha_zm, alpha_xm, alpha_ym
853  real(rp) :: ee1, ee2
854 
855  real(rp) :: alpha_nug ! grid nudging in inner domain
856 
857  integer :: i, j, k, iq
858  !---------------------------------------------------------------------------
859 
860  ! check invalid fraction
861  atmos_boundary_fracz = max( min( atmos_boundary_fracz, 1.0_rp ), eps )
862  atmos_boundary_fracx = max( min( atmos_boundary_fracx, 1.0_rp ), eps )
863  atmos_boundary_fracy = max( min( atmos_boundary_fracy, 1.0_rp ), eps )
864 
865  if ( atmos_boundary_tauz <= 0.0_rp ) then ! invalid tau
866  coef_z = 0.0_rp
867  else
868  coef_z = 1.0_rp / atmos_boundary_tauz
869  endif
870 
871  if ( atmos_boundary_taux <= 0.0_rp ) then ! invalid tau
872  coef_x = 0.0_rp
873  else
874  coef_x = 1.0_rp / atmos_boundary_taux
875  endif
876 
877  if ( atmos_boundary_tauy <= 0.0_rp ) then ! invalid tau
878  coef_y = 0.0_rp
879  else
880  coef_y = 1.0_rp / atmos_boundary_tauy
881  endif
882 
883  if ( atmos_grid_nudging_tau <= 0.0_rp ) then ! invalid tau
884  alpha_nug = 0.0_rp
885  else
886  alpha_nug = 1.0_rp / atmos_grid_nudging_tau
887  endif
888 
889  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
890  !$omp shared(JA,IA,KS,KE,CBFZ,ATMOS_BOUNDARY_FRACZ,FBFZ,ATMOS_BOUNDARY_LINEAR_V,coef_z,CBFX) &
891  !$omp shared(ATMOS_BOUNDARY_FRACX,PI,FBFX,ATMOS_BOUNDARY_LINEAR_H,coef_x) &
892  !$omp shared(ATMOS_BOUNDARY_EXP_H,CBFY,ATMOS_BOUNDARY_FRACY,FBFY,coef_y,l_bnd) &
893  !$omp shared(do_daughter_process) &
894  !$omp shared(ONLINE_USE_VELZ,ATMOS_BOUNDARY_USE_VELZ,ATMOS_BOUNDARY_alpha_VELZ,ATMOS_BOUNDARY_ALPHAFACT_VELZ) &
895  !$omp shared(ATMOS_BOUNDARY_USE_DENS,ATMOS_BOUNDARY_alpha_DENS,ATMOS_BOUNDARY_ALPHAFACT_DENS) &
896  !$omp shared(ATMOS_BOUNDARY_USE_VELX,ATMOS_BOUNDARY_alpha_VELX,ATMOS_BOUNDARY_ALPHAFACT_VELX) &
897  !$omp shared(ATMOS_BOUNDARY_USE_VELY,ATMOS_BOUNDARY_alpha_VELY,ATMOS_BOUNDARY_ALPHAFACT_VELY) &
898  !$omp shared(ATMOS_BOUNDARY_USE_PT,ATMOS_BOUNDARY_alpha_POTT,ATMOS_BOUNDARY_ALPHAFACT_PT) &
899  !$omp shared(ATMOS_BOUNDARY_USE_QV,ATMOS_BOUNDARY_alpha_QTRC,ATMOS_BOUNDARY_ALPHAFACT_QTRC) &
900  !$omp shared(ATMOS_BOUNDARY_DENS_ADJUST,ATMOS_BOUNDARY_DENS_ADJUST_tau,ATMOS_BOUNDARY_UPDATE_DT) &
901  !$omp shared(BND_QA,BND_IQ,I_QV) &
902  !$omp shared(ATMOS_GRID_NUDGING_uv,ATMOS_GRID_NUDGING_pt,ATMOS_GRID_NUDGING_qv) &
903  !$omp shared(alpha_nug) &
904  !$omp private(i,j,k,iq) &
905  !$omp private(ee1,ee2,alpha_z1,alpha_z2,alpha_x1,alpha_x2,alpha_y1,alpha_y2,alpha_zm,alpha_xm,alpha_ym)
906  do j = 1, ja
907  do i = 1, ia
908  do k = ks, ke
909  ee1 = cbfz(k)
910  if ( ee1 <= 1.0_rp - atmos_boundary_fracz ) then
911  ee1 = 0.0_rp
912  else
913  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracz ) / atmos_boundary_fracz
914  endif
915 
916  ee2 = fbfz(k)
917  if ( ee2 <= 1.0_rp - atmos_boundary_fracz ) then
918  ee2 = 0.0_rp
919  else
920  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracz ) / atmos_boundary_fracz
921  endif
922 
923  if ( .not. atmos_boundary_linear_v ) then
924  if ( ee1 > 0.0_rp .AND. ee1 <= 0.5_rp ) then
925  ee1 = 0.5_rp * ( 1.0_rp - cos( ee1*pi ) )
926  elseif( ee1 > 0.5_rp .AND. ee1 <= 1.0_rp ) then
927  ee1 = 0.5_rp * ( 1.0_rp + sin( (ee1-0.5_rp)*pi ) )
928  endif
929  if ( ee2 > 0.0_rp .AND. ee2 <= 0.5_rp ) then
930  ee2 = 0.5_rp * ( 1.0_rp - cos( ee2*pi ) )
931  elseif( ee2 > 0.5_rp .AND. ee2 <= 1.0_rp ) then
932  ee2 = 0.5_rp * ( 1.0_rp + sin( (ee2-0.5_rp)*pi ) )
933  endif
934  endif
935 
936  alpha_z1 = coef_z * ee1
937  alpha_z2 = coef_z * ee2
938  alpha_zm = ee1 / atmos_boundary_dens_adjust_tau
939 
940 
941  ee1 = cbfx(i)
942  if ( ee1 <= 1.0_rp - atmos_boundary_fracx ) then
943  ee1 = 0.0_rp
944  else
945  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracx ) / atmos_boundary_fracx
946  endif
947 
948  ee2 = fbfx(i)
949  if ( ee2 <= 1.0_rp - atmos_boundary_fracx ) then
950  ee2 = 0.0_rp
951  else
952  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracx ) / atmos_boundary_fracx
953  endif
954 
955  if ( .not. atmos_boundary_linear_h ) then
956  ee1 = ee1 * exp( -(1.0_rp-ee1) * atmos_boundary_exp_h )
957  ee2 = ee2 * exp( -(1.0_rp-ee2) * atmos_boundary_exp_h )
958  end if
959 
960  alpha_x1 = coef_x * ee1
961  alpha_x2 = coef_x * ee2
962  alpha_xm = ee1 / atmos_boundary_dens_adjust_tau
963 
964 
965  ee1 = cbfy(j)
966  if ( ee1 <= 1.0_rp - atmos_boundary_fracy ) then
967  ee1 = 0.0_rp
968  else
969  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracy ) / atmos_boundary_fracy
970  endif
971 
972  ee2 = fbfy(j)
973  if ( ee2 <= 1.0_rp - atmos_boundary_fracy ) then
974  ee2 = 0.0_rp
975  else
976  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracy ) / atmos_boundary_fracy
977  endif
978 
979  if ( .not. atmos_boundary_linear_h ) then
980  ee1 = ee1 * exp( -(1.0_rp-ee1) * atmos_boundary_exp_h )
981  ee2 = ee2 * exp( -(1.0_rp-ee2) * atmos_boundary_exp_h )
982  end if
983 
984  alpha_y1 = coef_y * ee1
985  alpha_y2 = coef_y * ee2
986  alpha_ym = ee1 / atmos_boundary_dens_adjust_tau
987 
988 
989  if ( l_bnd ) then
990  if ( ( do_daughter_process .and. online_use_velz ) &
991  .or. ( .not. do_daughter_process .and. atmos_boundary_use_velz ) ) then
992  atmos_boundary_alpha_velz(k,i,j) = max( alpha_z2, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
993  else
994  atmos_boundary_alpha_velz(k,i,j) = 0.0_rp
995  end if
996  if ( atmos_boundary_dens_adjust ) then
997  atmos_boundary_alpha_dens(k,i,j) = max( alpha_zm, alpha_xm, alpha_ym )
998  else
999  atmos_boundary_alpha_dens(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_dens
1000  end if
1001  atmos_boundary_alpha_velx(k,i,j) = max( alpha_z1, alpha_x2, alpha_y1 ) * atmos_boundary_alphafact_velx
1002  atmos_boundary_alpha_vely(k,i,j) = max( alpha_z1, alpha_x1, alpha_y2 ) * atmos_boundary_alphafact_vely
1003  atmos_boundary_alpha_pott(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_pt
1004  do iq = 1, bnd_qa
1005  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
1006  end do
1007  else
1008  if ( atmos_boundary_use_dens ) then
1009  atmos_boundary_alpha_dens(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_dens
1010  else
1011  atmos_boundary_alpha_dens(k,i,j) = 0.0_rp
1012  end if
1013  if ( atmos_boundary_use_velz ) then
1014  atmos_boundary_alpha_velz(k,i,j) = max( alpha_z2, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
1015  else
1016  atmos_boundary_alpha_velz(k,i,j) = 0.0_rp
1017  end if
1018  if ( atmos_boundary_use_velx ) then
1019  atmos_boundary_alpha_velx(k,i,j) = max( alpha_z1, alpha_x2, alpha_y1 ) * atmos_boundary_alphafact_velx
1020  else
1021  atmos_boundary_alpha_velx(k,i,j) = 0.0_rp
1022  end if
1023  if ( atmos_boundary_use_vely ) then
1024  atmos_boundary_alpha_vely(k,i,j) = max( alpha_z1, alpha_x1, alpha_y2 ) * atmos_boundary_alphafact_vely
1025  else
1026  atmos_boundary_alpha_vely(k,i,j) = 0.0_rp
1027  end if
1028  if ( atmos_boundary_use_pt ) then
1029  atmos_boundary_alpha_pott(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_pt
1030  else
1031  atmos_boundary_alpha_pott(k,i,j) = 0.0_rp
1032  end if
1033  do iq = 1, bnd_qa
1034  if ( i_qv > 0 .and. iq == bnd_iq(i_qv) ) then
1035  if ( atmos_boundary_use_qv ) then
1036  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
1037  else
1038  atmos_boundary_alpha_qtrc(k,i,j,iq) = 0.0_rp
1039  endif
1040  else
1041  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
1042  end if
1043  end do
1044 
1045  end if
1046 
1047  ! internal grid nudging
1048  if ( atmos_grid_nudging_uv ) then
1049  atmos_boundary_alpha_velx(k,i,j) = max( atmos_boundary_alpha_velx(k,i,j), alpha_nug )
1050  atmos_boundary_alpha_vely(k,i,j) = max( atmos_boundary_alpha_vely(k,i,j), alpha_nug )
1051  endif
1052  if ( atmos_grid_nudging_pt ) then
1053  atmos_boundary_alpha_pott(k,i,j) = max( atmos_boundary_alpha_pott(k,i,j), alpha_nug )
1054  endif
1055  if ( atmos_grid_nudging_qv ) then
1056  atmos_boundary_alpha_qtrc(k,i,j,1) = max( atmos_boundary_alpha_qtrc(k,i,j,1), alpha_nug )
1057  endif
1058 
1059  enddo
1060  enddo
1061  enddo
1062 
1063 
1064  call atmos_boundary_alpha_fillhalo
1065 
1066  return
1067  end subroutine atmos_boundary_setalpha
1068 
1069  !-----------------------------------------------------------------------------
1071  subroutine atmos_boundary_setinitval( &
1072  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC )
1073  implicit none
1074 
1075  real(rp), intent(in) :: dens(ka,ia,ja)
1076  real(rp), intent(in) :: momz(ka,ia,ja)
1077  real(rp), intent(in) :: momx(ka,ia,ja)
1078  real(rp), intent(in) :: momy(ka,ia,ja)
1079  real(rp), intent(in) :: rhot(ka,ia,ja)
1080  real(rp), intent(in) :: qtrc(ka,ia,ja,bnd_qa)
1081 
1082  integer :: i, j, k, iq, iqb
1083  !---------------------------------------------------------------------------
1084 
1085  do j = 1, ja
1086  do i = 1, ia
1087  do k = ks, ke
1088  atmos_boundary_dens(k,i,j) = dens(k,i,j)
1089  atmos_boundary_velz(k,i,j) = momz(k,i,j) / ( dens(k,i,j)+dens(k+1,i, j ) ) * 2.0_rp
1090  atmos_boundary_pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
1091  do iq = 1, qa
1092  iqb = bnd_iq(iq)
1093  if ( iqb > 0 ) atmos_boundary_qtrc(k,i,j,iqb) = qtrc(k,i,j,iq)
1094  end do
1095  enddo
1096  enddo
1097  enddo
1098 
1099  do j = 1, ja
1100  do i = 1, ia-1
1101  do k = ks, ke
1102  atmos_boundary_velx(k,i,j) = momx(k,i,j) / ( dens(k,i,j)+dens(k, i+1,j ) ) * 2.0_rp
1103  enddo
1104  enddo
1105  enddo
1106  do j = 1, ja
1107  do k = ks, ke
1108  atmos_boundary_velx(k,ia,j) = momx(k,ia,j) / dens(k,ia,j)
1109  enddo
1110  enddo
1111 
1112  do j = 1, ja-1
1113  do i = 1, ia
1114  do k = ks, ke
1115  atmos_boundary_vely(k,i,j) = momy(k,i,j) / ( dens(k,i,j)+dens(k, i, j+1) ) * 2.0_rp
1116  enddo
1117  enddo
1118  enddo
1119  do i = 1, ia
1120  do k = ks, ke
1121  atmos_boundary_vely(k,i,ja) = momy(k,i,ja) / dens(k,i,ja)
1122  enddo
1123  enddo
1124 
1125  call atmos_boundary_var_fillhalo
1126 
1127  return
1128  end subroutine atmos_boundary_setinitval
1129 
1130  !-----------------------------------------------------------------------------
1132  subroutine atmos_boundary_read
1133  use scale_prc, only: &
1134  prc_abort
1135  use scale_file_cartesc, only: &
1137  file_cartesc_check_coordinates, &
1138  file_cartesc_read, &
1140  implicit none
1141 
1142  integer :: fid
1143  integer :: iq, iqb
1144  !---------------------------------------------------------------------------
1145 
1146  call file_cartesc_open( atmos_boundary_in_basename, fid, aggregate=atmos_boundary_in_aggregate )
1147 
1148  if ( atmos_boundary_in_check_coordinates ) then
1149  call file_cartesc_check_coordinates( fid, atmos=.true. )
1150  end if
1151 
1152  if ( atmos_boundary_use_dens &
1153  .OR. atmos_boundary_use_velz &
1154  .OR. atmos_boundary_use_velx &
1155  .OR. atmos_boundary_use_vely &
1156  .OR. atmos_boundary_use_pt &
1157  ) then
1158  call file_cartesc_read( fid, 'DENS', 'ZXY', atmos_boundary_dens(:,:,:) )
1159  end if
1160  if ( atmos_boundary_use_dens ) then
1161  call file_cartesc_read( fid, 'ALPHA_DENS', 'ZXY', atmos_boundary_alpha_dens(:,:,:) )
1162  endif
1163 
1164  if ( atmos_boundary_use_velz ) then
1165  call file_cartesc_read( fid, 'VELZ', 'ZHXY', atmos_boundary_velz(:,:,:) )
1166  call file_cartesc_read( fid, 'ALPHA_VELZ', 'ZHXY', atmos_boundary_alpha_velz(:,:,:) )
1167  endif
1168 
1169  if ( atmos_boundary_use_velx ) then
1170  call file_cartesc_read( fid, 'VELX', 'ZXHY', atmos_boundary_velx(:,:,:) )
1171  call file_cartesc_read( fid, 'ALPHA_VELX', 'ZXHY', atmos_boundary_alpha_velx(:,:,:) )
1172  endif
1173 
1174  if ( atmos_boundary_use_vely ) then
1175  call file_cartesc_read( fid, 'VELY', 'ZXYH', atmos_boundary_vely(:,:,:) )
1176  call file_cartesc_read( fid, 'ALPHA_VELY', 'ZXYH', atmos_boundary_alpha_vely(:,:,:) )
1177  endif
1178 
1179  if ( atmos_boundary_use_pt ) then
1180  call file_cartesc_read( fid, 'PT', 'ZXY', atmos_boundary_pott(:,:,:) )
1181  call file_cartesc_read( fid, 'ALPHA_PT', 'ZXY', atmos_boundary_alpha_pott(:,:,:) )
1182  endif
1183 
1184  do iq = 1, qa
1185  iqb = bnd_iq(iq)
1186  if ( iqb > 0 ) then
1187  call file_cartesc_read( fid, tracer_name(iq), 'ZXY', atmos_boundary_qtrc(:,:,:,iqb) )
1188  call file_cartesc_read( fid, 'ALPHA_'//trim(tracer_name(iq)), 'ZXY', atmos_boundary_alpha_qtrc(:,:,:,iqb) )
1189  endif
1190  end do
1191 
1192  call file_cartesc_close( fid )
1193 
1194 
1195  call atmos_boundary_var_fillhalo
1196  call atmos_boundary_alpha_fillhalo
1197 
1198  return
1199  end subroutine atmos_boundary_read
1200 
1201  !-----------------------------------------------------------------------------
1203  subroutine atmos_boundary_write
1204  use scale_file_cartesc, only: &
1208  file_cartesc_write_var, &
1210  use scale_comm_cartesc_nest, only: &
1212  implicit none
1213 
1214  integer :: fid
1215  integer :: vid_dens, vid_a_dens
1216  integer :: vid_velz, vid_a_velz
1217  integer :: vid_velx, vid_a_velx
1218  integer :: vid_vely, vid_a_vely
1219  integer :: vid_pott, vid_a_pott
1220  integer :: vid_qtrc(bnd_qa), vid_a_qtrc(bnd_qa)
1221  integer :: iq, iqb
1222  !---------------------------------------------------------------------------
1223 
1224  call file_cartesc_create( atmos_boundary_out_basename, atmos_boundary_out_title, & ! [IN]
1225  atmos_boundary_out_dtype, & ! [IN]
1226  fid, & ! [OUT]
1227  aggregate = atmos_boundary_out_aggregate ) ! [IN]
1228 
1229  if ( atmos_boundary_use_dens &
1230  .OR. atmos_boundary_use_velz &
1231  .OR. atmos_boundary_use_velx &
1232  .OR. atmos_boundary_use_vely &
1233  .OR. atmos_boundary_use_pt &
1234  .OR. l_bnd ) then
1235  call file_cartesc_def_var( fid, 'DENS', 'Reference Density', 'kg/m3', 'ZXY', atmos_boundary_out_dtype, vid_dens )
1236  else
1237  vid_dens = -1
1238  end if
1239 
1240  if ( atmos_boundary_use_dens .OR. l_bnd ) then
1241  call file_cartesc_def_var( fid, 'ALPHA_DENS', 'Alpha for DENS', '1', 'ZXY', atmos_boundary_out_dtype, vid_a_dens )
1242  else
1243  vid_a_dens = -1
1244  end if
1245 
1246  if ( atmos_boundary_use_velz .OR. (l_bnd .AND. online_use_velz) ) then
1247  call file_cartesc_def_var( fid, 'VELZ', 'Reference Velocity w', 'm/s', 'ZHXY', atmos_boundary_out_dtype, vid_velz )
1248  call file_cartesc_def_var( fid, 'ALPHA_VELZ', 'Alpha for VELZ', '1', 'ZHXY', atmos_boundary_out_dtype, vid_a_velz )
1249  else
1250  vid_velz = -1
1251  end if
1252 
1253  if ( atmos_boundary_use_velx .OR. l_bnd ) then
1254  call file_cartesc_def_var( fid, 'VELX', 'Reference Velocity u', 'm/s', 'ZXHY', atmos_boundary_out_dtype, vid_velx )
1255  call file_cartesc_def_var( fid, 'ALPHA_VELX', 'Alpha for VELX', '1', 'ZXHY', atmos_boundary_out_dtype, vid_a_velx )
1256  else
1257  vid_velx = -1
1258  end if
1259 
1260  if ( atmos_boundary_use_vely .OR. l_bnd ) then
1261  call file_cartesc_def_var( fid, 'VELY', 'Reference Velocity y', 'm/s', 'ZXYH', atmos_boundary_out_dtype, vid_vely )
1262  call file_cartesc_def_var( fid, 'ALPHA_VELY', 'Alpha for VELY', '1', 'ZXYH', atmos_boundary_out_dtype, vid_a_vely )
1263  else
1264  vid_vely = -1
1265  end if
1266 
1267  if ( atmos_boundary_use_pt .OR. l_bnd ) then
1268  call file_cartesc_def_var( fid, 'PT', 'Reference PT', 'K', 'ZXY', atmos_boundary_out_dtype, vid_pott )
1269  call file_cartesc_def_var( fid, 'ALPHA_PT', 'Alpha for PT', '1', 'ZXY', atmos_boundary_out_dtype, vid_a_pott )
1270  else
1271  vid_pott = -1
1272  end if
1273 
1274  do iq = 1, qa
1275  iqb = bnd_iq(iq)
1276  if ( iqb > 0 ) then
1277  call file_cartesc_def_var( fid, tracer_name(iq), 'Reference '//trim(tracer_name(iq)), tracer_unit(iq), 'ZXY', atmos_boundary_out_dtype, vid_qtrc(iqb) )
1278  call file_cartesc_def_var( fid, 'ALPHA_'//trim(tracer_name(iq)), 'Alpha for '//trim(tracer_name(iq)), '1', 'ZXY', atmos_boundary_out_dtype, vid_a_qtrc(iqb) )
1279  end if
1280  end do
1281 
1282 
1283  call file_cartesc_enddef( fid )
1284 
1285 
1286  if ( vid_dens > 0 ) then
1287  call file_cartesc_write_var( fid, vid_dens, atmos_boundary_dens(:,:,:), 'DENS', 'ZXY' )
1288  end if
1289 
1290  if ( vid_a_dens > 0 ) then
1291  call file_cartesc_write_var( fid, vid_a_dens, atmos_boundary_alpha_dens(:,:,:), 'ALPHA_DENS', 'ZXY' )
1292  end if
1293 
1294  if ( vid_velz > 0 ) then
1295  call file_cartesc_write_var( fid, vid_velz, atmos_boundary_velz(:,:,:), 'VELZ', 'ZHXY' )
1296  call file_cartesc_write_var( fid, vid_a_velz, atmos_boundary_alpha_velz(:,:,:), 'ALPHA_VELZ', 'ZHXY' )
1297  end if
1298 
1299  if ( vid_velx > 0 ) then
1300  call file_cartesc_write_var( fid, vid_velx, atmos_boundary_velx(:,:,:), 'VELX', 'ZXHY' )
1301  call file_cartesc_write_var( fid, vid_a_velx, atmos_boundary_alpha_velx(:,:,:), 'ALPHA_VELX', 'ZXHY' )
1302  endif
1303 
1304  if ( vid_vely > 0 ) then
1305  call file_cartesc_write_var( fid, vid_vely, atmos_boundary_vely(:,:,:), 'VELY', 'ZXYH' )
1306  call file_cartesc_write_var( fid, vid_a_vely, atmos_boundary_alpha_vely(:,:,:), 'ALPHA_VELY', 'ZXYH' )
1307  end if
1308 
1309  if ( vid_pott > 0 ) then
1310  call file_cartesc_write_var( fid, vid_pott, atmos_boundary_pott(:,:,:), 'PT', 'ZXY' )
1311  call file_cartesc_write_var( fid, vid_a_pott, atmos_boundary_alpha_pott(:,:,:), 'ALPHA_PT', 'ZXY' )
1312  end if
1313 
1314  do iqb = 1, bnd_qa
1315  if ( vid_qtrc(iqb) > 0 ) then
1316  call file_cartesc_write_var( fid, vid_qtrc(iqb), atmos_boundary_qtrc(:,:,:,iqb), tracer_name(iqb), 'ZXY' )
1317  call file_cartesc_write_var( fid, vid_a_qtrc(iqb), atmos_boundary_alpha_qtrc(:,:,:,iqb), 'ALPHA_'//trim(tracer_name(iqb)), 'ZXY' )
1318  end if
1319  end do
1320 
1321 
1322  call file_cartesc_close( fid )
1323 
1324 
1325  return
1326  end subroutine atmos_boundary_write
1327 
1328  !-----------------------------------------------------------------------------
1330  subroutine atmos_boundary_generate
1331  use scale_atmos_refstate, only: &
1333  implicit none
1334 
1335  integer :: i, j, k, iq
1336  !---------------------------------------------------------------------------
1337 
1338  !$omp parallel do collapse(2)
1339  do j = 1, ja
1340  do i = 1, ia
1341  do k = ks, ke
1343  atmos_boundary_velz(k,i,j) = atmos_boundary_value_velz
1344  atmos_boundary_velx(k,i,j) = atmos_boundary_value_velx
1345  atmos_boundary_vely(k,i,j) = atmos_boundary_value_vely
1346  atmos_boundary_pott(k,i,j) = atmos_boundary_value_pt
1347  do iq = 1, bnd_qa
1348  atmos_boundary_qtrc(k,i,j,iq) = atmos_boundary_value_qtrc
1349  end do
1350  enddo
1351  enddo
1352  enddo
1353 
1354  call atmos_boundary_var_fillhalo
1355 
1356  return
1357  end subroutine atmos_boundary_generate
1358 
1359  !-----------------------------------------------------------------------------
1361  subroutine atmos_boundary_initialize_file
1362  use scale_calendar, only: &
1366  use scale_time, only: &
1367  time_nowdate
1368  use scale_file_cartesc, only: &
1370  file_cartesc_check_coordinates
1371  implicit none
1372 
1373  integer :: boundary_time_startday
1374  real(dp) :: boundary_time_startsec
1375  real(dp) :: boundary_time_startms
1376  integer :: boundary_time_offset_year
1377  character(len=27) :: boundary_chardate
1378 
1379  if ( atmos_boundary_start_date(1) == -9999 ) then
1380  atmos_boundary_start_date = time_nowdate
1381  end if
1382 
1383  !--- calculate time of the initial step in boundary file [no offset]
1384  boundary_time_startms = 0.0_dp
1385  boundary_time_offset_year = 0
1386  call calendar_date2char( boundary_chardate, & ! [OUT]
1387  atmos_boundary_start_date(:), & ! [IN]
1388  boundary_time_startms ) ! [IN]
1389 
1390  call calendar_date2daysec( boundary_time_startday, & ! [OUT]
1391  boundary_time_startsec, & ! [OUT]
1392  atmos_boundary_start_date(:), & ! [IN]
1393  boundary_time_startms, & ! [IN]
1394  boundary_time_offset_year ) ! [IN]
1395 
1396  boundary_time_initdaysec = calendar_combine_daysec( boundary_time_startday, boundary_time_startsec )
1397 
1398  log_info("ATMOS_BOUNDARY_initialize_file",'(1x,A,A)') 'BOUNDARY START Date : ', boundary_chardate
1399 
1400  call file_cartesc_open( atmos_boundary_in_basename, atmos_boundary_fid )
1401 
1402  if ( atmos_boundary_in_check_coordinates ) then
1403  call file_cartesc_check_coordinates( atmos_boundary_fid, atmos=.true. )
1404  end if
1405 
1406  return
1407  end subroutine atmos_boundary_initialize_file
1408 
1409  !-----------------------------------------------------------------------------
1411  subroutine atmos_boundary_set_file
1412  use scale_prc, only: &
1413  prc_abort
1414  use scale_const, only: &
1415  pi => const_pi
1416  use scale_time, only: &
1417  time_nowdate, &
1418  time_dtsec
1419  use scale_calendar, only: &
1422  implicit none
1423 
1424  real(RP) :: bnd_DENS(KA,IA,JA) ! damping coefficient for DENS (0-1)
1425  real(RP) :: bnd_VELZ(KA,IA,JA) ! damping coefficient for VELZ (0-1)
1426  real(RP) :: bnd_VELX(KA,IA,JA) ! damping coefficient for VELX (0-1)
1427  real(RP) :: bnd_VELY(KA,IA,JA) ! damping coefficient for VELY (0-1)
1428  real(RP) :: bnd_POTT(KA,IA,JA) ! damping coefficient for POTT (0-1)
1429  real(RP) :: bnd_QTRC(KA,IA,JA,BND_QA) ! damping coefficient for QTRC (0-1)
1430 
1431  integer :: run_time_startdate(6)
1432  integer :: run_time_startday
1433  real(DP) :: run_time_startsec
1434  real(DP) :: run_time_startms
1435  integer :: run_time_offset_year
1436  real(DP) :: run_time_nowdaysec
1437 
1438  real(DP) :: boundary_diff_daysec
1439  real(RP) :: boundary_inc_offset
1440  integer :: fillgaps_steps
1441 
1442  real(RP) :: total
1443 
1444  integer :: i, j, k, n, iq
1445  !---------------------------------------------------------------------------
1446 
1447  if ( atmos_boundary_update_dt <= 0.0_dp ) then
1448  log_error("ATMOS_BOUNDARY_set_file",*) 'You need specify ATMOS_BOUNDARY_UPDATE_DT as larger than 0.0'
1449  call prc_abort
1450  endif
1451  update_nstep = nint( atmos_boundary_update_dt / time_dtsec )
1452  if ( abs(update_nstep * time_dtsec - atmos_boundary_update_dt) > 1e-10_dp ) then
1453  log_error("ATMOS_BOUNDARY_set_file",*) 'ATMOS_BOUNDARY_UPDATE_DT is not multiple of DT'
1454  call prc_abort
1455  end if
1456 
1457  !--- recalculate time of the run [no offset]
1458  run_time_startdate(:) = time_nowdate(:)
1459  run_time_startms = 0.0_dp
1460  run_time_offset_year = 0
1461 
1462  call calendar_date2daysec( run_time_startday, & ! [OUT]
1463  run_time_startsec, & ! [OUT]
1464  run_time_startdate(:), & ! [IN]
1465  run_time_startms, & ! [IN]
1466  run_time_offset_year ) ! [IN]
1467 
1468  run_time_nowdaysec = calendar_combine_daysec( run_time_startday, run_time_startsec )
1469 
1470  boundary_diff_daysec = run_time_nowdaysec - boundary_time_initdaysec
1471  boundary_timestep = 1 + int( boundary_diff_daysec / atmos_boundary_update_dt )
1472  boundary_inc_offset = mod( boundary_diff_daysec, atmos_boundary_update_dt )
1473  fillgaps_steps = int( boundary_inc_offset / time_dtsec )
1474 
1475  log_info("ATMOS_BOUNDARY_set_file",*) 'BOUNDARY TIMESTEP NUMBER FOR INIT:', boundary_timestep
1476  log_info("ATMOS_BOUNDARY_set_file",*) 'BOUNDARY OFFSET:', boundary_inc_offset
1477  log_info("ATMOS_BOUNDARY_set_file",*) 'BOUNDARY FILLGAPS STEPS:', fillgaps_steps
1478 
1479 
1480  if ( atmos_boundary_dens_adjust ) then
1481  allocate( offset_time_fact(0:update_nstep) )
1482  total = 0.0_rp
1483  !$omp parallel do reduction(+:total)
1484  do n = 0, update_nstep
1485  offset_time_fact(n) = 1.0_rp - cos( 2.0_rp * pi * ( n - 1 ) / update_nstep )
1486  total = total + offset_time_fact(n)
1487  end do
1488  total = total / update_nstep
1489  !$omp parallel do
1490  do n = 0, update_nstep
1491  offset_time_fact(n) = offset_time_fact(n) / total
1492  end do
1493  end if
1494 
1495  ! read boundary data from input file
1496  call atmos_boundary_update_file( ref_now )
1497 
1498  if ( atmos_boundary_dens_adjust ) then
1499  call calc_mass( ref_now )
1500  end if
1501 
1502  boundary_timestep = boundary_timestep + 1
1503  call atmos_boundary_update_file( ref_new )
1504 
1505  if ( atmos_boundary_dens_adjust ) then
1506  call calc_mass( ref_new )
1507  call set_offset
1508  end if
1509 
1510  ! copy now to old
1511  !$omp parallel do default(none) private(i,j,k,iq) OMP_SCHEDULE_ collapse(2) &
1512  !$omp shared(JA,IA,KS,KE,ATMOS_BOUNDARY_ref_DENS,ref_old,ref_now,ATMOS_BOUNDARY_ref_VELX) &
1513  !$omp shared(ATMOS_BOUNDARY_ref_VELY,ATMOS_BOUNDARY_ref_POTT,BND_QA,ATMOS_BOUNDARY_ref_QTRC)
1514  do j = 1, ja
1515  do i = 1, ia
1516  do k = ks, ke
1517  atmos_boundary_ref_dens(k,i,j,ref_old) = atmos_boundary_ref_dens(k,i,j,ref_now)
1518  atmos_boundary_ref_velx(k,i,j,ref_old) = atmos_boundary_ref_velx(k,i,j,ref_now)
1519  atmos_boundary_ref_vely(k,i,j,ref_old) = atmos_boundary_ref_vely(k,i,j,ref_now)
1520  atmos_boundary_ref_pott(k,i,j,ref_old) = atmos_boundary_ref_pott(k,i,j,ref_now)
1521  do iq = 1, bnd_qa
1522  atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now)
1523  end do
1524  end do
1525  end do
1526  end do
1527 
1528  if ( atmos_boundary_use_velz ) then
1529  !$omp parallel do collapse(2)
1530  do j = 1, ja
1531  do i = 1, ia
1532  do k = ks, ke
1533  atmos_boundary_ref_velz(k,i,j,:) = atmos_boundary_value_velz
1534  end do
1535  end do
1536  end do
1537  end if
1538 
1539  now_step = fillgaps_steps
1540 
1541  ! set boundary data
1542  call set_boundary( atmos_boundary_use_velz )
1543 
1544  return
1545  end subroutine atmos_boundary_set_file
1546 
1547  !-----------------------------------------------------------------------------
1549  subroutine atmos_boundary_initialize_online
1550  use scale_prc, only: &
1551  prc_abort
1552  use scale_comm_cartesc_nest, only: &
1554  online_use_velz, &
1555  parent_dtsec, &
1556  nestqa => comm_cartesc_nest_bnd_qa
1557  implicit none
1558 
1559  ! parameters
1560  integer, parameter :: handle = 2
1561 
1562  atmos_boundary_update_dt = parent_dtsec(handle)
1563 
1564  if ( nestqa > bnd_qa ) then
1565  log_error("ATMOS_BOUNDARY_initialize_online",*) 'NEST_BND_QA exceeds BND_QA'
1566  log_error_cont(*) 'This must not be occur.'
1567  log_error_cont(*) 'Please send your configuration file to SCALE develop member.'
1568  call prc_abort
1569  end if
1570 
1571  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
1572 
1573  return
1574  end subroutine atmos_boundary_initialize_online
1575 
1576  !-----------------------------------------------------------------------------
1578  subroutine atmos_boundary_set_online
1579  use scale_prc, only: &
1580  prc_abort
1581  use scale_const, only: &
1582  pi => const_pi
1583  use scale_time, only: &
1584  time_dtsec, &
1585  time_nstep
1586  use scale_comm_cartesc_nest, only: &
1587  online_use_velz, &
1588  parent_nstep
1589  implicit none
1590 
1591  ! parameters
1592  integer, parameter :: handle = 2
1593 
1594  real(RP) :: total
1595 
1596  ! works
1597  integer :: i, j, k, n, iq
1598  !---------------------------------------------------------------------------
1599 
1600  ! import data from parent domain
1601  boundary_timestep = 1
1602  log_info("ATMOS_BOUNDARY_set_online",*) 'BOUNDARY TIMESTEP NUMBER FOR INIT:', boundary_timestep
1603 
1604  call atmos_boundary_update_online_daughter( ref_now )
1605 
1606  if ( atmos_boundary_dens_adjust ) then
1607  call calc_mass( ref_now )
1608  end if
1609 
1610  boundary_timestep = boundary_timestep + 1
1611  log_info("ATMOS_BOUNDARY_set_online",*) 'BOUNDARY TIMESTEP NUMBER FOR INIT:', boundary_timestep
1612 
1613  call atmos_boundary_update_online_daughter( ref_new )
1614 
1615  if ( atmos_boundary_dens_adjust ) then
1616  call calc_mass( ref_new )
1617  call set_offset
1618  end if
1619 
1620  ! copy now to old
1621  !$omp parallel do collapse(2)
1622  do j = 1, ja
1623  do i = 1, ia
1624  do k = ks, ke
1625  atmos_boundary_ref_dens(k,i,j,ref_old) = atmos_boundary_ref_dens(k,i,j,ref_now)
1626  atmos_boundary_ref_velx(k,i,j,ref_old) = atmos_boundary_ref_velx(k,i,j,ref_now)
1627  atmos_boundary_ref_vely(k,i,j,ref_old) = atmos_boundary_ref_vely(k,i,j,ref_now)
1628  atmos_boundary_ref_pott(k,i,j,ref_old) = atmos_boundary_ref_pott(k,i,j,ref_now)
1629  do iq = 1, bnd_qa
1630  atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now)
1631  end do
1632  end do
1633  end do
1634  end do
1635 
1636  update_nstep = nint( atmos_boundary_update_dt / time_dtsec )
1637  if ( update_nstep * time_dtsec /= atmos_boundary_update_dt ) then
1638  log_error("ATMOS_BOUNDARY_set_online",*) 'DT of the parent is not multiple of the DT'
1639  call prc_abort
1640  end if
1641  if ( update_nstep * parent_nstep(handle) /= time_nstep ) then
1642  log_error("ATMOS_BOUNDARY_set_online",*) 'DURATION must be the same as that of the parent'
1643  call prc_abort
1644  end if
1645 
1646  now_step = 0 ! should be set as zero in initialize process
1647 
1648  if ( atmos_boundary_dens_adjust ) then
1649  allocate( offset_time_fact(update_nstep) )
1650  total = 0.0_rp
1651  !$omp parallel do reduction(+:total)
1652  do n = 0, update_nstep
1653  offset_time_fact(n) = 1.0_rp - cos( 2.0_rp * pi * ( n - 1 ) / update_nstep )
1654  total = total + offset_time_fact(n)
1655  end do
1656  total = total / update_nstep
1657  !$omp parallel do
1658  do n = 0, update_nstep
1659  offset_time_fact(n) = offset_time_fact(n) / total
1660  end do
1661  end if
1662 
1663  ! set boundary data
1665 
1666  return
1667  end subroutine atmos_boundary_set_online
1668 
1669  !-----------------------------------------------------------------------------
1671  subroutine atmos_boundary_firstsend( &
1672  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, QV, Qe )
1673  use scale_atmos_hydrometeor, only: &
1674  n_hyd
1675  use mod_atmos_phy_mp_vars, only: &
1676  qa_mp
1677  implicit none
1678 
1679  ! arguments
1680  real(RP), intent(in) :: DENS(KA,IA,JA)
1681  real(RP), intent(in) :: MOMZ(KA,IA,JA)
1682  real(RP), intent(in) :: MOMX(KA,IA,JA)
1683  real(RP), intent(in) :: MOMY(KA,IA,JA)
1684  real(RP), intent(in) :: RHOT(KA,IA,JA)
1685  real(RP), intent(in) :: QTRC(KA,IA,JA,QA_MP)
1686  real(RP), intent(in) :: QV (KA,IA,JA)
1687  real(RP), intent(in) :: Qe (KA,IA,JA,N_HYD)
1688  !---------------------------------------------------------------------------
1689 
1690  ! send data at the first time
1691  if ( do_parent_process ) then !online [parent]
1692  ! issue send
1693  call atmos_boundary_send( dens, momz, momx, momy, rhot, qtrc, qv, qe )
1694  endif
1695 
1696  return
1697  end subroutine atmos_boundary_firstsend
1698 
1699  !-----------------------------------------------------------------------------
1705  nestqa => comm_cartesc_nest_bnd_qa
1706  use scale_file_cartesc, only: &
1708  implicit none
1709 
1710  ! works
1711  integer :: handle
1712  !---------------------------------------------------------------------------
1713 
1714  if ( do_parent_process ) then !online [parent]
1715  handle = 1
1716  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
1717  endif
1718 
1719  if ( do_daughter_process ) then !online [daughter]
1720  handle = 2
1721  call comm_cartesc_nest_recv_cancel( handle )
1722  endif
1723 
1724  if ( atmos_boundary_fid > 0 ) then
1725  call file_cartesc_close( atmos_boundary_fid )
1726  atmos_boundary_fid = -1
1727  end if
1728 
1729  return
1730  end subroutine atmos_boundary_driver_finalize
1731 
1732  !-----------------------------------------------------------------------------
1734  subroutine atmos_boundary_driver_update( &
1735  last_step )
1736  use scale_prc, only: &
1737  prc_abort
1738  use scale_comm_cartesc_nest, only: &
1739  online_use_velz, &
1741  implicit none
1742 
1743  logical, intent(in) :: last_step
1744 
1745  integer :: handle
1746  !---------------------------------------------------------------------------
1747 
1748  if ( l_bnd ) then
1749 
1750  ! step boundary
1751  now_step = now_step + 1
1752 
1753  ! update referce vars
1754  if ( last_step ) then
1755  now_step = min( now_step, update_nstep-1 )
1756  else if ( now_step >= update_nstep ) then
1757  now_step = 0
1758  boundary_timestep = boundary_timestep + 1
1759 
1760  call update_ref_index
1761 
1762  if ( do_daughter_process ) then !online [daughter]
1763  call atmos_boundary_update_online_daughter( ref_new )
1764  else
1765  call atmos_boundary_update_file( ref_new )
1766  end if
1767 
1768  if ( atmos_boundary_dens_adjust ) then
1769  call calc_mass( ref_new )
1770  end if
1771 
1772  end if
1773 
1775 
1776  if ( atmos_boundary_dens_adjust ) then
1777  call set_offset
1778  end if
1779 
1780  elseif ( do_parent_process ) then
1781  ! do nothing
1782  else
1783  log_error("ATMOS_BOUNDARY_update",*) '[BUG] invalid path'
1784  call prc_abort
1785  end if
1786 
1787  call history_bnd( atmos_boundary_dens, &
1793 
1794  ! To be enable to do asynchronous communicaton
1795  if ( do_daughter_process ) then !online [daughter]
1796  handle = 2
1797  call comm_cartesc_nest_test( handle )
1798  endif
1799 
1800  return
1801  end subroutine atmos_boundary_driver_update
1802 
1803  !-----------------------------------------------------------------------------
1805  subroutine atmos_boundary_driver_send
1808  use mod_atmos_vars, only: &
1809  dens, &
1810  momz, &
1811  momx, &
1812  momy, &
1813  rhot, &
1814  qtrc, &
1815  qv, &
1816  qe
1817  use mod_atmos_phy_mp_vars, only: &
1818  qs_mp, &
1819  qe_mp
1820  implicit none
1821 
1822  integer :: handle
1823  !---------------------------------------------------------------------------
1824 
1825  if ( do_parent_process ) then !online [parent]
1826  ! should be called every time step
1827  call atmos_boundary_update_online_parent( dens,momz,momx,momy,rhot,qtrc(:,:,:,qs_mp:qe_mp), qv, qe )
1828 
1829  ! To be enable to do asynchronous communicaton
1830  handle = 1
1831  call comm_cartesc_nest_test( handle )
1832  endif
1833 
1834  return
1835  end subroutine atmos_boundary_driver_send
1836 
1837  !-----------------------------------------------------------------------------
1839  subroutine atmos_boundary_update_file( ref )
1840  use scale_file_cartesc, only: &
1841  file_cartesc_read, &
1843  use mod_atmos_phy_mp_vars, only: &
1844  qs_mp
1845  implicit none
1846 
1847  integer, intent(in) :: ref
1848 
1849  integer :: fid, iq, iqb
1850  !---------------------------------------------------------------------------
1851 
1852  log_info("ATMOS_BOUNDARY_update_file",*) "Atmos Boundary: read from boundary file(timestep=", boundary_timestep, ")"
1853 
1854  fid = atmos_boundary_fid
1855 
1856  call file_cartesc_read( fid, 'DENS', 'ZXY', atmos_boundary_ref_dens(:,:,:,ref), step=boundary_timestep )
1857  call file_cartesc_read( fid, 'VELX', 'ZXHY', atmos_boundary_ref_velx(:,:,:,ref), step=boundary_timestep )
1858  call file_cartesc_read( fid, 'VELY', 'ZXYH', atmos_boundary_ref_vely(:,:,:,ref), step=boundary_timestep )
1859  call file_cartesc_read( fid, 'PT', 'ZXY', atmos_boundary_ref_pott(:,:,:,ref), step=boundary_timestep )
1860  do iq = 1, qa
1861  iqb = bnd_iq(iq)
1862  if ( iqb > 0 ) then
1863  call file_cartesc_read( fid, tracer_name(iq), 'ZXY', atmos_boundary_ref_qtrc(:,:,:,iqb,ref), step=boundary_timestep )
1864  end if
1865  end do
1866 
1867  call file_cartesc_flush( fid )
1868 
1869  ! fill HALO in reference
1870  call atmos_boundary_ref_fillhalo( ref )
1871 
1872  return
1873  end subroutine atmos_boundary_update_file
1874 
1875  !-----------------------------------------------------------------------------
1877  subroutine atmos_boundary_update_online_parent( &
1878  DENS, & ! [in]
1879  MOMZ, & ! [in]
1880  MOMX, & ! [in]
1881  MOMY, & ! [in]
1882  RHOT, & ! [in]
1883  QTRC, & ! [in]
1884  QV, & ! [in]
1885  Qe ) ! [in]
1886  use scale_comm_cartesc_nest, only: &
1888  nestqa => comm_cartesc_nest_bnd_qa
1889  use scale_atmos_hydrometeor, only: &
1890  n_hyd
1891  use mod_atmos_phy_mp_vars, only: &
1892  qa_mp
1893  implicit none
1894 
1895  ! arguments
1896  real(rp), intent(in) :: dens(ka,ia,ja)
1897  real(rp), intent(in) :: momz(ka,ia,ja)
1898  real(rp), intent(in) :: momx(ka,ia,ja)
1899  real(rp), intent(in) :: momy(ka,ia,ja)
1900  real(rp), intent(in) :: rhot(ka,ia,ja)
1901  real(rp), intent(in) :: qtrc(ka,ia,ja,qa_mp)
1902  real(rp), intent(in) :: qv (ka,ia,ja)
1903  real(rp), intent(in) :: qe (ka,ia,ja,n_hyd)
1904 
1905  integer, parameter :: handle = 1
1906  !---------------------------------------------------------------------------
1907 
1908  log_info("ATMOS_BOUNDARY_update_online_parent",*)"ATMOS BOUNDARY update online: PARENT"
1909 
1910  ! issue wait
1911  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
1912 
1913  ! issue send
1914  call atmos_boundary_send( dens, momz, momx, momy, rhot, qtrc, qv, qe )
1915 
1916  return
1917  end subroutine atmos_boundary_update_online_parent
1918 
1919  !-----------------------------------------------------------------------------
1921  subroutine atmos_boundary_update_online_daughter( &
1922  ref ) ! [in]
1923  use scale_comm_cartesc_nest, only: &
1925  nestqa => comm_cartesc_nest_bnd_qa
1926  implicit none
1927 
1928  ! arguments
1929  integer, intent(in) :: ref
1930 
1931  integer, parameter :: handle = 2
1932  !---------------------------------------------------------------------------
1933 
1934  log_info("ATMOS_BOUNDARY_update_online_daughter",'(1X,A,I5)') 'ATMOS BOUNDARY update online: DAUGHTER', boundary_timestep
1935 
1936  ! issue wait
1937  call atmos_boundary_recv( ref )
1938 
1939  ! fill HALO in reference
1940  call atmos_boundary_ref_fillhalo( ref )
1941 
1942  ! issue receive
1943  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
1944 
1945  return
1946  end subroutine atmos_boundary_update_online_daughter
1947 
1948  !-----------------------------------------------------------------------------
1950  subroutine atmos_boundary_send( &
1951  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, QV, Qe )
1952  use scale_comm_cartesc_nest, only: &
1955  daughter_ka, &
1956  daughter_ia, &
1957  daughter_ja, &
1958  nestqa => comm_cartesc_nest_bnd_qa
1959  use scale_atmos_hydrometeor, only: &
1960  n_hyd
1961  use mod_atmos_phy_mp_vars, only: &
1962  qa_mp
1963  implicit none
1964 
1965  ! parameters
1966  integer, parameter :: handle = 1
1967 
1968  ! arguments
1969  real(rp), intent(in) :: dens(ka,ia,ja)
1970  real(rp), intent(in) :: momz(ka,ia,ja)
1971  real(rp), intent(in) :: momx(ka,ia,ja)
1972  real(rp), intent(in) :: momy(ka,ia,ja)
1973  real(rp), intent(in) :: rhot(ka,ia,ja)
1974  real(rp), intent(in), target :: qtrc(ka,ia,ja,qa_mp)
1975  real(rp), intent(in) :: qv (ka,ia,ja)
1976  real(rp), intent(in) :: qe (ka,ia,ja,n_hyd)
1977 
1978  ! works
1979  real(rp), pointer :: q(:,:,:,:)
1980  real(rp) :: dummy_d( daughter_ka(handle), daughter_ia(handle), daughter_ja(handle), nestqa )
1981 
1982  integer :: iq
1983  !---------------------------------------------------------------------------
1984 
1985  if ( online_boundary_diagqhyd ) then
1986  q => q_work
1987  q(:,:,:,1) = qv(:,:,:)
1988  do iq = 2, nestqa
1989  q(:,:,:,iq) = qe(:,:,:,iq-1)
1990  end do
1991  else
1992  q => qtrc
1993  end if
1994 
1995  call comm_cartesc_nest_nestdown( handle, &
1996  nestqa, &
1997  dens(:,:,:), & !(KA,IA,JA)
1998  momz(:,:,:), & !(KA,IA,JA)
1999  momx(:,:,:), & !(KA,IA,JA)
2000  momy(:,:,:), & !(KA,IA,JA)
2001  rhot(:,:,:), & !(KA,IA,JA)
2002  q(:,:,:,:), & !(KA,IA,JA,NESTQA)
2003  dummy_d(:,:,:,1), & !(KA,IA,JA)
2004  dummy_d(:,:,:,1), & !(KA,IA,JA)
2005  dummy_d(:,:,:,1), & !(KA,IA,JA)
2006  dummy_d(:,:,:,1), & !(KA,IA,JA)
2007  dummy_d(:,:,:,1), & !(KA,IA,JA)
2008  dummy_d(:,:,:,:) ) !(KA,IA,JA,NESTQA)
2009 
2010  return
2011  end subroutine atmos_boundary_send
2012 
2013  !-----------------------------------------------------------------------------
2015  subroutine atmos_boundary_recv( &
2016  ref_idx )
2017  use scale_prc, only: &
2018  prc_abort
2019  use scale_comm_cartesc_nest, only: &
2022  parent_ka, &
2023  parent_ia, &
2024  parent_ja, &
2025  nestqa => comm_cartesc_nest_bnd_qa
2026  use mod_atmos_phy_mp_driver, only: &
2028  implicit none
2029 
2030  ! parameters
2031  integer, parameter :: handle = 2
2032 
2033  ! arguments
2034  integer, intent(in) :: ref_idx
2035 
2036  ! works
2037  real(rp), pointer :: q(:,:,:,:)
2038  real(rp) :: dummy_p( parent_ka(handle), parent_ia(handle), parent_ja(handle), nestqa )
2039  !---------------------------------------------------------------------------
2040 
2041 !OCL XFILL
2042  dummy_p(:,:,:,:) = 0.0_rp
2043 
2044  if ( online_boundary_diagqhyd ) then
2045  q => q_work
2046  else
2047  q => atmos_boundary_ref_qtrc(:,:,:,1:nestqa,ref_idx)
2048  end if
2049 
2050  call comm_cartesc_nest_nestdown( handle, &
2051  nestqa, &
2052  dummy_p(:,:,:,1), & !(KA,IA,JA)
2053  dummy_p(:,:,:,1), & !(KA,IA,JA)
2054  dummy_p(:,:,:,1), & !(KA,IA,JA)
2055  dummy_p(:,:,:,1), & !(KA,IA,JA)
2056  dummy_p(:,:,:,1), & !(KA,IA,JA)
2057  dummy_p(:,:,:,1:nestqa), & !(KA,IA,JA,NESTQA)
2058  atmos_boundary_ref_dens(:,:,:,ref_idx), & !(KA,IA,JA)
2059  atmos_boundary_ref_velz(:,:,:,ref_idx), & !(KA,IA,JA)
2060  atmos_boundary_ref_velx(:,:,:,ref_idx), & !(KA,IA,JA)
2061  atmos_boundary_ref_vely(:,:,:,ref_idx), & !(KA,IA,JA)
2062  atmos_boundary_ref_pott(:,:,:,ref_idx), & !(KA,IA,JA)
2063  q(:,:,:,:) ) !(KA,IA,JA,NESTQA)
2064 
2065 
2066  if ( online_boundary_diagqhyd ) then
2067  call atmos_phy_mp_driver_qhyd2qtrc( ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2068  q(:,:,:,1), q(:,:,:,2:), & ! (in)
2069  atmos_boundary_ref_qtrc(:,:,:,:,ref_idx) ) ! (out)
2070  end if
2071 
2072  return
2073  end subroutine atmos_boundary_recv
2074 
2075  !-----------------------------------------------------------------------------
2076  subroutine set_boundary( use_velz )
2077  use scale_prc_cartesc, only: &
2078  prc_has_w, &
2079  prc_has_e, &
2080  prc_has_s, &
2081  prc_has_n
2082  use mod_atmos_vars, only: &
2083  dens, &
2084  momz, &
2085  momx, &
2086  momy, &
2087  rhot, &
2088  qtrc
2089  implicit none
2090  logical, intent(in) :: use_velz
2091 
2092  real(RP) :: bnd_DENS(KA,IA,JA) ! damping coefficient for DENS (0-1)
2093  real(RP) :: bnd_VELZ(KA,IA,JA) ! damping coefficient for VELZ (0-1)
2094  real(RP) :: bnd_VELX(KA,IA,JA) ! damping coefficient for VELX (0-1)
2095  real(RP) :: bnd_VELY(KA,IA,JA) ! damping coefficient for VELY (0-1)
2096  real(RP) :: bnd_POTT(KA,IA,JA) ! damping coefficient for POTT (0-1)
2097  real(RP) :: bnd_QTRC(KA,IA,JA,BND_QA) ! damping coefficient for QTRC (0-1)
2098 
2099  integer :: i, j, k, iq, iqb
2100 
2101 
2102  ! get boundaryal coefficients
2103  call get_boundary( bnd_dens(:,:,:), & ! [OUT]
2104  bnd_velz(:,:,:), & ! [OUT]
2105  bnd_velx(:,:,:), & ! [OUT]
2106  bnd_vely(:,:,:), & ! [OUT]
2107  bnd_pott(:,:,:), & ! [OUT]
2108  bnd_qtrc(:,:,:,:), & ! [OUT]
2109  now_step, & ! [IN]
2110  update_nstep ) ! [IN]
2111 
2112  ! update boundary vars
2113  !$omp parallel do collapse(2)
2114  do j = 1, ja
2115  do i = 1, ia
2116  do k = ks, ke
2117  atmos_boundary_dens(k,i,j) = bnd_dens(k,i,j)
2118  atmos_boundary_velx(k,i,j) = bnd_velx(k,i,j)
2119  atmos_boundary_vely(k,i,j) = bnd_vely(k,i,j)
2120  atmos_boundary_pott(k,i,j) = bnd_pott(k,i,j)
2121  do iq = 1, bnd_qa
2122  atmos_boundary_qtrc(k,i,j,iq) = bnd_qtrc(k,i,j,iq)
2123  end do
2124  end do
2125  end do
2126  end do
2127  if ( use_velz ) then
2128  !$omp parallel do collapse(2)
2129  do j = 1, ja
2130  do i = 1, ia
2131  do k = ks, ke
2132  atmos_boundary_velz(k,i,j) = bnd_velz(k,i,j)
2133  end do
2134  end do
2135  end do
2136  end if
2137 
2138  ! fill HALO in western region
2139  if ( .NOT. prc_has_w ) then
2140  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2141  !$omp shared(JA,JS,IS,KS,KE,QA,DENS,MOMX,RHOT,QTRC) &
2142  !$omp shared(ATMOS_BOUNDARY_DENS,ATMOS_BOUNDARY_VELX) &
2143  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC) &
2144  !$omp shared(BND_QA,BND_IQ) &
2145  !$omp private(i,j,k,iq,iqb)
2146  do j = 1, ja
2147  do i = 1, is-1
2148  do k = ks, ke
2149  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2150  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2151  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2152  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2153  do iq = 1, qa
2154  iqb = bnd_iq(iq)
2155  if ( iqb > 0 ) then
2156  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2157  else
2158  qtrc(k,i,j,iq) = qtrc(k,is,j,iq)
2159  end if
2160  end do
2161  end do
2162  end do
2163  end do
2164  !$omp parallel do
2165  do j = 1, ja-1
2166  do i = 1, is-1
2167  do k = ks, ke
2168  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2169  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2170  end do
2171  end do
2172  end do
2173  do i = 1, is-1
2174  do k = ks, ke
2175  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) &
2176  * atmos_boundary_dens(k,i,ja)
2177  end do
2178  end do
2179  if ( use_velz ) then
2180  !$omp parallel do
2181  do j = 1, ja
2182  do i = 1, is-1
2183  do k = ks, ke-1
2184  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2185  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2186  end do
2187  end do
2188  end do
2189  else
2190  !$omp parallel do
2191  do j = 1, ja
2192  do i = 1, is-1
2193  do k = ks, ke-1
2194  momz(k,i,j) = momz(k,is,j)
2195  end do
2196  end do
2197  end do
2198  end if
2199  end if
2200 
2201  ! fill HALO in eastern region
2202  if ( .NOT. prc_has_e ) then
2203  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2204  !$omp shared(JA,IE,IA,KS,KE,QA) &
2205  !$omp shared(DENS,RHOT,QTRC) &
2206  !$omp shared(ATMOS_BOUNDARY_DENS,ATMOS_BOUNDARY_VELX) &
2207  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC) &
2208  !$omp shared(BND_QA,BND_IQ) &
2209  !$omp private(i,j,k,iq,iqb)
2210  do j = 1, ja
2211  do i = ie+1, ia
2212  do k = ks, ke
2213  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2214  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2215  do iq = 1, qa
2216  iqb = bnd_iq(iq)
2217  if ( iqb > 0 ) then
2218  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2219  else
2220  qtrc(k,i,j,iq) = qtrc(k,ie,j,iq)
2221  end if
2222  end do
2223  end do
2224  end do
2225  end do
2226  !$omp parallel do
2227  do j = 1, ja
2228  do i = ie, ia-1
2229  do k = ks, ke
2230  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2231  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2232  end do
2233  end do
2234  end do
2235  !$omp parallel do
2236  do j = 1, ja
2237  do k = ks, ke
2238  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) * atmos_boundary_dens(k,ia,j)
2239  end do
2240  end do
2241  !$omp parallel do
2242  do j = 1, ja-1
2243  do i = ie+1, ia
2244  do k = ks, ke
2245  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2246  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2247  end do
2248  end do
2249  end do
2250  do i = ie+1, ia
2251  do k = ks, ke
2252  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) &
2253  * atmos_boundary_dens(k,i,ja)
2254  end do
2255  end do
2256  if ( use_velz ) then
2257  !$omp parallel do
2258  do j = 1, ja
2259  do i = ie+1, ia
2260  do k = ks, ke-1
2261  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2262  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2263  end do
2264  end do
2265  end do
2266  else
2267  !$omp parallel do
2268  do j = 1, ja
2269  do i = ie+1, ia
2270  do k = ks, ke-1
2271  momz(k,i,j) = momz(k,ie,j)
2272  end do
2273  end do
2274  end do
2275  end if
2276  end if
2277 
2278  ! fill HALO in southern region
2279  if ( .NOT. prc_has_s ) then
2280  do j = 1, js-1
2281  do i = 1, ia
2282  do k = ks, ke
2283  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2284  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2285  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2286  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2287  do iq = 1, qa
2288  iqb = bnd_iq(iq)
2289  if ( iqb > 0 ) then
2290  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2291  else
2292  qtrc(k,i,j,iq) = qtrc(k,i,js,iq)
2293  end if
2294  end do
2295  end do
2296  end do
2297  end do
2298  do j = 1, js-1
2299  do i = 1, ia-1
2300  do k = ks, ke
2301  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2302  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2303  end do
2304  end do
2305  end do
2306  do j = 1, js-1
2307  do k = ks, ke
2308  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) &
2309  * atmos_boundary_dens(k,ia,j)
2310  end do
2311  end do
2312  if ( use_velz ) then
2313  do j = 1, js-1
2314  do i = 1, ia
2315  do k = ks, ke-1
2316  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2317  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2318  end do
2319  end do
2320  end do
2321  else
2322  do j = 1, js-1
2323  do i = 1, ia
2324  do k = ks, ke-1
2325  momz(k,i,j) = momz(k,i,js)
2326  end do
2327  end do
2328  end do
2329  end if
2330  end if
2331 
2332  ! fill HALO in northern region
2333  if ( .NOT. prc_has_n ) then
2334  do j = je+1, ja
2335  do i = 1, ia
2336  do k = ks, ke
2337  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2338  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2339  do iq = 1, qa
2340  iqb = bnd_iq(iq)
2341  if ( iqb > 0 ) then
2342  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2343  else
2344  qtrc(k,i,j,iq) = qtrc(k,i,je,iq)
2345  end if
2346  end do
2347  end do
2348  end do
2349  end do
2350  do j = je, ja-1
2351  do i = 1, ia
2352  do k = ks, ke
2353  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2354  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2355  end do
2356  end do
2357  end do
2358  do i = 1, ia
2359  do k = ks, ke
2360  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) * atmos_boundary_dens(k,i,ja)
2361  end do
2362  end do
2363  do j = je+1, ja
2364  do i = 1, ia-1
2365  do k = ks, ke
2366  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2367  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2368  end do
2369  end do
2370  end do
2371  do j = je+1, ja
2372  do k = ks, ke
2373  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) &
2374  * atmos_boundary_dens(k,ia,j)
2375  end do
2376  end do
2377  if ( use_velz ) then
2378  do j = je+1, ja
2379  do i = 1, ia
2380  do k = ks, ke-1
2381  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2382  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2383  end do
2384  end do
2385  end do
2386  else
2387  do j = je+1, ja
2388  do i = 1, ia
2389  do k = ks, ke-1
2390  momz(k,i,j) = momz(k,i,je)
2391  end do
2392  end do
2393  end do
2394  end if
2395  end if
2396 
2397  return
2398  end subroutine set_boundary
2399  !-----------------------------------------------------------------------------
2401  subroutine get_boundary_same_parent( &
2402  bnd_DENS, &
2403  bnd_VELZ, &
2404  bnd_VELX, &
2405  bnd_VELY, &
2406  bnd_POTT, &
2407  bnd_QTRC, &
2408  now_step, &
2409  update_step )
2410  implicit none
2411 
2412  ! arguments
2413  real(RP), intent(out) :: bnd_DENS(:,:,:)
2414  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2415  real(RP), intent(out) :: bnd_VELX(:,:,:)
2416  real(RP), intent(out) :: bnd_VELY(:,:,:)
2417  real(RP), intent(out) :: bnd_POTT(:,:,:)
2418  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2419  integer, intent(in) :: now_step
2420  integer, intent(in) :: update_step
2421 
2422  ! works
2423  integer :: i, j, k, iq
2424  integer :: ref
2425  !---------------------------------------------------------------------------
2426 
2427  if ( now_step == update_step ) then
2428  ref = ref_new
2429  else
2430  ref = ref_now
2431  end if
2432 
2433  !$omp parallel do collapse(2)
2434  do j = 1, ja
2435  do i = 1, ia
2436  do k = ks, ke
2437  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref)
2438  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref)
2439  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref)
2440  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref)
2441  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref)
2442  do iq = 1, bnd_qa
2443  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref)
2444  end do
2445  end do
2446  end do
2447  end do
2448 
2449  return
2450  end subroutine get_boundary_same_parent
2451 
2452  !-----------------------------------------------------------------------------
2454  subroutine get_boundary_nearest_neighbor( &
2455  bnd_DENS, &
2456  bnd_VELZ, &
2457  bnd_VELX, &
2458  bnd_VELY, &
2459  bnd_POTT, &
2460  bnd_QTRC, &
2461  now_step, &
2462  update_step )
2463  implicit none
2464 
2465  ! parameters
2466  real(RP) :: EPS = 1.0e-4_rp
2467 
2468  ! arguments
2469  real(RP), intent(out) :: bnd_DENS(:,:,:)
2470  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2471  real(RP), intent(out) :: bnd_VELX(:,:,:)
2472  real(RP), intent(out) :: bnd_VELY(:,:,:)
2473  real(RP), intent(out) :: bnd_POTT(:,:,:)
2474  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2475  integer, intent(in) :: now_step
2476  integer, intent(in) :: update_step
2477 
2478  ! works
2479  integer :: i, j, k, iq
2480  integer :: ref_idx
2481 
2482  real(RP) :: real_nstep
2483  real(RP) :: half_nstep
2484  !---------------------------------------------------------------------------
2485 
2486  real_nstep = real( now_step, kind=rp )
2487  half_nstep = real( update_nstep, kind=rp ) * 0.5_rp
2488 
2489  ! this step before half of the parent step
2490  if( ( real_nstep - eps ) < half_nstep ) then
2491  ref_idx = ref_now
2492 
2493  ! this step after half of the parent step
2494  else if( ( real_nstep - 1.0_rp + eps ) > half_nstep ) then
2495  ref_idx = ref_new
2496 
2497  ! this step across half of the parent step
2498  else
2499  ref_idx = ref_now
2500 
2501  end if
2502 
2503  !$omp parallel do collapse(2)
2504  do j = 1, ja
2505  do i = 1, ia
2506  do k = ks, ke
2507  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_idx)
2508  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_idx)
2509  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_idx)
2510  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_idx)
2511  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_idx)
2512  do iq = 1, bnd_qa
2513  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_idx)
2514  end do
2515  end do
2516  end do
2517  end do
2518 
2519  return
2520  end subroutine get_boundary_nearest_neighbor
2521 
2522  !-----------------------------------------------------------------------------
2524  subroutine get_boundary_lerp_initpoint( &
2525  bnd_DENS, &
2526  bnd_VELZ, &
2527  bnd_VELX, &
2528  bnd_VELY, &
2529  bnd_POTT, &
2530  bnd_QTRC, &
2531  now_step, &
2532  update_step )
2533  implicit none
2534 
2535  ! arguments
2536  real(RP), intent(out) :: bnd_DENS(:,:,:)
2537  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2538  real(RP), intent(out) :: bnd_VELX(:,:,:)
2539  real(RP), intent(out) :: bnd_VELY(:,:,:)
2540  real(RP), intent(out) :: bnd_POTT(:,:,:)
2541  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2542  integer, intent(in) :: now_step
2543  integer, intent(in) :: update_step
2544 
2545  ! works
2546  integer :: i, j, k, iq
2547 
2548  real(RP) :: fact
2549  !---------------------------------------------------------------------------
2550 
2551  fact = ( now_step + 0.5_rp ) / update_step
2552 
2553  !$omp parallel do default(none) private(i,j,k,iq) OMP_SCHEDULE_ collapse(2) &
2554  !$omp shared(JA,IA,KS,KE,bnd_DENS,ATMOS_BOUNDARY_ref_DENS,ref_now,fact,ref_new,bnd_VELZ) &
2555  !$omp shared(ATMOS_BOUNDARY_ref_VELZ,bnd_VELX,ATMOS_BOUNDARY_ref_VELX,bnd_VELY) &
2556  !$omp shared(ATMOS_BOUNDARY_ref_VELY,bnd_POTT,ATMOS_BOUNDARY_ref_POTT,BND_QA,bnd_QTRC,ATMOS_BOUNDARY_ref_QTRC)
2557  do j = 1, ja
2558  do i = 1, ia
2559  do k = ks, ke
2560  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_now) * ( 1.0_rp-fact ) &
2561  + atmos_boundary_ref_dens(k,i,j,ref_new) * fact
2562  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_now) * ( 1.0_rp-fact ) &
2563  + atmos_boundary_ref_velz(k,i,j,ref_new) * fact
2564  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_now) * ( 1.0_rp-fact) &
2565  + atmos_boundary_ref_velx(k,i,j,ref_new) * fact
2566  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_now) * ( 1.0_rp-fact) &
2567  + atmos_boundary_ref_vely(k,i,j,ref_new) * fact
2568  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_now) * ( 1.0_rp-fact ) &
2569  + atmos_boundary_ref_pott(k,i,j,ref_new) * fact
2570  do iq = 1, bnd_qa
2571  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * ( 1.0_rp-fact ) &
2572  + atmos_boundary_ref_qtrc(k,i,j,iq,ref_new) * fact
2573  end do
2574  end do
2575  end do
2576  end do
2577 
2578  return
2579  end subroutine get_boundary_lerp_initpoint
2580 
2581  !-----------------------------------------------------------------------------
2583  subroutine get_boundary_lerp_midpoint( &
2584  bnd_DENS, &
2585  bnd_VELZ, &
2586  bnd_VELX, &
2587  bnd_VELY, &
2588  bnd_POTT, &
2589  bnd_QTRC, &
2590  now_step, &
2591  update_step )
2592  use scale_time, only: &
2593  time_dtsec
2594  implicit none
2595 
2596  ! parameters
2597  real(RP) :: EPS = 1.0e-4_rp
2598 
2599  ! arguments
2600  real(RP), intent(out) :: bnd_DENS(:,:,:)
2601  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2602  real(RP), intent(out) :: bnd_VELX(:,:,:)
2603  real(RP), intent(out) :: bnd_VELY(:,:,:)
2604  real(RP), intent(out) :: bnd_POTT(:,:,:)
2605  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2606  integer, intent(in) :: now_step
2607  integer, intent(in) :: update_step
2608 
2609  ! works
2610  integer :: i, j, k, iq
2611 
2612  real(RP) :: real_nstep
2613  real(RP) :: half_nstep
2614  real(RP) :: t1, t2
2615  !---------------------------------------------------------------------------
2616 
2617  real_nstep = real( now_step, kind=rp )
2618  half_nstep = real( update_nstep, kind=rp ) * 0.5_rp
2619 
2620  ! this step before half of the parent step
2621  if( ( real_nstep - eps ) < half_nstep ) then
2622 
2623  t1 = time_dtsec / atmos_boundary_update_dt * ( real_nstep + half_nstep - 0.5_rp )
2624 
2625  !$omp parallel do collapse(2)
2626  do j = 1, ja
2627  do i = 1, ia
2628  do k = ks, ke
2629  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_now) * t1 &
2630  - atmos_boundary_ref_dens(k,i,j,ref_old) * ( t1 - 1.0_rp )
2631  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_now) * t1 &
2632  - atmos_boundary_ref_velz(k,i,j,ref_old) * ( t1 - 1.0_rp )
2633  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_now) * t1 &
2634  - atmos_boundary_ref_velx(k,i,j,ref_old) * ( t1 - 1.0_rp )
2635  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_now) * t1 &
2636  - atmos_boundary_ref_vely(k,i,j,ref_old) * ( t1 - 1.0_rp )
2637  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_now) * t1 &
2638  - atmos_boundary_ref_pott(k,i,j,ref_old) * ( t1 - 1.0_rp )
2639  do iq = 1, bnd_qa
2640  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * t1 &
2641  - atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) * ( t1 - 1.0_rp )
2642  end do
2643  end do
2644  end do
2645  end do
2646 
2647  ! this step after half of the parent step
2648  else if( ( real_nstep - 1.0_rp + eps ) > half_nstep ) then
2649 
2650  t1 = time_dtsec / atmos_boundary_update_dt * ( real_nstep - half_nstep - 0.5_rp )
2651 
2652  !$omp parallel do collapse(2)
2653  do j = 1, ja
2654  do i = 1, ia
2655  do k = ks, ke
2656  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_new) * t1 &
2657  - atmos_boundary_ref_dens(k,i,j,ref_now) * ( t1 - 1.0_rp )
2658  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_new) * t1 &
2659  - atmos_boundary_ref_velz(k,i,j,ref_now) * ( t1 - 1.0_rp )
2660  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_new) * t1 &
2661  - atmos_boundary_ref_velx(k,i,j,ref_now) * ( t1 - 1.0_rp )
2662  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_new) * t1 &
2663  - atmos_boundary_ref_vely(k,i,j,ref_now) * ( t1 - 1.0_rp )
2664  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_new) * t1 &
2665  - atmos_boundary_ref_pott(k,i,j,ref_now) * ( t1 - 1.0_rp )
2666  do iq = 1, bnd_qa
2667  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_new) * t1 &
2668  - atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * ( t1 - 1.0_rp )
2669  end do
2670  end do
2671  end do
2672  end do
2673 
2674  ! this step across half of the parent step
2675  else
2676 
2677  t1 = time_dtsec / atmos_boundary_update_dt * ( real_nstep + half_nstep - 1.0_rp )
2678  t2 = time_dtsec / atmos_boundary_update_dt * ( real_nstep - half_nstep )
2679 
2680  !$omp parallel do collapse(2)
2681  do j = 1, ja
2682  do i = 1, ia
2683  do k = ks, ke
2684  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_new) * t2 * 0.25_rp &
2685  + atmos_boundary_ref_dens(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2686  - atmos_boundary_ref_dens(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2687  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_new) * t2 * 0.25_rp &
2688  + atmos_boundary_ref_velz(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2689  - atmos_boundary_ref_velz(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2690  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_new) * t2 * 0.25_rp &
2691  + atmos_boundary_ref_velx(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2692  - atmos_boundary_ref_velx(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2693  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_new) * t2 * 0.25_rp &
2694  + atmos_boundary_ref_vely(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2695  - atmos_boundary_ref_vely(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2696  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_new) * t2 * 0.25_rp &
2697  + atmos_boundary_ref_pott(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2698  - atmos_boundary_ref_pott(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2699  do iq = 1, bnd_qa
2700  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_new) * t2 * 0.25_rp &
2701  + atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2702  - atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) * ( t1 - 1.0_rp )
2703  end do
2704  end do
2705  end do
2706  end do
2707 
2708  end if
2709 
2710  return
2711  end subroutine get_boundary_lerp_midpoint
2712 
2713  !-----------------------------------------------------------------------------
2715  subroutine update_ref_index
2716  implicit none
2717 
2718  ! works
2719  integer :: ref_tmp
2720  !---------------------------------------------------------------------------
2721 
2722  ref_tmp = ref_old
2723  ref_old = ref_now
2724  ref_now = ref_new
2725  ref_new = ref_tmp
2726 
2727  return
2728  end subroutine update_ref_index
2729 
2730  subroutine history_bnd( &
2731  ATMOS_BOUNDARY_DENS, &
2732  ATMOS_BOUNDARY_VELZ, &
2733  ATMOS_BOUNDARY_VELX, &
2734  ATMOS_BOUNDARY_VELY, &
2735  ATMOS_BOUNDARY_POTT, &
2736  ATMOS_BOUNDARY_QTRC )
2737  use scale_file_history, only: &
2738  file_history_in
2739  use mod_atmos_phy_mp_vars, only: &
2740  qs_mp
2741  implicit none
2742  real(RP), intent(in) :: ATMOS_BOUNDARY_DENS(KA,IA,JA)
2743  real(RP), intent(in) :: ATMOS_BOUNDARY_VELZ(KA,IA,JA)
2744  real(RP), intent(in) :: ATMOS_BOUNDARY_VELX(KA,IA,JA)
2745  real(RP), intent(in) :: ATMOS_BOUNDARY_VELY(KA,IA,JA)
2746  real(RP), intent(in) :: ATMOS_BOUNDARY_POTT(KA,IA,JA)
2747  real(RP), intent(in) :: ATMOS_BOUNDARY_QTRC(KA,IA,JA,BND_QA)
2748 
2749  integer :: iq, iqb
2750 
2751  call file_history_in( atmos_boundary_dens(:,:,:), 'DENS_BND', 'Boundary Density', 'kg/m3' )
2752  call file_history_in( atmos_boundary_velz(:,:,:), 'VELZ_BND', 'Boundary velocity z-direction', 'm/s', dim_type='ZHXY' )
2753  call file_history_in( atmos_boundary_velx(:,:,:), 'VELX_BND', 'Boundary velocity x-direction', 'm/s', dim_type='ZXHY' )
2754  call file_history_in( atmos_boundary_vely(:,:,:), 'VELY_BND', 'Boundary velocity y-direction', 'm/s', dim_type='ZXYH' )
2755  call file_history_in( atmos_boundary_pott(:,:,:), 'PT_BND', 'Boundary potential temperature', 'K' )
2756  do iq = 1, qa
2757  iqb = bnd_iq(iq)
2758  if ( iqb > 0 ) then
2759  call file_history_in( atmos_boundary_qtrc(:,:,:,iqb), trim(tracer_name(iq))//'_BND', &
2760  trim(tracer_name(iq))//' in boundary', 'kg/kg' )
2761  end if
2762  enddo
2763 
2764  return
2765  end subroutine history_bnd
2766 
2767  subroutine calc_mass( ref )
2768  use scale_prc_cartesc, only: &
2769  prc_has_w, &
2770  prc_has_e, &
2771  prc_has_s, &
2772  prc_has_n
2773  use scale_statistics, only: &
2774  statistics_total
2775  use scale_atmos_grid_cartesc_real, only: &
2776  vol => atmos_grid_cartesc_real_vol, &
2777  totvol => atmos_grid_cartesc_real_totvol, &
2778  areazxv_y => atmos_grid_cartesc_real_areazxv_y, &
2779  totareazuy_x => atmos_grid_cartesc_real_totareazuy_x, &
2780  totareazxv_y => atmos_grid_cartesc_real_totareazxv_y
2781  use scale_atmos_refstate, only: &
2782  dens_ref => atmos_refstate_dens
2783  use mod_atmos_vars, only: &
2784  dens
2785  implicit none
2786  integer, intent(in) :: ref
2787 
2788  real(DP) :: masstot, masstot_current
2789  real(DP) :: massflx
2790  real(DP) :: offset_band, offset_bias
2791  real(DP) :: ref_tot
2792  real(DP) :: flx_w, flx_e, flx_s, flx_n
2793  real(DP) :: ref_w, ref_e, ref_s, ref_n
2794 
2795  real(RP), target :: work_x(KA,JA), work_y(KA,IA)
2796  real(RP), pointer :: ptr(:,:)
2797 
2798  integer :: k, i, j
2799 
2800  ! total mass
2801  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
2802  atmos_boundary_ref_dens(:,:,:,ref), & ! (in)
2803  "DENS_bnd", & ! (in)
2804  vol(:,:,:), totvol, & ! (in)
2805  log_suppress = .true., global = .true., & ! (in)
2806  sum = masstot ) ! (out)
2807 
2808 !!$ call STATISTICS_total( KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2809 !!$ DENS(:,:,:), "DENS_bnd_update", & ! (in)
2810 !!$ VOL(:,:,:), TOTVOL, & ! (in)
2811 !!$ log_suppress = .true., global = .true., & ! (in)
2812 !!$ sum = masstot_current ) ! (out)
2813 
2814 
2815  ! West
2816  if ( .NOT. prc_has_w ) then
2817  !$omp parallel do
2818  do j = js, je
2819  do k = ks, ke
2820  work_x(k,j) = atmos_boundary_ref_velx(k,is-1,j,ref) &
2821  * ( atmos_boundary_ref_dens(k,is-1,j,ref) + atmos_boundary_ref_dens(k,is,j,ref) ) * 0.5_rp
2822  end do
2823  end do
2824  ptr => work_x
2825  else
2826  ptr => zero_x
2827  end if
2828  call statistics_total( ka, ks, ke, ja, js, je, &
2829  ptr(:,:), "MFLUX_bnd_w", & ! (in)
2830  areazuy_w(:,:), totareazuy_x(is-1), & ! (in)
2831  log_suppress = .true., global = .true., & ! (in)
2832  sum = flx_w ) ! (out)
2833  if ( .NOT. prc_has_w ) then
2834  !$omp parallel do
2835  do j = js, je
2836  do k = ks, ke
2837  work_x(k,j) = dens_ref(k,is,j)
2838  end do
2839  end do
2840  end if
2841  call statistics_total( ka, ks, ke, ja, js, je, &
2842  ptr(:,:), "DENS_ref_w", & ! (in)
2843  areazuy_w(:,:), totareazuy_x(is-1), & ! (in)
2844  log_suppress = .true., global = .true., & ! (in)
2845  sum = ref_w ) ! (out)
2846 
2847  ! East
2848  if ( .NOT. prc_has_e ) then
2849  !$omp parallel do
2850  do j = js, je
2851  do k = ks, ke
2852  work_x(k,j) = atmos_boundary_ref_velx(k,ie,j,ref) &
2853  * ( atmos_boundary_ref_dens(k,ie,j,ref) + atmos_boundary_ref_dens(k,ie+1,j,ref) ) * 0.5_rp
2854  end do
2855  end do
2856  ptr => work_x
2857  else
2858  ptr => zero_x
2859  end if
2860  call statistics_total( ka, ks, ke, ja, js, je, &
2861  ptr(:,:), "MFLUX_bnd_e", & ! (in)
2862  areazuy_e(:,:), totareazuy_x(ie), & ! (in)
2863  log_suppress = .true., global = .true., & ! (in)
2864  sum = flx_e ) ! (out)
2865  if ( .NOT. prc_has_e ) then
2866  !$omp parallel do
2867  do j = js, je
2868  do k = ks, ke
2869  work_x(k,j) = dens_ref(k,ie,j)
2870  end do
2871  end do
2872  end if
2873  call statistics_total( ka, ks, ke, ja, js, je, &
2874  ptr(:,:), "DENS_ref_e", & ! (in)
2875  areazuy_e(:,:), totareazuy_x(ie), & ! (in)
2876  log_suppress = .true., global = .true., & ! (in)
2877  sum = ref_e ) ! (out)
2878 
2879  ! South
2880  if ( .NOT. prc_has_s ) then
2881  !$omp parallel do
2882  do i = is, ie
2883  do k = ks, ke
2884  work_y(k,i) = atmos_boundary_ref_vely(k,i,js-1,ref) &
2885  * ( atmos_boundary_ref_dens(k,i,js-1,ref) + atmos_boundary_ref_dens(k,i,js,ref) ) * 0.5_rp
2886  end do
2887  end do
2888  ptr => work_y
2889  else
2890  ptr => zero_y
2891  end if
2892  call statistics_total( ka, ks, ke, ia, is, ie, &
2893  ptr(:,:), "MFLUX_bnd_s", & ! (in)
2894  areazxv_y(:,:,js-1), totareazxv_y(js-1), & ! (in)
2895  log_suppress = .true., global = .true., & ! (in)
2896  sum = flx_s ) ! (out)
2897  if ( .NOT. prc_has_s ) then
2898  !$omp parallel do
2899  do i = is, ie
2900  do k = ks, ke
2901  work_y(k,i) = dens_ref(k,i,js)
2902  end do
2903  end do
2904  end if
2905  call statistics_total( ka, ks, ke, ia, is, ie, &
2906  ptr(:,:), "DENS_ref_s", & ! (in)
2907  areazxv_y(:,:,js-1), totareazxv_y(js-1), & ! (in)
2908  log_suppress = .true., global = .true., & ! (in)
2909  sum = ref_s ) ! (out)
2910 
2911  ! North
2912  if ( .NOT. prc_has_n ) then
2913  !$omp parallel do
2914  do i = is, ie
2915  do k = ks, ke
2916  work_y(k,i) = atmos_boundary_ref_vely(k,i,je,ref) &
2917  * ( atmos_boundary_ref_dens(k,i,je,ref) + atmos_boundary_ref_dens(k,i,je+1,ref) ) * 0.5_rp
2918  end do
2919  end do
2920  ptr => work_y
2921  else
2922  ptr => zero_y
2923  end if
2924  call statistics_total( ka, ks, ke, ia, is, ie, &
2925  ptr(:,:), "MFLX_bnd_n", & ! (in)
2926  areazxv_y(:,:,je), totareazxv_y(je), & ! (in)
2927  log_suppress = .true., global = .true., & ! (in)
2928  sum = flx_n ) ! (out)
2929  if ( .NOT. prc_has_n ) then
2930  !$omp parallel do
2931  do i = is, ie
2932  do k = ks, ke
2933  work_y(k,i) = dens_ref(k,i,je)
2934  end do
2935  end do
2936  end if
2937  call statistics_total( ka, ks, ke, ia, is, ie, &
2938  ptr(:,:), "DENS_ref_n", & ! (in)
2939  areazxv_y(:,:,je), totareazxv_y(je), & ! (in)
2940  log_suppress = .true., global = .true., & ! (in)
2941  sum = ref_n ) ! (out)
2942 
2943  massflx = flx_w - flx_e + flx_s - flx_n
2944 
2945  offset_band = ( masstot - masstot_now ) / atmos_boundary_update_dt &
2946  - ( massflx + massflx_now ) * 0.5_dp
2947 ! offset_bias = ( MASSTOT_now - masstot_current ) / ATMOS_BOUNDARY_UPDATE_DT
2948 
2949  log_info("ATMOS_BOUNDARY_calc_mass",*) "Offset_band is: ", offset_band, "(", masstot, masstot_now, massflx, massflx_now, ")"
2950 
2951  ref_tot = ref_w + ref_e + ref_s + ref_n
2952  offset_band = offset_band / ref_tot
2953 ! offset_bias = offset_bias / ref_tot
2954 
2955  log_info_cont(*) " per dens ", offset_band
2956 
2957  ! density of the reference state is used as weight
2958  if ( .not. prc_has_w ) then
2959  !$omp parallel do
2960  do j = js, je
2961  do k = ks, ke
2962  mflux_offset_x(k,j,1,1) = offset_band * dens_ref(k,is,j)
2963 ! MFLUX_OFFSET_X(k,j,1,2) = offset_bias * DENS_ref(k,IS,j)
2964  end do
2965  end do
2966  end if
2967  if ( .not. prc_has_e ) then
2968  !$omp parallel do
2969  do j = js, je
2970  do k = ks, ke
2971  mflux_offset_x(k,j,2,1) = - offset_band * dens_ref(k,ie,j)
2972 ! MFLUX_OFFSET_X(k,j,2,2) = - offset_bias * DENS_ref(k,IE,j)
2973  end do
2974  end do
2975  end if
2976  if ( .not. prc_has_s ) then
2977  !$omp parallel do
2978  do i = is, ie
2979  do k = ks, ke
2980  mflux_offset_y(k,i,1,1) = offset_band * dens_ref(k,i,js)
2981 ! MFLUX_OFFSET_Y(k,i,1,2) = offset_bias * DENS_ref(k,i,JS)
2982  end do
2983  end do
2984  end if
2985  if ( .not. prc_has_n ) then
2986  !$omp parallel do
2987  do i = is, ie
2988  do k = ks, ke
2989  mflux_offset_y(k,i,2,1) = - offset_band * dens_ref(k,i,je)
2990 ! MFLUX_OFFSET_Y(k,i,2,2) = - offset_bias * DENS_ref(k,i,JE)
2991  end do
2992  end do
2993  end if
2994 
2995  masstot_now = masstot
2996  massflx_now = massflx
2997 
2998  return
2999  end subroutine calc_mass
3000 
3001 
3002  subroutine set_offset
3003  integer :: k, i, j, n
3004 
3005  !$omp parallel do
3006  do j = js, je
3007  do k = ks, ke
3008  do n = 1, 2
3009  atmos_boundary_mflux_offset_x(k,j,n) = mflux_offset_x(k,j,n,1) * offset_time_fact(now_step) !&
3010 ! + MFLUX_OFFSET_X(k,j,n,2)
3011  end do
3012  end do
3013  end do
3014 
3015  !$omp parallel do
3016  do i = is, ie
3017  do k = ks, ke
3018  do n = 1, 2
3019  atmos_boundary_mflux_offset_y(k,i,n) = mflux_offset_y(k,i,n,1) * offset_time_fact(now_step) !&
3020 ! + MFLUX_OFFSET_Y(k,i,n,2)
3021  end do
3022  end do
3023  end do
3024 
3025  return
3026  end subroutine set_offset
3027 
3028 end module mod_atmos_bnd_driver
scale_atmos_grid_cartesc_index::isb
integer, public isb
Definition: scale_atmos_grid_cartesC_index.F90:63
scale_comm_cartesc_nest::comm_cartesc_nest_recv_cancel
subroutine, public comm_cartesc_nest_recv_cancel(HANDLE)
Sub-command for data transfer from parent to daughter: nestdown.
Definition: scale_comm_cartesC_nest.F90:2407
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
mod_atmos_phy_mp_vars
module Atmosphere / Physics Cloud Microphysics
Definition: mod_atmos_phy_mp_vars.F90:12
mod_atmos_vars::qe
real(rp), dimension(:,:,:,:), allocatable, target, public qe
Definition: mod_atmos_vars.F90:104
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
mod_atmos_bnd_driver::atmos_boundary_alpha_pott
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_pott
Definition: mod_atmos_bnd_driver.F90:56
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
mod_atmos_bnd_driver::atmos_boundary_update_flag
logical, public atmos_boundary_update_flag
Definition: mod_atmos_bnd_driver.F90:64
mod_atmos_phy_mp_vars::qa_mp
integer, public qa_mp
Definition: mod_atmos_phy_mp_vars.F90:77
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:34
mod_atmos_bnd_driver::atmos_boundary_set_online
subroutine atmos_boundary_set_online
Set boundary value for real case experiment [online daughter].
Definition: mod_atmos_bnd_driver.F90:1579
scale_tracer::tracer_unit
character(len=h_short), dimension(qa_max), public tracer_unit
Definition: scale_tracer.F90:40
scale_index
module Index
Definition: scale_index.F90:11
scale_prc_cartesc::prc_has_s
logical, public prc_has_s
Definition: scale_prc_cartesC.F90:50
scale_file_cartesc::file_cartesc_enddef
subroutine, public file_cartesc_enddef(fid)
Exit netCDF file define mode.
Definition: scale_file_cartesC.F90:943
scale_file_cartesc::file_cartesc_def_var
subroutine, public file_cartesc_def_var(fid, varname, desc, unit, dim_type, datatype, vid, standard_name, timeintv, nsteps, cell_measures)
Define a variable to file.
Definition: scale_file_cartesC.F90:3307
mod_atmos_bnd_driver::atmos_boundary_driver_finalize
subroutine, public atmos_boundary_driver_finalize
Finalize boundary value.
Definition: mod_atmos_bnd_driver.F90:1702
scale_calendar::calendar_combine_daysec
real(dp) function, public calendar_combine_daysec(absday, abssec)
Combine day and second.
Definition: scale_calendar.F90:405
scale_atmos_refstate::atmos_refstate_dens
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_dens
refernce density [kg/m3]
Definition: scale_atmos_refstate.F90:40
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
mod_atmos_phy_mp_vars::qs_mp
integer, public qs_mp
Definition: mod_atmos_phy_mp_vars.F90:78
scale_comm_cartesc_nest
module Communication CartesianC nesting
Definition: scale_comm_cartesC_nest.F90:12
scale_comm_cartesc_nest::parent_ka
integer, dimension(2), public parent_ka
parent max number in z-direction (with halo)
Definition: scale_comm_cartesC_nest.F90:57
mod_atmos_bnd_driver::atmos_boundary_driver_send
subroutine, public atmos_boundary_driver_send
Send data to child domain.
Definition: mod_atmos_bnd_driver.F90:1806
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_refstate
module atmosphere / reference state
Definition: scale_atmos_refstate.F90:12
mod_atmos_bnd_driver::atmos_boundary_dens
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_dens
Definition: mod_atmos_bnd_driver.F90:45
mod_atmos_bnd_driver::atmos_boundary_alpha_dens
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_dens
Definition: mod_atmos_bnd_driver.F90:52
mod_atmos_bnd_driver::atmos_boundary_qtrc
real(rp), dimension(:,:,:,:), allocatable, public atmos_boundary_qtrc
Definition: mod_atmos_bnd_driver.F90:50
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfx
face buffer factor (0-1): x
Definition: scale_atmos_grid_cartesC.F90:72
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:97
scale_calendar
module CALENDAR
Definition: scale_calendar.F90:13
mod_atmos_bnd_driver::atmos_boundary_velx
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_velx
Definition: mod_atmos_bnd_driver.F90:47
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:79
mod_atmos_bnd_driver::atmos_boundary_driver_setup
subroutine, public atmos_boundary_driver_setup
Setup.
Definition: mod_atmos_bnd_driver.F90:223
mod_atmos_phy_ch_vars::qs_ch
integer, public qs_ch
Definition: mod_atmos_phy_ch_vars.F90:61
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_file
module file
Definition: scale_file.F90:15
mod_atmos_bnd_driver::qtrc
real(rp), dimension(0-1), allocatable, public qtrc
Definition: mod_atmos_bnd_driver.F90:50
mod_atmos_bnd_driver::set_offset
subroutine set_offset
Definition: mod_atmos_bnd_driver.F90:3003
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:48
scale_atmos_grid_cartesc_index::jeb
integer, public jeb
Definition: scale_atmos_grid_cartesC_index.F90:66
mod_atmos_phy_mp_driver
module atmosphere / physics / cloud microphysics
Definition: mod_atmos_phy_mp_driver.F90:12
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:49
scale_comm_cartesc_nest::comm_cartesc_nest_nestdown
subroutine, public comm_cartesc_nest_nestdown(HANDLE, BND_QA, DENS_send, MOMZ_send, MOMX_send, MOMY_send, RHOT_send, QTRC_send, DENS_recv, VELZ_recv, VELX_recv, VELY_recv, POTT_recv, QTRC_recv)
Boundary data transfer from parent to daughter: nestdown.
Definition: scale_comm_cartesC_nest.F90:1868
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfx
center buffer factor (0-1): x
Definition: scale_atmos_grid_cartesC.F90:70
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_atmos_bnd_driver::atmos_boundary_pott
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_pott
Definition: mod_atmos_bnd_driver.F90:49
scale_const::const_pi
real(rp), public const_pi
pi
Definition: scale_const.F90:31
scale_comm_cartesc_nest::online_iam_parent
logical, public online_iam_parent
a flag to say "I am a parent"
Definition: scale_comm_cartesC_nest.F90:100
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
mod_atmos_bnd_driver::calc_mass
subroutine calc_mass(ref)
Definition: mod_atmos_bnd_driver.F90:2768
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_io
module STDIO
Definition: scale_io.F90:10
scale_comm_cartesc_nest::parent_ia
integer, dimension(2), public parent_ia
parent max number in x-direction (with halo)
Definition: scale_comm_cartesC_nest.F90:58
scale_comm_cartesc_nest::online_boundary_use_qhyd
logical, public online_boundary_use_qhyd
Definition: scale_comm_cartesC_nest.F90:105
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfz
center buffer factor (0-1): z
Definition: scale_atmos_grid_cartesC.F90:46
mod_atmos_bnd_driver::atmos_boundary_alpha_vely
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_vely
Definition: mod_atmos_bnd_driver.F90:55
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:44
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:83
mod_atmos_bnd_driver::atmos_boundary_driver_set
subroutine, public atmos_boundary_driver_set
set
Definition: mod_atmos_bnd_driver.F90:600
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfy
face buffer factor (0-1): y
Definition: scale_atmos_grid_cartesC.F90:73
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:76
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazuy_x
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazuy_x
virtical area (zuy, normal x) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:66
mod_atmos_bnd_driver::bnd_iq
integer, dimension(:), allocatable, public bnd_iq
Definition: mod_atmos_bnd_driver.F90:43
mod_atmos_bnd_driver::bnd_qa
integer, public bnd_qa
Definition: mod_atmos_bnd_driver.F90:42
scale_file_cartesc::file_cartesc_close
subroutine, public file_cartesc_close(fid)
Close a netCDF file.
Definition: scale_file_cartesC.F90:1023
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfz
face buffer factor (0-1): z
Definition: scale_atmos_grid_cartesC.F90:47
scale_calendar::calendar_date2daysec
subroutine, public calendar_date2daysec(absday, abssec, ymdhms, subsec, offset_year)
Convert from gregorian date to absolute day/second.
Definition: scale_calendar.F90:159
mod_atmos_bnd_driver
module ATMOSPHERE / Boundary treatment
Definition: mod_atmos_bnd_driver.F90:13
mod_atmos_bnd_driver::atmos_boundary_smoother_fact
real(rp), public atmos_boundary_smoother_fact
Definition: mod_atmos_bnd_driver.F90:62
scale_tracer::tracer_name
character(len=h_short), dimension(qa_max), public tracer_name
Definition: scale_tracer.F90:38
scale_prof
module profiler
Definition: scale_prof.F90:11
mod_atmos_phy_ch_vars
module Atmosphere / Physics Chemistry
Definition: mod_atmos_phy_ch_vars.F90:12
mod_atmos_vars::momx
real(rp), dimension(:,:,:), allocatable, target, public momx
Definition: mod_atmos_vars.F90:77
mod_atmos_bnd_driver::atmos_boundary_alpha_velz
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_velz
Definition: mod_atmos_bnd_driver.F90:53
scale_comm_cartesc_nest::online_iam_daughter
logical, public online_iam_daughter
a flag to say "I am a daughter"
Definition: scale_comm_cartesC_nest.F90:101
mod_atmos_bnd_driver::set_boundary
subroutine set_boundary(use_velz)
Definition: mod_atmos_bnd_driver.F90:2077
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_atmos_vars::momy
real(rp), dimension(:,:,:), allocatable, target, public momy
Definition: mod_atmos_vars.F90:78
scale_time
module TIME
Definition: scale_time.F90:11
mod_atmos_vars::qv
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
Definition: mod_atmos_vars.F90:96
mod_atmos_bnd_driver::atmos_boundary_set_file
subroutine atmos_boundary_set_file
Set boundary value for real case experiment.
Definition: mod_atmos_bnd_driver.F90:1412
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
mod_atmos_bnd_driver::atmos_boundary_alpha_qtrc
real(rp), dimension(:,:,:,:), allocatable, public atmos_boundary_alpha_qtrc
Definition: mod_atmos_bnd_driver.F90:57
mod_atmos_bnd_driver::atmos_boundary_vely
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_vely
Definition: mod_atmos_bnd_driver.F90:48
scale_comm_cartesc_nest::parent_nstep
integer, dimension(2), public parent_nstep
parent step [number]
Definition: scale_comm_cartesC_nest.F90:63
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:77
mod_atmos_bnd_driver::atmos_boundary_velz
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_velz
Definition: mod_atmos_bnd_driver.F90:46
scale_comm_cartesc_nest::daughter_ia
integer, dimension(2), public daughter_ia
daughter max number in x-direction (with halo)
Definition: scale_comm_cartesC_nest.F90:69
scale_comm_cartesc_nest::parent_dtsec
real(dp), dimension(2), public parent_dtsec
parent DT [sec]
Definition: scale_comm_cartesC_nest.F90:62
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
mod_atmos_bnd_driver::dens
real(rp), dimension(0-1), allocatable, public dens
Definition: mod_atmos_bnd_driver.F90:45
scale_file::file_aggregate
logical, public file_aggregate
Definition: scale_file.F90:182
scale_file_cartesc::file_cartesc_create
subroutine, public file_cartesc_create(basename, title, datatype, fid, date, subsec, haszcoord, append, aggregate, single)
Create/open a netCDF file.
Definition: scale_file_cartesC.F90:780
scale_comm_cartesc_nest::online_use_velz
logical, public online_use_velz
Definition: scale_comm_cartesC_nest.F90:103
scale_file_cartesc::file_cartesc_open
subroutine, public file_cartesc_open(basename, fid, aggregate)
open a netCDF file for read
Definition: scale_file_cartesC.F90:746
scale_file_cartesc::file_cartesc_flush
subroutine, public file_cartesc_flush(fid)
Flush all pending requests to a netCDF file (PnetCDF only)
Definition: scale_file_cartesC.F90:997
scale_time::time_dtsec
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:33
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_comm_cartesc_nest::use_nesting
logical, public use_nesting
Definition: scale_comm_cartesC_nest.F90:98
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.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_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totareazuy_x
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_real_totareazuy_x
total area (zuy, normal x) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:80
mod_atmos_bnd_driver::real
logical, public real
Definition: mod_atmos_bnd_driver.F90:64
scale_comm_cartesc_nest::online_boundary_diagqhyd
logical, public online_boundary_diagqhyd
Definition: scale_comm_cartesC_nest.F90:106
scale_atmos_grid_cartesc_index::ieb
integer, public ieb
Definition: scale_atmos_grid_cartesC_index.F90:64
mod_atmos_phy_mp_vars::qe_mp
integer, public qe_mp
Definition: mod_atmos_phy_mp_vars.F90:79
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfy
center buffer factor (0-1): y
Definition: scale_atmos_grid_cartesC.F90:71
scale_comm_cartesc_nest::daughter_ja
integer, dimension(2), public daughter_ja
daughter max number in y-direction (with halo)
Definition: scale_comm_cartesC_nest.F90:70
mod_atmos_bnd_driver::atmos_boundary_mflux_offset_y
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_mflux_offset_y
Definition: mod_atmos_bnd_driver.F90:60
mod_atmos_phy_ch_vars::qe_ch
integer, public qe_ch
Definition: mod_atmos_phy_ch_vars.F90:62
mod_atmos_bnd_driver::atmos_boundary_mflux_offset_x
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_mflux_offset_x
Definition: mod_atmos_bnd_driver.F90:59
scale_time::time_nstep
integer, public time_nstep
total steps [number]
Definition: scale_time.F90:72
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol
real(rp), public atmos_grid_cartesc_real_totvol
total volume (zxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:87
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_comm_cartesc_nest::parent_ja
integer, dimension(2), public parent_ja
parent max number in y-direction (with halo)
Definition: scale_comm_cartesC_nest.F90:59
mod_atmos_bnd_driver::atmos_boundary_alpha_velx
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_velx
Definition: mod_atmos_bnd_driver.F90:54
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_comm_cartesc_nest::comm_cartesc_nest_test
subroutine, public comm_cartesc_nest_test(HANDLE)
[check communication status] Inter-communication
Definition: scale_comm_cartesC_nest.F90:2813
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazxv_y
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazxv_y
virtical area (zxv, normal y) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:67
scale_comm_cartesc_nest::daughter_ka
integer, dimension(2), public daughter_ka
daughter max number in z-direction (with halo)
Definition: scale_comm_cartesC_nest.F90:68
scale_atmos_grid_cartesc_index::jsb
integer, public jsb
Definition: scale_atmos_grid_cartesC_index.F90:65
mod_atmos_bnd_driver::atmos_boundary_driver_update
subroutine, public atmos_boundary_driver_update(last_step)
Update boundary value with a constant time boundary.
Definition: mod_atmos_bnd_driver.F90:1736
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:47
mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc
subroutine, public atmos_phy_mp_driver_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
Definition: mod_atmos_phy_mp_driver.F90:1361
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:56
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_file_cartesc
module file / cartesianC
Definition: scale_file_cartesC.F90:11
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:79
scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue
subroutine, public comm_cartesc_nest_recvwait_issue(HANDLE, BND_QA)
Sub-command for data transfer from parent to daughter: nestdown.
Definition: scale_comm_cartesC_nest.F90:2303
scale_comm_cartesc_nest::comm_cartesc_nest_bnd_qa
integer, public comm_cartesc_nest_bnd_qa
number of tracer treated in nesting system
Definition: scale_comm_cartesC_nest.F90:94
scale_calendar::calendar_date2char
subroutine, public calendar_date2char(chardate, ymdhms, subsec)
Convert from gregorian date to absolute day/second.
Definition: scale_calendar.F90:603
scale_comm_cartesc_nest::offline
logical, public offline
Definition: scale_comm_cartesC_nest.F90:99
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totareazxv_y
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_real_totareazxv_y
total area (zxv, normal y) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:81
mod_atmos_bnd_driver::update_ref_index
subroutine update_ref_index
Update indices of array of boundary references.
Definition: mod_atmos_bnd_driver.F90:2716