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_setalpha
73  private :: atmos_boundary_setinitval
74  private :: atmos_boundary_read
75  private :: atmos_boundary_write
76  private :: atmos_boundary_generate
77  private :: atmos_boundary_initialize_file
78  private :: atmos_boundary_initialize_online
79  private :: atmos_boundary_update
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  !
88  !-----------------------------------------------------------------------------
89  !
90  !++ Private parameters & variables
91  !
92  character(len=H_SHORT), private :: atmos_boundary_type = 'NONE'
93  character(len=H_LONG), private :: atmos_boundary_in_basename = ''
94  logical, private :: atmos_boundary_in_basename_add_num = .false.
95  integer, private :: atmos_boundary_in_number_of_files = 1
96  logical, private :: atmos_boundary_in_check_coordinates = .true.
97  integer, private :: atmos_boundary_in_step_limit = 1000
98  logical, private :: atmos_boundary_in_aggregate
99  character(len=H_LONG), private :: atmos_boundary_out_basename = ''
100  character(len=H_MID), private :: atmos_boundary_out_title = 'SCALE-RM BOUNDARY CONDITION'
101  character(len=H_SHORT), private :: atmos_boundary_out_dtype = 'DEFAULT'
102  logical, private :: atmos_boundary_out_aggregate
103 
104  logical, private :: atmos_boundary_use_dens = .false. ! read from file?
105  logical, private :: atmos_boundary_use_velz = .false. ! read from file?
106  logical, private :: atmos_boundary_use_velx = .false. ! read from file?
107  logical, private :: atmos_boundary_use_vely = .false. ! read from file?
108  logical, private :: atmos_boundary_use_pt = .false. ! read from file?
109  logical, private :: atmos_boundary_use_qv = .false. ! read from file?
110  logical, private :: atmos_boundary_use_qhyd = .false. ! read from file?
111  logical, private :: atmos_boundary_use_chem = .false. ! read from file?
112 
113  real(rp), private :: atmos_boundary_value_velz = 0.0_rp ! velocity w at boundary, 0 [m/s]
114  real(rp), private :: atmos_boundary_value_velx = 0.0_rp ! velocity u at boundary, 0 [m/s]
115  real(rp), private :: atmos_boundary_value_vely = 0.0_rp ! velocity v at boundary, 0 [m/s]
116  real(rp), private :: atmos_boundary_value_pt = 300.0_rp ! potential temp. at boundary, 300 [K]
117  real(rp), private :: atmos_boundary_value_qtrc = 0.0_rp ! tracer at boundary, 0 [kg/kg]
118 
119  real(rp), private :: atmos_boundary_alphafact_dens = 1.0_rp ! alpha factor again default
120  real(rp), private :: atmos_boundary_alphafact_velz = 1.0_rp ! alpha factor again default
121  real(rp), private :: atmos_boundary_alphafact_velx = 1.0_rp ! alpha factor again default
122  real(rp), private :: atmos_boundary_alphafact_vely = 1.0_rp ! alpha factor again default
123  real(rp), private :: atmos_boundary_alphafact_pt = 1.0_rp ! alpha factor again default
124  real(rp), private :: atmos_boundary_alphafact_qtrc = 1.0_rp ! alpha factor again default
125 
126  real(rp), private :: atmos_boundary_fracz = 1.0_rp ! fraction of boundary region for dumping (z) (0-1)
127  real(rp), private :: atmos_boundary_fracx = 1.0_rp ! fraction of boundary region for dumping (x) (0-1)
128  real(rp), private :: atmos_boundary_fracy = 1.0_rp ! fraction of boundary region for dumping (y) (0-1)
129  real(rp), private :: atmos_boundary_tauz ! maximum value for damping tau (z) [s]
130  real(rp), private :: atmos_boundary_taux ! maximum value for damping tau (x) [s]
131  real(rp), private :: atmos_boundary_tauy ! maximum value for damping tau (y) [s]
132 
133  real(dp), private :: update_dt = 0.0_dp ! inteval time of boudary data update [s]
134  integer, private :: update_nstep
135 
136  logical, private :: atmos_grid_nudging_uv = .false. ! grid nudging
137  logical, private :: atmos_grid_nudging_pt = .false. ! grid nudging
138  logical, private :: atmos_grid_nudging_qv = .false. ! grid nudging
139  real(rp), private :: atmos_grid_nudging_tau ! Damping tau for grid nudging [s]
140 
141  ! work
142  real(rp), private, allocatable :: dens_ref(:,:,:)
143  real(rp), private, allocatable :: velz_ref(:,:,:)
144  real(rp), private, allocatable :: velx_ref(:,:,:)
145  real(rp), private, allocatable :: vely_ref(:,:,:)
146  real(rp), private, allocatable :: pott_ref(:,:,:)
147  real(rp), private, allocatable, target :: qtrc_ref(:,:,:,:)
148  real(rp), private, allocatable, target :: q_send_work(:,:,:,:) ! QV + Qe
149  real(rp), private, allocatable, target :: q_recv_work(:,:,:,:) ! QV + Qe
150 
151 
152  integer, private :: now_step
153  logical, private :: atmos_boundary_linear_v = .false. ! linear or non-linear profile of relax region
154  logical, private :: atmos_boundary_linear_h = .true. ! linear or non-linear profile of relax region
155  real(rp), private :: atmos_boundary_exp_h = 2.0_rp ! factor of non-linear profile of relax region
156  logical, private :: atmos_boundary_online = .false. ! boundary online update by communicate inter-domain
157  logical, private :: atmos_boundary_dens_adjust = .false.
158  real(rp), private :: atmos_boundary_dens_adjust_tau
159 
160  logical, private :: do_parent_process = .false.
161  logical, private :: do_daughter_process = .false.
162  logical, private :: l_bnd = .false.
163 
164  ! for mass flux offset
165  real(dp), private :: masstot_now = 0.0_dp
166  real(dp), private :: massflx_now = 0.0_dp
167  real(rp), allocatable, private :: areazuy_w(:,:), areazuy_e(:,:)
168  real(rp), allocatable, private :: offset_time_fact(:)
169  real(rp), allocatable, private :: mflux_offset_x(:,:,:,:)
170  real(rp), allocatable, private :: mflux_offset_y(:,:,:,:)
171  real(rp), allocatable, private, target :: zero_x(:,:), zero_y(:,:)
172 
173 
174  !-----------------------------------------------------------------------------
175 contains
176  !-----------------------------------------------------------------------------
178  subroutine atmos_boundary_driver_setup
179  use scale_prc, only: &
180  prc_abort
181  use scale_const, only: &
182  undef => const_undef
183  use scale_time, only: &
184  dt => time_dtsec
185  use scale_file, only: &
187  use scale_comm_cartesc_nest, only: &
188  online_use_velz, &
190  use_nesting, &
195  online_send_qa, &
197  use scale_atmos_hydrometeor, only: &
199  i_qv
200  use mod_atmos_phy_mp_vars, only: &
201  qs_mp, &
202  qe_mp
203  use mod_atmos_phy_ch_vars, only: &
204  qs_ch, &
205  qe_ch
206  use scale_atmos_grid_cartesc_real, only: &
208  implicit none
209 
210  namelist / param_atmos_boundary / &
211  atmos_boundary_type, &
212  atmos_boundary_in_basename, &
213  atmos_boundary_in_basename_add_num, &
214  atmos_boundary_in_number_of_files, &
215  atmos_boundary_in_check_coordinates, &
216  atmos_boundary_in_step_limit, &
217  atmos_boundary_in_aggregate, &
218  atmos_boundary_out_basename, &
219  atmos_boundary_out_title, &
220  atmos_boundary_out_dtype, &
221  atmos_boundary_out_aggregate, &
222  atmos_boundary_use_velz, &
223  atmos_boundary_use_velx, &
224  atmos_boundary_use_vely, &
225  atmos_boundary_use_pt, &
226  atmos_boundary_use_dens, &
227  atmos_boundary_use_qv, &
228  atmos_boundary_use_qhyd, &
229  atmos_boundary_use_chem, &
230  atmos_boundary_dens_adjust, &
231  atmos_boundary_dens_adjust_tau, &
232  atmos_boundary_value_velz, &
233  atmos_boundary_value_velx, &
234  atmos_boundary_value_vely, &
235  atmos_boundary_value_pt, &
236  atmos_boundary_value_qtrc, &
237  atmos_boundary_alphafact_dens, &
238  atmos_boundary_alphafact_velz, &
239  atmos_boundary_alphafact_velx, &
240  atmos_boundary_alphafact_vely, &
241  atmos_boundary_alphafact_pt, &
242  atmos_boundary_alphafact_qtrc, &
244  atmos_boundary_fracz, &
245  atmos_boundary_fracx, &
246  atmos_boundary_fracy, &
247  atmos_boundary_tauz, &
248  atmos_boundary_taux, &
249  atmos_boundary_tauy, &
250  atmos_boundary_linear_v, &
251  atmos_boundary_linear_h, &
252  atmos_boundary_exp_h, &
253  atmos_grid_nudging_uv, &
254  atmos_grid_nudging_pt, &
255  atmos_grid_nudging_qv, &
256  atmos_grid_nudging_tau
257 
258  integer :: k, i, j, iq
259  integer :: ierr
260  !---------------------------------------------------------------------------
261 
262  log_newline
263  log_info("ATMOS_BOUNDARY_setup",*) 'Setup'
264 
265 
266  atmos_boundary_tauz = dt * 10.0_rp
267  atmos_boundary_taux = dt * 10.0_rp
268  atmos_boundary_tauy = dt * 10.0_rp
269 
270  atmos_boundary_in_aggregate = file_aggregate
271  atmos_boundary_out_aggregate = file_aggregate
272 
273  atmos_boundary_dens_adjust = .false.
274  atmos_boundary_dens_adjust_tau = -1.0_rp
275 
276  atmos_grid_nudging_tau = 10.0_rp * 24.0_rp * 3600.0_rp ! 10days [s]
277 
278 
279  !--- read namelist
280  rewind(io_fid_conf)
281  read(io_fid_conf,nml=param_atmos_boundary,iostat=ierr)
282  if( ierr < 0 ) then !--- missing
283  log_info("ATMOS_BOUNDARY_setup",*) 'Not found namelist. Default used.'
284  elseif( ierr > 0 ) then !--- fatal error
285  log_error("ATMOS_BOUNDARY_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_BOUNDARY. Check!'
286  call prc_abort
287  endif
288  log_nml(param_atmos_boundary)
289 
290  ! setting switches
291  if( .NOT. use_nesting ) then
292  atmos_boundary_online = .false.
293  else
294  atmos_boundary_online = .true.
295  endif
296  do_parent_process = .false.
297  do_daughter_process = .false.
298  if ( atmos_boundary_online ) then
299  if ( online_iam_parent ) then
300  do_parent_process = .true.
301  endif
302  if ( online_iam_daughter ) then
303  do_daughter_process = .true.
304  atmos_boundary_use_velz = online_use_velz
305  atmos_boundary_use_qhyd = online_boundary_use_qhyd
306  endif
307  endif
308 
309  allocate( bnd_iq(qa) )
310  bnd_iq(:) = -1
311  bnd_qa = 0
312  if ( .not. atmos_hydrometeor_dry ) then
313  bnd_qa = bnd_qa + 1
314  bnd_iq(i_qv) = bnd_qa
315  if( atmos_boundary_use_qhyd ) then
316  do iq = qs_mp+1, qe_mp
317  bnd_qa = bnd_qa + 1
318  bnd_iq(iq) = bnd_qa
319  end do
320  end if
321  end if
322  if( atmos_boundary_use_chem ) then
323  do iq = qs_ch, qe_ch
324  bnd_qa = bnd_qa + 1
325  bnd_iq(iq) = bnd_qa
326  end do
327  endif
328  !$acc enter data copyin(BND_IQ)
329 
330  allocate( atmos_boundary_dens(ka,ia,ja) )
331  allocate( atmos_boundary_velz(ka,ia,ja) )
332  allocate( atmos_boundary_velx(ka,ia,ja) )
333  allocate( atmos_boundary_vely(ka,ia,ja) )
334  allocate( atmos_boundary_pott(ka,ia,ja) )
335  allocate( atmos_boundary_qtrc(ka,ia,ja,max(bnd_qa,1)) )
336  atmos_boundary_dens(:,:,:) = undef
337  atmos_boundary_velz(:,:,:) = undef
338  atmos_boundary_velx(:,:,:) = undef
339  atmos_boundary_vely(:,:,:) = undef
340  atmos_boundary_pott(:,:,:) = undef
341  atmos_boundary_qtrc(:,:,:,:) = undef
342  !$acc enter data create(ATMOS_BOUNDARY_DENS, ATMOS_BOUNDARY_VELZ, ATMOS_BOUNDARY_VELX, ATMOS_BOUNDARY_VELY, ATMOS_BOUNDARY_POTT, ATMOS_BOUNDARY_QTRC)
343 
344  allocate( atmos_boundary_alpha_dens(ka,ia,ja) )
345  allocate( atmos_boundary_alpha_velz(ka,ia,ja) )
346  allocate( atmos_boundary_alpha_velx(ka,ia,ja) )
347  allocate( atmos_boundary_alpha_vely(ka,ia,ja) )
348  allocate( atmos_boundary_alpha_pott(ka,ia,ja) )
349  allocate( atmos_boundary_alpha_qtrc(ka,ia,ja,max(bnd_qa,1)) )
350  atmos_boundary_alpha_dens(:,:,:) = 0.0_rp
351  atmos_boundary_alpha_velz(:,:,:) = 0.0_rp
352  atmos_boundary_alpha_velx(:,:,:) = 0.0_rp
353  atmos_boundary_alpha_vely(:,:,:) = 0.0_rp
354  atmos_boundary_alpha_pott(:,:,:) = 0.0_rp
355  atmos_boundary_alpha_qtrc(:,:,:,:) = 0.0_rp
356  !$acc enter data copyin(ATMOS_BOUNDARY_alpha_DENS, ATMOS_BOUNDARY_alpha_VELZ, ATMOS_BOUNDARY_alpha_VELX, ATMOS_BOUNDARY_alpha_VELY, ATMOS_BOUNDARY_alpha_POTT, ATMOS_BOUNDARY_alpha_QTRC)
357 
358  allocate( atmos_boundary_mflux_offset_x(ka,ja,2) )
359  allocate( atmos_boundary_mflux_offset_y(ka,ia,2) )
360  atmos_boundary_mflux_offset_x(:,:,:) = 0.0_rp
361  atmos_boundary_mflux_offset_y(:,:,:) = 0.0_rp
362  !$acc enter data copyin(ATMOS_BOUNDARY_MFLUX_OFFSET_X, ATMOS_BOUNDARY_MFLUX_OFFSET_Y)
363 
364  if ( atmos_boundary_type == 'REAL' .OR. do_daughter_process ) then
365  l_bnd = .true.
366  else
367  l_bnd = .false.
368  end if
369 
370  if ( l_bnd ) then
371 
372  allocate( dens_ref(ka,ia,ja) )
373  allocate( velx_ref(ka,ia,ja) )
374  allocate( vely_ref(ka,ia,ja) )
375  allocate( pott_ref(ka,ia,ja) )
376  allocate( qtrc_ref(ka,ia,ja,max(bnd_qa,1)) )
377  dens_ref(:,:,:) = 0.0_rp
378  velx_ref(:,:,:) = 0.0_rp
379  vely_ref(:,:,:) = 0.0_rp
380  pott_ref(:,:,:) = 0.0_rp
381  qtrc_ref(:,:,:,:) = 0.0_rp
382  !$acc enter data copyin(DENS_ref, VELX_ref, VELY_ref, POTT_ref, QTRC_ref)
383  if ( atmos_boundary_use_velz ) then
384  allocate( velz_ref(ka,ia,ja) )
385  velz_ref(:,:,:) = 0.0_rp
386  !$acc enter data copyin(VELZ_ref)
387  end if
388 
389  ! initialize boundary value (reading file or waiting parent domain)
390  if ( do_daughter_process ) then
391  call atmos_boundary_initialize_online
392  else
393  if ( atmos_boundary_in_basename /= '' ) then
394  call atmos_boundary_initialize_file
395  else
396  log_error("ATMOS_BOUNDARY_setup",*) 'You need specify ATMOS_BOUNDARY_IN_BASENAME'
397  call prc_abort
398  endif
399  endif
400 
401  if ( atmos_boundary_dens_adjust_tau <= 0.0_rp ) then
402  atmos_boundary_dens_adjust_tau = max( real(update_dt,kind=rp) / 6.0_rp, &
403  atmos_boundary_taux, atmos_boundary_tauy )
404  end if
405 
406  call atmos_boundary_setalpha
407 
409 
410  ! for mass flux offset
411  if ( atmos_boundary_dens_adjust ) then
412  allocate( areazuy_w(ka,ja), areazuy_e(ka,ja) )
413  allocate( mflux_offset_x(ka,ja,2,2), mflux_offset_y(ka,ia,2,2) )
414  allocate( zero_x(ka,ja), zero_y(ka,ia) )
415 
416  !$omp parallel do
417  do j = js, je
418  do k = ks, ke
419  areazuy_w(k,j) = atmos_grid_cartesc_real_areazuy_x(k,is-1,j)
420  areazuy_e(k,j) = atmos_grid_cartesc_real_areazuy_x(k,ie ,j)
421  end do
422  end do
423  !$omp parallel do
424  do j = js, je
425  do k = ks, ke
426  mflux_offset_x(k,j,:,:) = 0.0_rp
427  end do
428  end do
429  !$omp parallel do
430  do i = is, ie
431  do k = ks, ke
432  mflux_offset_y(k,i,:,:) = 0.0_rp
433  end do
434  end do
435 
436  !$omp parallel do
437  do j = js, je
438  do k = ks, ke
439  zero_x(k,j) = 0.0_rp
440  end do
441  end do
442  !$omp parallel do
443  do i = is, ie
444  do k = ks, ke
445  zero_y(k,i) = 0.0_rp
446  end do
447  end do
448 
449  !$acc enter data copyin(AREAZUY_W, AREAZUY_E, MFLUX_OFFSET_X, MFLUX_OFFSET_Y, zero_x, zero_y)
450  end if
451 
452  elseif ( atmos_boundary_type == 'NONE' ) then
453 
455 
456  elseif ( atmos_boundary_type == 'CONST' ) then
457 
458  call atmos_boundary_setalpha
459 
461 
462  elseif ( atmos_boundary_type == 'INIT' ) then
463 
464  call atmos_boundary_setalpha
465 
467 
468  elseif ( atmos_boundary_type == 'OFFLINE' ) then
469 
470  if ( atmos_boundary_in_basename /= '' ) then
471  call atmos_boundary_read
472  else
473  log_error("ATMOS_BOUNDARY_setup",*) 'You need specify ATMOS_BOUNDARY_IN_BASENAME'
474  call prc_abort
475  endif
476 
478 
479  else
480  log_error("ATMOS_BOUNDARY_setup",*) 'unsupported ATMOS_BOUNDARY_TYPE. Check!', trim(atmos_boundary_type)
481  call prc_abort
482  endif
483 
485 
486 
487  !----- report data -----
488  log_newline
489  log_info("ATMOS_BOUNDARY_setup",*) 'Atmospheric boundary parameters '
490  log_info_cont(*) 'Atmospheric boundary type : ', atmos_boundary_type
491  log_newline
492  log_info_cont(*) 'Is VELZ used in atmospheric boundary? : ', atmos_boundary_use_velz
493  log_info_cont(*) 'Is VELX used in atmospheric boundary? : ', atmos_boundary_use_velx
494  log_info_cont(*) 'Is VELY used in atmospheric boundary? : ', atmos_boundary_use_vely
495  log_info_cont(*) 'Is PT used in atmospheric boundary? : ', atmos_boundary_use_pt
496  log_info_cont(*) 'Is DENS used in atmospheric boundary? : ', atmos_boundary_use_dens
497  log_info_cont(*) 'Is QV used in atmospheric boundary? : ', atmos_boundary_use_qv
498  log_info_cont(*) 'Is QHYD used in atmospheric boundary? : ', atmos_boundary_use_qhyd
499  log_info_cont(*) 'Is CHEM used in atmospheric boundary? : ', atmos_boundary_use_chem
500  log_newline
501  log_info_cont(*) 'Atmospheric boundary VELZ values : ', atmos_boundary_value_velz
502  log_info_cont(*) 'Atmospheric boundary VELX values : ', atmos_boundary_value_velx
503  log_info_cont(*) 'Atmospheric boundary VELY values : ', atmos_boundary_value_vely
504  log_info_cont(*) 'Atmospheric boundary PT values : ', atmos_boundary_value_pt
505  log_info_cont(*) 'Atmospheric boundary QTRC values : ', atmos_boundary_value_qtrc
506  log_newline
507  log_info_cont(*) 'Atmospheric boundary smoother factor : ', atmos_boundary_smoother_fact
508  log_info_cont(*) 'Atmospheric boundary z-fraction : ', atmos_boundary_fracz
509  log_info_cont(*) 'Atmospheric boundary x-fraction : ', atmos_boundary_fracx
510  log_info_cont(*) 'Atmospheric boundary y-fraction : ', atmos_boundary_fracy
511  log_info_cont(*) 'Atmospheric boundary z-relaxation time : ', atmos_boundary_tauz
512  log_info_cont(*) 'Atmospheric boundary x-relaxation time : ', atmos_boundary_taux
513  log_info_cont(*) 'Atmospheric boundary y-relaxation time : ', atmos_boundary_tauy
514  log_newline
515  log_newline
516  log_info_cont(*) 'Linear profile in vertically relax region : ', atmos_boundary_linear_v
517  log_info_cont(*) 'Linear profile in horizontally relax region : ', atmos_boundary_linear_h
518  log_info_cont(*) 'Non-linear factor in horizontally relax region : ', atmos_boundary_exp_h
519  log_newline
520  log_info_cont(*) 'Online nesting for lateral boundary : ', atmos_boundary_online
521  log_info_cont(*) 'Does lateral boundary exist in this domain? : ', l_bnd
522  log_newline
523  log_info_cont(*) 'Is grid nudging used for VELX & VELY? : ', atmos_grid_nudging_uv
524  log_info_cont(*) 'Is grid nudging used for POTT? : ', atmos_grid_nudging_pt
525  log_info_cont(*) 'Is grid nudging used for QV? : ', atmos_grid_nudging_qv
526  log_info_cont(*) 'Relaxation time for grid nudging : ', atmos_grid_nudging_tau
527  log_info_cont(*) 'Density adjustment : ', atmos_boundary_dens_adjust
528  if ( atmos_boundary_dens_adjust ) then
529  log_info_cont(*) 'Density relaxation time : ', atmos_boundary_dens_adjust_tau
530  end if
531 
532  if ( online_send_diagqhyd ) then
533  allocate( q_send_work(ka,ia,ja,online_send_qa) )
534  !$acc enter data create( Q_SEND_WORK )
535  end if
536  if ( online_recv_diagqhyd ) then
537  allocate( q_recv_work(ka,ia,ja,online_recv_qa) )
538  !$acc enter data create( Q_RECV_WORK )
539  end if
540 
541  return
542  end subroutine atmos_boundary_driver_setup
543 
544  !-----------------------------------------------------------------------------
546  subroutine atmos_boundary_driver_set( &
547  time )
548  use scale_const, only: &
549  pi => const_pi
550  use mod_atmos_vars, only: &
551  dens, &
552  momz, &
553  momx, &
554  momy, &
555  rhot, &
556  qtrc, &
557  qv, &
558  qe
559  use mod_atmos_phy_mp_vars, only: &
560  qs_mp, &
561  qe_mp
562  implicit none
563  real(dp), intent(in) :: time
564 
565  real(rp) :: total
566  integer :: n
567 
568  if ( do_parent_process ) then !online [parent]
569  call atmos_boundary_firstsend( &
570  dens, momz, momx, momy, rhot, qtrc(:,:,:,qs_mp:qe_mp), qv, qe )
571  end if
572 
573  if ( l_bnd ) then
574 
575  if ( atmos_boundary_dens_adjust ) then
576  now_step = 0
577  allocate( offset_time_fact(0:update_nstep) )
578  total = 0.0_rp
579  !$omp parallel do reduction(+:total)
580  do n = 0, update_nstep
581  offset_time_fact(n) = 1.0_rp - cos( 2.0_rp * pi * ( n - 1 ) / update_nstep )
582  total = total + offset_time_fact(n)
583  end do
584  total = total / update_nstep
585  !$omp parallel do
586  do n = 0, update_nstep
587  offset_time_fact(n) = offset_time_fact(n) / total
588  end do
589  end if
590 
591  ! initialize boundary value (reading file or waiting parent domain)
592  if ( do_daughter_process ) then
593  call atmos_boundary_set_online( time )
594  else
595  if ( atmos_boundary_in_basename /= '' ) then
596  call atmos_boundary_set_file( time )
597  endif
598  endif
599 
600  elseif ( atmos_boundary_type == 'CONST' ) then
601 
602  call atmos_boundary_generate
603 
604  elseif ( atmos_boundary_type == 'INIT' ) then
605 
606  call atmos_boundary_setinitval( dens, & ! [IN]
607  momz, & ! [IN]
608  momx, & ! [IN]
609  momy, & ! [IN]
610  rhot, & ! [IN]
611  qtrc ) ! [IN]
612  endif
613 
614  if( atmos_boundary_out_basename /= '' ) then
615  call atmos_boundary_write
616  endif
617 
618  if ( atmos_boundary_update_flag ) then
619 
620  call history_bnd( &
627  end if
628 
629  return
630  end subroutine atmos_boundary_driver_set
631 
632  !-----------------------------------------------------------------------------
634  subroutine atmos_boundary_var_fillhalo
635  use scale_comm_cartesc, only: &
636  comm_vars8, &
637  comm_wait
638  implicit none
639 
640  integer :: i, j, iq
641  !---------------------------------------------------------------------------
642 
643 !!$ do j = 1, JA
644 !!$ do i = 1, IA
645 !!$ ATMOS_BOUNDARY_DENS( 1:KS-1,i,j) = ATMOS_BOUNDARY_DENS(KS,i,j)
646 !!$ ATMOS_BOUNDARY_VELZ( 1:KS-1,i,j) = ATMOS_BOUNDARY_VELZ(KS,i,j)
647 !!$ ATMOS_BOUNDARY_VELX( 1:KS-1,i,j) = ATMOS_BOUNDARY_VELX(KS,i,j)
648 !!$ ATMOS_BOUNDARY_VELY( 1:KS-1,i,j) = ATMOS_BOUNDARY_VELY(KS,i,j)
649 !!$ ATMOS_BOUNDARY_POTT( 1:KS-1,i,j) = ATMOS_BOUNDARY_POTT(KS,i,j)
650 !!$
651 !!$ ATMOS_BOUNDARY_DENS(KE+1:KA, i,j) = ATMOS_BOUNDARY_DENS(KE,i,j)
652 !!$ ATMOS_BOUNDARY_VELZ(KE+1:KA, i,j) = ATMOS_BOUNDARY_VELZ(KE,i,j)
653 !!$ ATMOS_BOUNDARY_VELX(KE+1:KA, i,j) = ATMOS_BOUNDARY_VELX(KE,i,j)
654 !!$ ATMOS_BOUNDARY_VELY(KE+1:KA, i,j) = ATMOS_BOUNDARY_VELY(KE,i,j)
655 !!$ ATMOS_BOUNDARY_POTT(KE+1:KA, i,j) = ATMOS_BOUNDARY_POTT(KE,i,j)
656 !!$
657 !!$ do iq = 1, BND_QA
658 !!$ ATMOS_BOUNDARY_QTRC( 1:KS-1,i,j,iq) = ATMOS_BOUNDARY_QTRC(KS,i,j,iq)
659 !!$ ATMOS_BOUNDARY_QTRC(KE+1:KA, i,j,iq) = ATMOS_BOUNDARY_QTRC(KE,i,j,iq)
660 !!$ end do
661 !!$ end do
662 !!$ end do
663 
664  call comm_vars8( atmos_boundary_dens(:,:,:), 1 )
665  call comm_vars8( atmos_boundary_velz(:,:,:), 2 )
666  call comm_vars8( atmos_boundary_velx(:,:,:), 3 )
667  call comm_vars8( atmos_boundary_vely(:,:,:), 4 )
668  call comm_vars8( atmos_boundary_pott(:,:,:), 5 )
669  do iq = 1, bnd_qa
670  call comm_vars8( atmos_boundary_qtrc(:,:,:,iq), 5+iq )
671  end do
672 
673  call comm_wait ( atmos_boundary_dens(:,:,:), 1, .false. )
674  call comm_wait ( atmos_boundary_velz(:,:,:), 2, .false. )
675  call comm_wait ( atmos_boundary_velx(:,:,:), 3, .false. )
676  call comm_wait ( atmos_boundary_vely(:,:,:), 4, .false. )
677  call comm_wait ( atmos_boundary_pott(:,:,:), 5, .false. )
678  do iq = 1, bnd_qa
679  call comm_wait ( atmos_boundary_qtrc(:,:,:,iq), 5+iq, .false. )
680  end do
681 
682  return
683  end subroutine atmos_boundary_var_fillhalo
684 
685  !-----------------------------------------------------------------------------
687  subroutine atmos_boundary_alpha_fillhalo
688  use scale_comm_cartesc, only: &
689  comm_vars8, &
690  comm_wait
691  implicit none
692 
693  integer :: i, j, iq
694  !---------------------------------------------------------------------------
695 
696 !!$ do j = 1, JA
697 !!$ do i = 1, IA
698 !!$ ATMOS_BOUNDARY_alpha_DENS( 1:KS-1,i,j) = ATMOS_BOUNDARY_alpha_DENS(KS,i,j)
699 !!$ ATMOS_BOUNDARY_alpha_VELZ( 1:KS-1,i,j) = ATMOS_BOUNDARY_alpha_VELZ(KS,i,j)
700 !!$ ATMOS_BOUNDARY_alpha_VELX( 1:KS-1,i,j) = ATMOS_BOUNDARY_alpha_VELX(KS,i,j)
701 !!$ ATMOS_BOUNDARY_alpha_VELY( 1:KS-1,i,j) = ATMOS_BOUNDARY_alpha_VELY(KS,i,j)
702 !!$ ATMOS_BOUNDARY_alpha_POTT( 1:KS-1,i,j) = ATMOS_BOUNDARY_alpha_POTT(KS,i,j)
703 !!$
704 !!$ ATMOS_BOUNDARY_alpha_DENS(KE+1:KA, i,j) = ATMOS_BOUNDARY_alpha_DENS(KE,i,j)
705 !!$ ATMOS_BOUNDARY_alpha_VELZ(KE+1:KA, i,j) = ATMOS_BOUNDARY_alpha_VELZ(KE,i,j)
706 !!$ ATMOS_BOUNDARY_alpha_VELX(KE+1:KA, i,j) = ATMOS_BOUNDARY_alpha_VELX(KE,i,j)
707 !!$ ATMOS_BOUNDARY_alpha_VELY(KE+1:KA, i,j) = ATMOS_BOUNDARY_alpha_VELY(KE,i,j)
708 !!$ ATMOS_BOUNDARY_alpha_POTT(KE+1:KA, i,j) = ATMOS_BOUNDARY_alpha_POTT(KE,i,j)
709 !!$
710 !!$ do iq = 1, BND_QA
711 !!$ ATMOS_BOUNDARY_alpha_QTRC( 1:KS-1,i,j,iq) = ATMOS_BOUNDARY_alpha_QTRC(KS,i,j,iq)
712 !!$ ATMOS_BOUNDARY_alpha_QTRC(KE+1:KA, i,j,iq) = ATMOS_BOUNDARY_alpha_QTRC(KE,i,j,iq)
713 !!$ end do
714 !!$ enddo
715 !!$ enddo
716 
717  call comm_vars8( atmos_boundary_alpha_dens(:,:,:), 1 )
718  call comm_vars8( atmos_boundary_alpha_velz(:,:,:), 2 )
719  call comm_vars8( atmos_boundary_alpha_velx(:,:,:), 3 )
720  call comm_vars8( atmos_boundary_alpha_vely(:,:,:), 4 )
721  call comm_vars8( atmos_boundary_alpha_pott(:,:,:), 5 )
722  do iq = 1, bnd_qa
723  call comm_vars8( atmos_boundary_alpha_qtrc(:,:,:,iq), 5+iq )
724  end do
725 
726  call comm_wait ( atmos_boundary_alpha_dens(:,:,:), 1, .false. )
727  call comm_wait ( atmos_boundary_alpha_velz(:,:,:), 2, .false. )
728  call comm_wait ( atmos_boundary_alpha_velx(:,:,:), 3, .false. )
729  call comm_wait ( atmos_boundary_alpha_vely(:,:,:), 4, .false. )
730  call comm_wait ( atmos_boundary_alpha_pott(:,:,:), 5, .false. )
731  do iq = 1, bnd_qa
732  call comm_wait ( atmos_boundary_alpha_qtrc(:,:,:,iq), 5+iq, .false. )
733  end do
734 
735  !$acc update host(ATMOS_BOUNDARY_alpha_DENS,ATMOS_BOUNDARY_alpha_VELZ,ATMOS_BOUNDARY_alpha_VELX,ATMOS_BOUNDARY_alpha_VELY,ATMOS_BOUNDARY_alpha_POTT)
736 
737  return
738  end subroutine atmos_boundary_alpha_fillhalo
739 
740  !-----------------------------------------------------------------------------
742  subroutine atmos_boundary_setalpha
743  use scale_const, only: &
744  eps => const_eps, &
745  pi => const_pi
746  use scale_atmos_grid_cartesc, only: &
747  cbfz => atmos_grid_cartesc_cbfz, &
748  cbfx => atmos_grid_cartesc_cbfx, &
749  cbfy => atmos_grid_cartesc_cbfy, &
750  fbfz => atmos_grid_cartesc_fbfz, &
751  fbfx => atmos_grid_cartesc_fbfx, &
753  use scale_atmos_hydrometeor, only: &
754  i_qv
755  real(rp) :: coef_z, alpha_z1, alpha_z2
756  real(rp) :: coef_x, alpha_x1, alpha_x2
757  real(rp) :: coef_y, alpha_y1, alpha_y2
758  real(rp) :: alpha_zm, alpha_xm, alpha_ym
759  real(rp) :: ee1, ee2
760 
761  real(rp) :: alpha_nug ! grid nudging in inner domain
762 
763  integer :: i, j, k, iq
764  !---------------------------------------------------------------------------
765 
766  ! check invalid fraction
767  atmos_boundary_fracz = max( min( atmos_boundary_fracz, 1.0_rp ), eps )
768  atmos_boundary_fracx = max( min( atmos_boundary_fracx, 1.0_rp ), eps )
769  atmos_boundary_fracy = max( min( atmos_boundary_fracy, 1.0_rp ), eps )
770 
771  if ( atmos_boundary_tauz <= 0.0_rp ) then ! invalid tau
772  coef_z = 0.0_rp
773  else
774  coef_z = 1.0_rp / atmos_boundary_tauz
775  endif
776 
777  if ( atmos_boundary_taux <= 0.0_rp ) then ! invalid tau
778  coef_x = 0.0_rp
779  else
780  coef_x = 1.0_rp / atmos_boundary_taux
781  endif
782 
783  if ( atmos_boundary_tauy <= 0.0_rp ) then ! invalid tau
784  coef_y = 0.0_rp
785  else
786  coef_y = 1.0_rp / atmos_boundary_tauy
787  endif
788 
789  if ( atmos_grid_nudging_tau <= 0.0_rp ) then ! invalid tau
790  alpha_nug = 0.0_rp
791  else
792  alpha_nug = 1.0_rp / atmos_grid_nudging_tau
793  endif
794 
795  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
796  !$omp shared(JA,IA,KS,KE,CBFZ,ATMOS_BOUNDARY_FRACZ,FBFZ,ATMOS_BOUNDARY_LINEAR_V,coef_z,CBFX) &
797  !$omp shared(ATMOS_BOUNDARY_FRACX,FBFX,ATMOS_BOUNDARY_LINEAR_H,coef_x) &
798  !$omp shared(ATMOS_BOUNDARY_EXP_H,CBFY,ATMOS_BOUNDARY_FRACY,FBFY,coef_y,l_bnd) &
799  !$omp shared(do_daughter_process) &
800  !$omp shared(ATMOS_BOUNDARY_USE_VELZ,ATMOS_BOUNDARY_alpha_VELZ,ATMOS_BOUNDARY_ALPHAFACT_VELZ) &
801  !$omp shared(ATMOS_BOUNDARY_USE_DENS,ATMOS_BOUNDARY_alpha_DENS,ATMOS_BOUNDARY_ALPHAFACT_DENS) &
802  !$omp shared(ATMOS_BOUNDARY_USE_VELX,ATMOS_BOUNDARY_alpha_VELX,ATMOS_BOUNDARY_ALPHAFACT_VELX) &
803  !$omp shared(ATMOS_BOUNDARY_USE_VELY,ATMOS_BOUNDARY_alpha_VELY,ATMOS_BOUNDARY_ALPHAFACT_VELY) &
804  !$omp shared(ATMOS_BOUNDARY_USE_PT,ATMOS_BOUNDARY_alpha_POTT,ATMOS_BOUNDARY_ALPHAFACT_PT) &
805  !$omp shared(ATMOS_BOUNDARY_USE_QV,ATMOS_BOUNDARY_alpha_QTRC,ATMOS_BOUNDARY_ALPHAFACT_QTRC) &
806  !$omp shared(ATMOS_BOUNDARY_DENS_ADJUST,ATMOS_BOUNDARY_DENS_ADJUST_tau) &
807  !$omp shared(BND_QA,BND_IQ,I_QV) &
808  !$omp shared(ATMOS_GRID_NUDGING_uv,ATMOS_GRID_NUDGING_pt,ATMOS_GRID_NUDGING_qv) &
809  !$omp shared(alpha_nug) &
810  !$omp private(i,j,k,iq) &
811  !$omp private(ee1,ee2,alpha_z1,alpha_z2,alpha_x1,alpha_x2,alpha_y1,alpha_y2,alpha_zm,alpha_xm,alpha_ym)
812  !$acc kernels
813  do j = 1, ja
814  do i = 1, ia
815  do k = ks, ke
816  ee1 = cbfz(k)
817  if ( ee1 <= 1.0_rp - atmos_boundary_fracz ) then
818  ee1 = 0.0_rp
819  else
820  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracz ) / atmos_boundary_fracz
821  endif
822 
823  ee2 = fbfz(k)
824  if ( ee2 <= 1.0_rp - atmos_boundary_fracz ) then
825  ee2 = 0.0_rp
826  else
827  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracz ) / atmos_boundary_fracz
828  endif
829 
830  if ( .not. atmos_boundary_linear_v ) then
831  if ( ee1 > 0.0_rp .AND. ee1 <= 0.5_rp ) then
832  ee1 = 0.5_rp * ( 1.0_rp - cos( ee1*pi ) )
833  elseif( ee1 > 0.5_rp .AND. ee1 <= 1.0_rp ) then
834  ee1 = 0.5_rp * ( 1.0_rp + sin( (ee1-0.5_rp)*pi ) )
835  endif
836  if ( ee2 > 0.0_rp .AND. ee2 <= 0.5_rp ) then
837  ee2 = 0.5_rp * ( 1.0_rp - cos( ee2*pi ) )
838  elseif( ee2 > 0.5_rp .AND. ee2 <= 1.0_rp ) then
839  ee2 = 0.5_rp * ( 1.0_rp + sin( (ee2-0.5_rp)*pi ) )
840  endif
841  endif
842 
843  alpha_z1 = coef_z * ee1
844  alpha_z2 = coef_z * ee2
845  alpha_zm = ee1 / atmos_boundary_dens_adjust_tau
846 
847 
848  ee1 = cbfx(i)
849  if ( ee1 <= 1.0_rp - atmos_boundary_fracx ) then
850  ee1 = 0.0_rp
851  else
852  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracx ) / atmos_boundary_fracx
853  endif
854 
855  ee2 = fbfx(i)
856  if ( ee2 <= 1.0_rp - atmos_boundary_fracx ) then
857  ee2 = 0.0_rp
858  else
859  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracx ) / atmos_boundary_fracx
860  endif
861 
862  if ( .not. atmos_boundary_linear_h ) then
863  ee1 = ee1 * exp( -(1.0_rp-ee1) * atmos_boundary_exp_h )
864  ee2 = ee2 * exp( -(1.0_rp-ee2) * atmos_boundary_exp_h )
865  end if
866 
867  alpha_x1 = coef_x * ee1
868  alpha_x2 = coef_x * ee2
869  alpha_xm = ee1 / atmos_boundary_dens_adjust_tau
870 
871 
872  ee1 = cbfy(j)
873  if ( ee1 <= 1.0_rp - atmos_boundary_fracy ) then
874  ee1 = 0.0_rp
875  else
876  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracy ) / atmos_boundary_fracy
877  endif
878 
879  ee2 = fbfy(j)
880  if ( ee2 <= 1.0_rp - atmos_boundary_fracy ) then
881  ee2 = 0.0_rp
882  else
883  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracy ) / atmos_boundary_fracy
884  endif
885 
886  if ( .not. atmos_boundary_linear_h ) then
887  ee1 = ee1 * exp( -(1.0_rp-ee1) * atmos_boundary_exp_h )
888  ee2 = ee2 * exp( -(1.0_rp-ee2) * atmos_boundary_exp_h )
889  end if
890 
891  alpha_y1 = coef_y * ee1
892  alpha_y2 = coef_y * ee2
893  alpha_ym = ee1 / atmos_boundary_dens_adjust_tau
894 
895 
896  if ( l_bnd ) then
897  if ( atmos_boundary_use_velz ) then
898  atmos_boundary_alpha_velz(k,i,j) = max( alpha_z2, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
899  else
900  atmos_boundary_alpha_velz(k,i,j) = 0.0_rp
901  end if
902  if ( atmos_boundary_dens_adjust ) then
903  atmos_boundary_alpha_dens(k,i,j) = max( alpha_zm, alpha_xm, alpha_ym )
904  else
905  atmos_boundary_alpha_dens(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_dens
906  end if
907  atmos_boundary_alpha_velx(k,i,j) = max( alpha_z1, alpha_x2, alpha_y1 ) * atmos_boundary_alphafact_velx
908  atmos_boundary_alpha_vely(k,i,j) = max( alpha_z1, alpha_x1, alpha_y2 ) * atmos_boundary_alphafact_vely
909  atmos_boundary_alpha_pott(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_pt
910  !$acc loop seq
911  do iq = 1, bnd_qa
912  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
913  end do
914  else
915  if ( atmos_boundary_use_dens ) then
916  atmos_boundary_alpha_dens(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_dens
917  else
918  atmos_boundary_alpha_dens(k,i,j) = 0.0_rp
919  end if
920  if ( atmos_boundary_use_velz ) then
921  atmos_boundary_alpha_velz(k,i,j) = max( alpha_z2, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
922  else
923  atmos_boundary_alpha_velz(k,i,j) = 0.0_rp
924  end if
925  if ( atmos_boundary_use_velx ) then
926  atmos_boundary_alpha_velx(k,i,j) = max( alpha_z1, alpha_x2, alpha_y1 ) * atmos_boundary_alphafact_velx
927  else
928  atmos_boundary_alpha_velx(k,i,j) = 0.0_rp
929  end if
930  if ( atmos_boundary_use_vely ) then
931  atmos_boundary_alpha_vely(k,i,j) = max( alpha_z1, alpha_x1, alpha_y2 ) * atmos_boundary_alphafact_vely
932  else
933  atmos_boundary_alpha_vely(k,i,j) = 0.0_rp
934  end if
935  if ( atmos_boundary_use_pt ) then
936  atmos_boundary_alpha_pott(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_pt
937  else
938  atmos_boundary_alpha_pott(k,i,j) = 0.0_rp
939  end if
940  !$acc loop seq
941  do iq = 1, bnd_qa
942  if ( i_qv > 0 .and. iq == bnd_iq(i_qv) ) then
943  if ( atmos_boundary_use_qv ) then
944  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
945  else
946  atmos_boundary_alpha_qtrc(k,i,j,iq) = 0.0_rp
947  endif
948  else
949  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
950  end if
951  end do
952 
953  end if
954 
955  ! internal grid nudging
956  if ( atmos_grid_nudging_uv ) then
957  atmos_boundary_alpha_velx(k,i,j) = max( atmos_boundary_alpha_velx(k,i,j), alpha_nug )
958  atmos_boundary_alpha_vely(k,i,j) = max( atmos_boundary_alpha_vely(k,i,j), alpha_nug )
959  endif
960  if ( atmos_grid_nudging_pt ) then
961  atmos_boundary_alpha_pott(k,i,j) = max( atmos_boundary_alpha_pott(k,i,j), alpha_nug )
962  endif
963  if ( atmos_grid_nudging_qv ) then
964  atmos_boundary_alpha_qtrc(k,i,j,1) = max( atmos_boundary_alpha_qtrc(k,i,j,1), alpha_nug )
965  endif
966 
967  enddo
968  enddo
969  enddo
970  !$acc end kernels
971 
972  call atmos_boundary_alpha_fillhalo
973 
974  return
975  end subroutine atmos_boundary_setalpha
976 
977  !-----------------------------------------------------------------------------
979  subroutine atmos_boundary_setinitval( &
980  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC )
981  implicit none
982 
983  real(rp), intent(in) :: dens(ka,ia,ja)
984  real(rp), intent(in) :: momz(ka,ia,ja)
985  real(rp), intent(in) :: momx(ka,ia,ja)
986  real(rp), intent(in) :: momy(ka,ia,ja)
987  real(rp), intent(in) :: rhot(ka,ia,ja)
988  real(rp), intent(in) :: qtrc(ka,ia,ja,bnd_qa)
989 
990  integer :: i, j, k, iq, iqb
991  !---------------------------------------------------------------------------
992 
993  !$acc data copyin(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC)
994 
995  !$acc kernels
996  do j = 1, ja
997  do i = 1, ia
998  do k = ks, ke
999  atmos_boundary_dens(k,i,j) = dens(k,i,j)
1000  atmos_boundary_velz(k,i,j) = momz(k,i,j) / ( dens(k,i,j)+dens(k+1,i, j ) ) * 2.0_rp
1001  atmos_boundary_pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
1002  !$acc loop seq
1003  do iq = 1, qa
1004  iqb = bnd_iq(iq)
1005  if ( iqb > 0 ) atmos_boundary_qtrc(k,i,j,iqb) = qtrc(k,i,j,iq)
1006  end do
1007  enddo
1008  enddo
1009  enddo
1010  !$acc end kernels
1011 
1012  !$acc kernels
1013  do j = 1, ja
1014  do i = 1, ia-1
1015  do k = ks, ke
1016  atmos_boundary_velx(k,i,j) = momx(k,i,j) / ( dens(k,i,j)+dens(k, i+1,j ) ) * 2.0_rp
1017  enddo
1018  enddo
1019  enddo
1020  !$acc end kernels
1021  !$acc kernels
1022  do j = 1, ja
1023  do k = ks, ke
1024  atmos_boundary_velx(k,ia,j) = momx(k,ia,j) / dens(k,ia,j)
1025  enddo
1026  enddo
1027  !$acc end kernels
1028 
1029  !$acc kernels
1030  do j = 1, ja-1
1031  do i = 1, ia
1032  do k = ks, ke
1033  atmos_boundary_vely(k,i,j) = momy(k,i,j) / ( dens(k,i,j)+dens(k, i, j+1) ) * 2.0_rp
1034  enddo
1035  enddo
1036  enddo
1037  !$acc end kernels
1038  !$acc kernels
1039  do i = 1, ia
1040  do k = ks, ke
1041  atmos_boundary_vely(k,i,ja) = momy(k,i,ja) / dens(k,i,ja)
1042  enddo
1043  enddo
1044  !$acc end kernels
1045 
1046  call atmos_boundary_var_fillhalo
1047 
1048  !$acc end data
1049 
1050  return
1051  end subroutine atmos_boundary_setinitval
1052 
1053  !-----------------------------------------------------------------------------
1055  subroutine atmos_boundary_read
1056  use scale_prc, only: &
1057  prc_abort
1058  use scale_file_cartesc, only: &
1060  file_cartesc_check_coordinates, &
1061  file_cartesc_read, &
1064  implicit none
1065 
1066  integer :: fid
1067  integer :: iq, iqb
1068  !---------------------------------------------------------------------------
1069 
1070  call file_cartesc_open( atmos_boundary_in_basename, fid, aggregate=atmos_boundary_in_aggregate )
1071 
1072  if ( atmos_boundary_in_check_coordinates ) then
1073  call file_cartesc_check_coordinates( fid, atmos=.true. )
1074  end if
1075 
1076  if ( atmos_boundary_use_dens &
1077  .OR. atmos_boundary_use_velz &
1078  .OR. atmos_boundary_use_velx &
1079  .OR. atmos_boundary_use_vely &
1080  .OR. atmos_boundary_use_pt &
1081  ) then
1082  call file_cartesc_read( fid, 'DENS', 'ZXY', atmos_boundary_dens(:,:,:) )
1083 #ifdef _OPENACC
1084  call file_cartesc_flush(fid)
1085 #endif
1086  !$acc update device(ATMOS_BOUNDARY_DENS)
1087  end if
1088  if ( atmos_boundary_use_dens ) then
1089  call file_cartesc_read( fid, 'ALPHA_DENS', 'ZXY', atmos_boundary_alpha_dens(:,:,:) )
1090 #ifdef _OPENACC
1091  call file_cartesc_flush(fid)
1092 #endif
1093  !$acc update device(ATMOS_BOUNDARY_alpha_DENS)
1094  endif
1095 
1096  if ( atmos_boundary_use_velz ) then
1097  call file_cartesc_read( fid, 'VELZ', 'ZHXY', atmos_boundary_velz(:,:,:) )
1098  call file_cartesc_read( fid, 'ALPHA_VELZ', 'ZHXY', atmos_boundary_alpha_velz(:,:,:) )
1099 #ifdef _OPENACC
1100  call file_cartesc_flush(fid)
1101 #endif
1102  !$acc update device(ATMOS_BOUNDARY_VELZ, ATMOS_BOUNDARY_alpha_VELZ)
1103  endif
1104 
1105  if ( atmos_boundary_use_velx ) then
1106  call file_cartesc_read( fid, 'VELX', 'ZXHY', atmos_boundary_velx(:,:,:) )
1107  call file_cartesc_read( fid, 'ALPHA_VELX', 'ZXHY', atmos_boundary_alpha_velx(:,:,:) )
1108 #ifdef _OPENACC
1109  call file_cartesc_flush(fid)
1110 #endif
1111  !$acc update device(ATMOS_BOUNDARY_VELX, ATMOS_BOUNDARY_alpha_VELX)
1112  endif
1113 
1114  if ( atmos_boundary_use_vely ) then
1115  call file_cartesc_read( fid, 'VELY', 'ZXYH', atmos_boundary_vely(:,:,:) )
1116  call file_cartesc_read( fid, 'ALPHA_VELY', 'ZXYH', atmos_boundary_alpha_vely(:,:,:) )
1117 #ifdef _OPENACC
1118  call file_cartesc_flush(fid)
1119 #endif
1120  !$acc update device(ATMOS_BOUNDARY_VELY, ATMOS_BOUNDARY_alpha_VELY)
1121  endif
1122 
1123  if ( atmos_boundary_use_pt ) then
1124  call file_cartesc_read( fid, 'PT', 'ZXY', atmos_boundary_pott(:,:,:) )
1125  call file_cartesc_read( fid, 'ALPHA_PT', 'ZXY', atmos_boundary_alpha_pott(:,:,:) )
1126 #ifdef _OPENACC
1127  call file_cartesc_flush(fid)
1128 #endif
1129  !$acc update device(ATMOS_BOUNDARY_POTT, ATMOS_BOUNDARY_alpha_POTT)
1130  endif
1131 
1132  do iq = 1, qa
1133  iqb = bnd_iq(iq)
1134  if ( iqb > 0 ) then
1135  call file_cartesc_read( fid, tracer_name(iq), 'ZXY', atmos_boundary_qtrc(:,:,:,iqb) )
1136  call file_cartesc_read( fid, 'ALPHA_'//trim(tracer_name(iq)), 'ZXY', atmos_boundary_alpha_qtrc(:,:,:,iqb) )
1137  endif
1138 #ifdef _OPENACC
1139  call file_cartesc_flush(fid)
1140  !$acc update device(ATMOS_BOUNDARY_QTRC, ATMOS_BOUNDARY_alpha_QTRC)
1141 #endif
1142  end do
1143 
1144  call file_cartesc_close( fid )
1145 
1146 
1147  call atmos_boundary_var_fillhalo
1148  call atmos_boundary_alpha_fillhalo
1149 
1150  return
1151  end subroutine atmos_boundary_read
1152 
1153  !-----------------------------------------------------------------------------
1155  subroutine atmos_boundary_write
1156  use scale_file_cartesc, only: &
1160  file_cartesc_write_var, &
1162  implicit none
1163 
1164  integer :: fid
1165  integer :: vid_dens, vid_a_dens
1166  integer :: vid_velz, vid_a_velz
1167  integer :: vid_velx, vid_a_velx
1168  integer :: vid_vely, vid_a_vely
1169  integer :: vid_pott, vid_a_pott
1170  integer :: vid_qtrc(bnd_qa), vid_a_qtrc(bnd_qa)
1171  integer :: iq, iqb
1172  !---------------------------------------------------------------------------
1173 
1174  call file_cartesc_create( atmos_boundary_out_basename, atmos_boundary_out_title, & ! [IN]
1175  atmos_boundary_out_dtype, & ! [IN]
1176  fid, & ! [OUT]
1177  aggregate = atmos_boundary_out_aggregate ) ! [IN]
1178 
1179  if ( atmos_boundary_use_dens &
1180  .OR. atmos_boundary_use_velz &
1181  .OR. atmos_boundary_use_velx &
1182  .OR. atmos_boundary_use_vely &
1183  .OR. atmos_boundary_use_pt &
1184  .OR. l_bnd ) then
1185  call file_cartesc_def_var( fid, 'DENS', 'Reference Density', 'kg/m3', 'ZXY', atmos_boundary_out_dtype, vid_dens )
1186  else
1187  vid_dens = -1
1188  end if
1189 
1190  if ( atmos_boundary_use_dens .OR. l_bnd ) then
1191  call file_cartesc_def_var( fid, 'ALPHA_DENS', 'Alpha for DENS', '1', 'ZXY', atmos_boundary_out_dtype, vid_a_dens )
1192  else
1193  vid_a_dens = -1
1194  end if
1195 
1196  if ( atmos_boundary_use_velz ) then
1197  call file_cartesc_def_var( fid, 'VELZ', 'Reference Velocity w', 'm/s', 'ZHXY', atmos_boundary_out_dtype, vid_velz )
1198  call file_cartesc_def_var( fid, 'ALPHA_VELZ', 'Alpha for VELZ', '1', 'ZHXY', atmos_boundary_out_dtype, vid_a_velz )
1199  else
1200  vid_velz = -1
1201  end if
1202 
1203  if ( atmos_boundary_use_velx .OR. l_bnd ) then
1204  call file_cartesc_def_var( fid, 'VELX', 'Reference Velocity u', 'm/s', 'ZXHY', atmos_boundary_out_dtype, vid_velx )
1205  call file_cartesc_def_var( fid, 'ALPHA_VELX', 'Alpha for VELX', '1', 'ZXHY', atmos_boundary_out_dtype, vid_a_velx )
1206  else
1207  vid_velx = -1
1208  end if
1209 
1210  if ( atmos_boundary_use_vely .OR. l_bnd ) then
1211  call file_cartesc_def_var( fid, 'VELY', 'Reference Velocity y', 'm/s', 'ZXYH', atmos_boundary_out_dtype, vid_vely )
1212  call file_cartesc_def_var( fid, 'ALPHA_VELY', 'Alpha for VELY', '1', 'ZXYH', atmos_boundary_out_dtype, vid_a_vely )
1213  else
1214  vid_vely = -1
1215  end if
1216 
1217  if ( atmos_boundary_use_pt .OR. l_bnd ) then
1218  call file_cartesc_def_var( fid, 'PT', 'Reference PT', 'K', 'ZXY', atmos_boundary_out_dtype, vid_pott )
1219  call file_cartesc_def_var( fid, 'ALPHA_PT', 'Alpha for PT', '1', 'ZXY', atmos_boundary_out_dtype, vid_a_pott )
1220  else
1221  vid_pott = -1
1222  end if
1223 
1224  do iq = 1, qa
1225  iqb = bnd_iq(iq)
1226  if ( iqb > 0 ) then
1227  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) )
1228  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) )
1229  end if
1230  end do
1231 
1232 
1233  call file_cartesc_enddef( fid )
1234 
1235 
1236  if ( vid_dens > 0 ) then
1237  call file_cartesc_write_var( fid, vid_dens, atmos_boundary_dens(:,:,:), 'DENS', 'ZXY' )
1238  end if
1239 
1240  if ( vid_a_dens > 0 ) then
1241  call file_cartesc_write_var( fid, vid_a_dens, atmos_boundary_alpha_dens(:,:,:), 'ALPHA_DENS', 'ZXY' )
1242  end if
1243 
1244  if ( vid_velz > 0 ) then
1245  call file_cartesc_write_var( fid, vid_velz, atmos_boundary_velz(:,:,:), 'VELZ', 'ZHXY' )
1246  call file_cartesc_write_var( fid, vid_a_velz, atmos_boundary_alpha_velz(:,:,:), 'ALPHA_VELZ', 'ZHXY' )
1247  end if
1248 
1249  if ( vid_velx > 0 ) then
1250  call file_cartesc_write_var( fid, vid_velx, atmos_boundary_velx(:,:,:), 'VELX', 'ZXHY' )
1251  call file_cartesc_write_var( fid, vid_a_velx, atmos_boundary_alpha_velx(:,:,:), 'ALPHA_VELX', 'ZXHY' )
1252  endif
1253 
1254  if ( vid_vely > 0 ) then
1255  call file_cartesc_write_var( fid, vid_vely, atmos_boundary_vely(:,:,:), 'VELY', 'ZXYH' )
1256  call file_cartesc_write_var( fid, vid_a_vely, atmos_boundary_alpha_vely(:,:,:), 'ALPHA_VELY', 'ZXYH' )
1257  end if
1258 
1259  if ( vid_pott > 0 ) then
1260  call file_cartesc_write_var( fid, vid_pott, atmos_boundary_pott(:,:,:), 'PT', 'ZXY' )
1261  call file_cartesc_write_var( fid, vid_a_pott, atmos_boundary_alpha_pott(:,:,:), 'ALPHA_PT', 'ZXY' )
1262  end if
1263 
1264  do iqb = 1, bnd_qa
1265  if ( vid_qtrc(iqb) > 0 ) then
1266  call file_cartesc_write_var( fid, vid_qtrc(iqb), atmos_boundary_qtrc(:,:,:,iqb), tracer_name(iqb), 'ZXY' )
1267  call file_cartesc_write_var( fid, vid_a_qtrc(iqb), atmos_boundary_alpha_qtrc(:,:,:,iqb), 'ALPHA_'//trim(tracer_name(iqb)), 'ZXY' )
1268  end if
1269  end do
1270 
1271 
1272  call file_cartesc_close( fid )
1273 
1274 
1275  return
1276  end subroutine atmos_boundary_write
1277 
1278  !-----------------------------------------------------------------------------
1280  subroutine atmos_boundary_generate
1281  use scale_atmos_refstate, only: &
1283  implicit none
1284 
1285  integer :: i, j, k, iq
1286  !---------------------------------------------------------------------------
1287 
1288  !$omp parallel do collapse(2)
1289  !$acc kernels
1290  do j = 1, ja
1291  do i = 1, ia
1292  do k = ks, ke
1294  atmos_boundary_velz(k,i,j) = atmos_boundary_value_velz
1295  atmos_boundary_velx(k,i,j) = atmos_boundary_value_velx
1296  atmos_boundary_vely(k,i,j) = atmos_boundary_value_vely
1297  atmos_boundary_pott(k,i,j) = atmos_boundary_value_pt
1298  !$acc loop seq
1299  do iq = 1, bnd_qa
1300  atmos_boundary_qtrc(k,i,j,iq) = atmos_boundary_value_qtrc
1301  end do
1302  enddo
1303  enddo
1304  enddo
1305  !$acc end kernels
1306 
1307  call atmos_boundary_var_fillhalo
1308 
1309  return
1310  end subroutine atmos_boundary_generate
1311 
1312  !-----------------------------------------------------------------------------
1314  subroutine atmos_boundary_initialize_file
1315  use scale_prc, only: &
1316  prc_abort
1317  use scale_file_external_input, only: &
1318  file_external_input_regist
1319  use scale_time, only: &
1320  time_dtsec
1321  implicit none
1322 
1323  character(len=4), parameter :: vars_name(4) = (/ 'DENS', 'VELX', 'VELY', 'PT ' /)
1324 
1325  integer :: k, i, j
1326  integer :: n, iq, iqb
1327 
1328  do n = 1, 4
1329  call file_external_input_regist( atmos_boundary_in_basename, & ! [IN]
1330  atmos_boundary_in_basename_add_num, & ! [IN]
1331  atmos_boundary_in_number_of_files, & ! [IN]
1332  vars_name(n), & ! [IN]
1333  'ZXY', & ! [IN]
1334  .false., .false., .false., 0, 0.0_rp, & ! [IN]
1335  check_coordinates = atmos_boundary_in_check_coordinates, & ! [IN]
1336  step_limit = atmos_boundary_in_step_limit, & ! [IN]
1337  update_dt = update_dt ) ! [OUT]
1338  end do
1339 
1340  if ( update_dt <= 0.0_dp ) then
1341  log_error("ATMOS_BOUNDARY_initialize_file",*) 'Time in the file is invalid'
1342  call prc_abort
1343  endif
1344  update_nstep = nint( update_dt / time_dtsec )
1345  if ( abs(update_nstep * time_dtsec - update_dt) > 1e-10_dp ) then
1346  log_error("ATMOS_BOUNDARY_initialize_file",*) 'Time step in the file is not multiple of DT'
1347  call prc_abort
1348  end if
1349 
1350  do iq = 1, qa
1351  iqb = bnd_iq(iq)
1352  if ( iqb > 0 ) &
1353  call file_external_input_regist( atmos_boundary_in_basename, & ! [IN]
1354  atmos_boundary_in_basename_add_num, & ! [IN]
1355  atmos_boundary_in_number_of_files, & ! [IN]
1356  tracer_name(iq), & ! [IN]
1357  'ZXY', & ! [IN]
1358  .false., .false., .false., 0, 0.0_rp, & ! [IN]
1359  check_coordinates = atmos_boundary_in_check_coordinates, & ! [IN]
1360  step_limit = atmos_boundary_in_step_limit ) ! [IN]
1361  end do
1362 
1363  if ( atmos_boundary_use_velz ) then
1364  !$omp parallel do collapse(2)
1365  !$acc kernels
1366  do j = 1, ja
1367  do i = 1, ia
1368  do k = ks, ke
1369  atmos_boundary_velz(k,i,j) = atmos_boundary_value_velz
1370  end do
1371  end do
1372  end do
1373  !$acc end kernels
1374  end if
1375 
1376  return
1377  end subroutine atmos_boundary_initialize_file
1378 
1379  !-----------------------------------------------------------------------------
1381  subroutine atmos_boundary_set_file( &
1382  time )
1384  file_external_input_get_ref, &
1385  i_prev
1386  implicit none
1387  real(DP), intent(in) :: time
1388 
1389  logical :: error
1390  !---------------------------------------------------------------------------
1391 
1392  if ( atmos_boundary_dens_adjust ) then
1393  call file_external_input_get_ref( 'DENS', dens_ref(:,:,:), error, i_prev )
1394  call file_external_input_get_ref( 'VELX', velx_ref(:,:,:), error, i_prev )
1395  call file_external_input_get_ref( 'VELY', vely_ref(:,:,:), error, i_prev )
1396  call calc_mass
1397  end if
1398 
1399  ! read boundary data from input file
1400  call atmos_boundary_update_file( time )
1401 
1402  return
1403  end subroutine atmos_boundary_set_file
1404 
1405  !-----------------------------------------------------------------------------
1407  subroutine atmos_boundary_initialize_online
1408  use scale_prc, only: &
1409  prc_abort
1410  use scale_comm_cartesc_nest, only: &
1415  use scale_time, only: &
1416  time_dtsec, &
1417  time_nstep
1418  implicit none
1419 
1420  update_dt = online_parent_dtsec
1421  update_nstep = nint( update_dt / time_dtsec )
1422  if ( update_nstep * time_dtsec /= update_dt ) then
1423  log_error("ATMOS_BOUNDARY_initialize_online",*) 'DT of the parent is not multiple of the DT'
1424  call prc_abort
1425  end if
1426  if ( update_nstep * online_parent_nstep /= time_nstep ) then
1427  log_error("ATMOS_BOUNDARY_initialize_online",*) 'DURATION must be the same as that of the parent'
1428  call prc_abort
1429  end if
1430 
1431  if ( online_recv_qa > bnd_qa ) then
1432  log_error("ATMOS_BOUNDARY_initialize_online",*) 'ONLINE_RECV_QA exceeds BND_QA'
1433  log_error_cont(*) 'This must not be occur.'
1434  log_error_cont(*) 'Please send your configuration file to SCALE develop member.'
1435  call prc_abort
1436  end if
1437 
1439 
1440  return
1441  end subroutine atmos_boundary_initialize_online
1442 
1443  !-----------------------------------------------------------------------------
1445  subroutine atmos_boundary_set_online( &
1446  time )
1447  use scale_time, only: &
1449  use scale_comm_cartesc_nest, only: &
1450  parent_dtsec => online_parent_dtsec, &
1451  parent_nstep => online_parent_nstep
1452  use scale_file_external_input, only: &
1453  file_external_input_regist
1454  implicit none
1455  real(DP), intent(in) :: time
1456 
1457  integer :: nstep
1458  integer :: iq, iqb
1459  !---------------------------------------------------------------------------
1460 
1461  ! import data from parent domain
1462  call atmos_boundary_update_online_daughter( time, .true., .true. )
1463 
1464  nstep = parent_nstep + 1
1465 
1466  call file_external_input_regist( 'DENS', & ! [IN]
1467  dens_ref(:,:,:), & ! [IN]
1468  'ZXY', & ! [IN]
1469  nstep, & ! [IN]
1470  time, & ! [IN]
1471  parent_dtsec ) ! [IN]
1472 
1473  call file_external_input_regist( 'VELX', & ! [IN]
1474  velx_ref(:,:,:), & ! [IN]
1475  'ZXY', & ! [IN]
1476  nstep, & ! [IN]
1477  time, & ! [IN]
1478  parent_dtsec ) ! [IN]
1479 
1480  call file_external_input_regist( 'VELY', & ! [IN]
1481  vely_ref(:,:,:), & ! [IN]
1482  'ZXY', & ! [IN]
1483  nstep, & ! [IN]
1484  time, & ! [IN]
1485  parent_dtsec ) ! [IN]
1486 
1487  call file_external_input_regist( 'PT', & ! [IN]
1488  pott_ref(:,:,:), & ! [IN]
1489  'ZXY', & ! [IN]
1490  nstep, & ! [IN]
1491  time, & ! [IN]
1492  parent_dtsec ) ! [IN]
1493 
1494  do iq = 1, qa
1495  iqb = bnd_iq(iq)
1496  if ( iqb > 0 ) &
1497  call file_external_input_regist( tracer_name(iq), & ! [IN]
1498  qtrc_ref(:,:,:,iqb), & ! [IN]
1499  'ZXY', & ! [IN]
1500  nstep, & ! [IN]
1501  time, & ! [IN]
1502  parent_dtsec ) ! [IN]
1503  end do
1504 
1505  if ( atmos_boundary_use_velz ) then
1506  call file_external_input_regist( 'VELZ', & ! [IN]
1507  velz_ref(:,:,:), & ! [IN]
1508  'ZXY', & ! [IN]
1509  nstep, & ! [IN]
1510  time, & ! [IN]
1511  parent_dtsec ) ! [IN]
1512 
1513  end if
1514 
1515  call atmos_boundary_update_online_daughter( time, .false., .true. )
1516 
1517  return
1518  end subroutine atmos_boundary_set_online
1519 
1520  !-----------------------------------------------------------------------------
1522  subroutine atmos_boundary_firstsend( &
1523  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, QV, Qe )
1524  use scale_atmos_hydrometeor, only: &
1525  n_hyd
1526  use mod_atmos_phy_mp_vars, only: &
1527  qa_mp
1528  implicit none
1529 
1530  ! arguments
1531  real(RP), intent(in) :: DENS(KA,IA,JA)
1532  real(RP), intent(in) :: MOMZ(KA,IA,JA)
1533  real(RP), intent(in) :: MOMX(KA,IA,JA)
1534  real(RP), intent(in) :: MOMY(KA,IA,JA)
1535  real(RP), intent(in) :: RHOT(KA,IA,JA)
1536  real(RP), intent(in) :: QTRC(KA,IA,JA,QA_MP)
1537  real(RP), intent(in) :: QV (KA,IA,JA)
1538  real(RP), intent(in) :: Qe (KA,IA,JA,N_HYD)
1539  !---------------------------------------------------------------------------
1540 
1541  ! send data at the first time
1542  if ( do_parent_process ) then !online [parent]
1543  ! issue send
1544  call atmos_boundary_send( dens, momz, momx, momy, rhot, qtrc, qv, qe )
1545  endif
1546 
1547  return
1548  end subroutine atmos_boundary_firstsend
1549 
1550  !-----------------------------------------------------------------------------
1556  use scale_file_cartesc, only: &
1558  implicit none
1559 
1560  !---------------------------------------------------------------------------
1561 
1562  if ( do_parent_process ) then !online [parent]
1564  endif
1565 
1566  if ( do_daughter_process ) then !online [daughter]
1568  endif
1569 
1570  !$acc exit data delete(BND_IQ)
1571  deallocate( bnd_iq )
1572 
1573  !$acc exit data delete(ATMOS_BOUNDARY_DENS, ATMOS_BOUNDARY_VELZ, ATMOS_BOUNDARY_VELX, ATMOS_BOUNDARY_VELY, ATMOS_BOUNDARY_POTT, ATMOS_BOUNDARY_QTRC)
1574  deallocate( atmos_boundary_dens )
1575  deallocate( atmos_boundary_velz )
1576  deallocate( atmos_boundary_velx )
1577  deallocate( atmos_boundary_vely )
1578  deallocate( atmos_boundary_pott )
1579  deallocate( atmos_boundary_qtrc )
1580 
1581  !$acc exit data delete(ATMOS_BOUNDARY_alpha_DENS, ATMOS_BOUNDARY_alpha_VELZ, ATMOS_BOUNDARY_alpha_VELX, ATMOS_BOUNDARY_alpha_VELY, ATMOS_BOUNDARY_alpha_POTT, ATMOS_BOUNDARY_alpha_QTRC)
1582  deallocate( atmos_boundary_alpha_dens )
1583  deallocate( atmos_boundary_alpha_velz )
1584  deallocate( atmos_boundary_alpha_velx )
1585  deallocate( atmos_boundary_alpha_vely )
1586  deallocate( atmos_boundary_alpha_pott )
1587  deallocate( atmos_boundary_alpha_qtrc )
1588 
1589  !$acc exit data delete(ATMOS_BOUNDARY_MFLUX_OFFSET_X, ATMOS_BOUNDARY_MFLUX_OFFSET_Y)
1590  deallocate( atmos_boundary_mflux_offset_x )
1591  deallocate( atmos_boundary_mflux_offset_y )
1592 
1593  if ( l_bnd ) then
1594  !$acc exit data delete(DENS_ref, VELX_ref, VELY_ref, POTT_ref, QTRC_ref)
1595  deallocate( dens_ref )
1596  deallocate( velx_ref )
1597  deallocate( vely_ref )
1598  deallocate( pott_ref )
1599  deallocate( qtrc_ref )
1600 
1601  if ( atmos_boundary_use_velz ) then
1602  !$acc exit data delete(VELZ_ref)
1603  deallocate( velz_ref )
1604  end if
1605 
1606  if ( atmos_boundary_dens_adjust ) then
1607  !$acc exit data delete(AREAZUY_W, AREAZUY_E, MFLUX_OFFSET_X, MFLUX_OFFSET_Y, zero_x, zero_y)
1608  deallocate( areazuy_w, areazuy_e )
1609  deallocate( mflux_offset_x, mflux_offset_y )
1610  deallocate( zero_x )
1611  deallocate( zero_y )
1612  end if
1613  end if
1614 
1615  if ( allocated( q_send_work ) ) then
1616  !$acc exit data delete(Q_SEND_WORK)
1617  deallocate( q_send_work )
1618  end if
1619  if ( allocated( q_recv_work ) ) then
1620  !$acc exit data delete(Q_RECV_WORK)
1621  deallocate( q_recv_work )
1622  end if
1623 
1624 
1625  return
1626  end subroutine atmos_boundary_driver_finalize
1627 
1628  !-----------------------------------------------------------------------------
1630  subroutine atmos_boundary_driver_update( &
1631  time )
1632  use scale_prc, only: &
1633  prc_abort
1634  use mod_atmos_phy_mp_vars, only: &
1635  qs_mp
1636  use scale_comm_cartesc_nest, only: &
1638  implicit none
1639  real(dp), intent(in) :: time
1640 
1641  !---------------------------------------------------------------------------
1642 
1643  if ( l_bnd ) then
1644 
1645  if ( do_daughter_process ) then !online [daughter]
1646  call atmos_boundary_update_online_daughter( time, .false., .false. )
1647  else
1648  call atmos_boundary_update_file( time )
1649  end if
1650 
1651  call set_boundary
1652 
1653  elseif ( do_parent_process ) then
1654  ! do nothing
1655  else
1656  log_error("ATMOS_BOUNDARY_driver_update",*) '[BUG] invalid path'
1657  call prc_abort
1658  end if
1659 
1666 
1667  ! To be enable to do asynchronous communicaton
1668  if ( do_daughter_process ) then !online [daughter]
1670  endif
1671 
1672  return
1673  end subroutine atmos_boundary_driver_update
1674 
1675  !-----------------------------------------------------------------------------
1677  subroutine atmos_boundary_driver_send
1680  use mod_atmos_vars, only: &
1681  dens, &
1682  momz, &
1683  momx, &
1684  momy, &
1685  rhot, &
1686  qtrc, &
1687  qv, &
1688  qe
1689  use mod_atmos_phy_mp_vars, only: &
1690  qs_mp, &
1691  qe_mp
1692  implicit none
1693 
1694  !---------------------------------------------------------------------------
1695 
1696  if ( do_parent_process ) then !online [parent]
1697  ! should be called every time step
1698  call atmos_boundary_update_online_parent( dens,momz,momx,momy,rhot,qtrc(:,:,:,qs_mp:qe_mp), qv, qe )
1699 
1700  ! To be enable to do asynchronous communicaton
1702  endif
1703 
1704  return
1705  end subroutine atmos_boundary_driver_send
1706 
1707  !-----------------------------------------------------------------------------
1709  subroutine atmos_boundary_update( &
1710  time )
1711  use scale_prc, only: &
1712  prc_abort
1713  use scale_file_external_input, only: &
1714  file_external_input_update
1715  implicit none
1716  real(dp), intent(in) :: time
1717 
1718  logical :: error
1719  integer :: iq, iqb
1720 
1721  call file_external_input_update( 'DENS', time, atmos_boundary_dens(:,:,:), error )
1722  if ( error ) then
1723  log_error("ATMOS_BOUNDARY_update",*) 'Failed to read data: DENS'
1724  call prc_abort
1725  end if
1726 
1727  call file_external_input_update( 'VELX', time, atmos_boundary_velx(:,:,:), error )
1728  if ( error ) then
1729  log_error("ATMOS_BOUNDARY_update",*) 'Failed to read data: VELX'
1730  call prc_abort
1731  end if
1732 
1733  call file_external_input_update( 'VELY', time, atmos_boundary_vely(:,:,:), error )
1734  if ( error ) then
1735  log_error("ATMOS_BOUNDARY_update",*) 'Failed to read data: VELY'
1736  call prc_abort
1737  end if
1738 
1739  call file_external_input_update( 'PT', time, atmos_boundary_pott(:,:,:), error )
1740  if ( error ) then
1741  log_error("ATMOS_BOUNDARY_update",*) 'Failed to read data: PT'
1742  call prc_abort
1743  end if
1744 
1745  if ( do_daughter_process .and. atmos_boundary_use_velz ) then
1746  call file_external_input_update( 'VELZ', time, atmos_boundary_velz(:,:,:), error )
1747  if ( error ) then
1748  log_error("ATMOS_BOUNDARY_update",*) 'Failed to read data: VELZ'
1749  call prc_abort
1750  end if
1751  end if
1752 
1753  do iq = 1, qa
1754  iqb = bnd_iq(iq)
1755  if ( iqb > 0 ) then
1756  call file_external_input_update( tracer_name(iq), time, atmos_boundary_qtrc(:,:,:,iqb), error )
1757  if ( error ) then
1758  log_error("ATMOS_BOUNDARY_update",*) 'Failed to read data: ', trim(tracer_name(iq))
1759  call prc_abort
1760  end if
1761  end if
1762  end do
1763 
1764  call atmos_boundary_var_fillhalo
1765 
1766  return
1767  end subroutine atmos_boundary_update
1768 
1769  !-----------------------------------------------------------------------------
1771  subroutine atmos_boundary_update_file( &
1772  time )
1773  use scale_file_external_input, only: &
1774  file_external_input_get_ref
1775  implicit none
1776  real(dp), intent(in) :: time
1777 
1778  logical :: error
1779  !---------------------------------------------------------------------------
1780 
1781  call atmos_boundary_update( time )
1782 
1783  if ( atmos_boundary_dens_adjust ) then
1784  call file_external_input_get_ref( 'DENS', dens_ref(:,:,:), error )
1785  call file_external_input_get_ref( 'VELX', velx_ref(:,:,:), error )
1786  call file_external_input_get_ref( 'VELY', vely_ref(:,:,:), error )
1787  call calc_mass
1788  call set_offset
1789  end if
1790 
1791  return
1792  end subroutine atmos_boundary_update_file
1793 
1794  !-----------------------------------------------------------------------------
1796  subroutine atmos_boundary_update_online_parent( &
1797  DENS, & ! [in]
1798  MOMZ, & ! [in]
1799  MOMX, & ! [in]
1800  MOMY, & ! [in]
1801  RHOT, & ! [in]
1802  QTRC, & ! [in]
1803  QV, & ! [in]
1804  Qe ) ! [in]
1805  use scale_comm_cartesc_nest, only: &
1807  use scale_atmos_hydrometeor, only: &
1808  n_hyd
1809  use mod_atmos_phy_mp_vars, only: &
1810  qa_mp
1811  implicit none
1812 
1813  ! arguments
1814  real(rp), intent(in) :: dens(ka,ia,ja)
1815  real(rp), intent(in) :: momz(ka,ia,ja)
1816  real(rp), intent(in) :: momx(ka,ia,ja)
1817  real(rp), intent(in) :: momy(ka,ia,ja)
1818  real(rp), intent(in) :: rhot(ka,ia,ja)
1819  real(rp), intent(in) :: qtrc(ka,ia,ja,qa_mp)
1820  real(rp), intent(in) :: qv (ka,ia,ja)
1821  real(rp), intent(in) :: qe (ka,ia,ja,n_hyd)
1822 
1823  !---------------------------------------------------------------------------
1824 
1825  log_info("ATMOS_BOUNDARY_update_online_parent",*)"ATMOS BOUNDARY update online: PARENT"
1826 
1827  ! issue wait
1829 
1830  ! issue send
1831  call atmos_boundary_send( dens, momz, momx, momy, rhot, qtrc, qv, qe )
1832 
1833  return
1834  end subroutine atmos_boundary_update_online_parent
1835 
1836  !-----------------------------------------------------------------------------
1838  subroutine atmos_boundary_update_online_daughter( &
1839  time, &
1840  init, &
1841  force )
1842  use scale_comm_cartesc_nest, only: &
1844  use scale_file_external_input, only: &
1846  file_external_input_put_ref
1847  implicit none
1848  real(dp), intent(in) :: time
1849  logical, intent(in) :: init
1850  logical, intent(in) :: force
1851 
1852  logical :: do_read
1853  logical :: error
1854  integer :: iq, iqb
1855  !---------------------------------------------------------------------------
1856 
1857  log_info("ATMOS_BOUNDARY_update_online_daughter",'(1X,A,I5)') 'ATMOS BOUNDARY update online: DAUGHTER'
1858 
1859  if ( .not. force ) then
1860  call file_external_input_query( "DENS", & ! [IN]
1861  time, & ! [IN]
1862  do_read ) ! [OUT]
1863  else
1864  do_read = .true.
1865  end if
1866 
1867  if ( do_read ) then
1868  ! issue wait
1869  call atmos_boundary_recv
1870 
1871  ! issue receive
1873 
1874  if ( atmos_boundary_dens_adjust ) then
1875  call calc_mass
1876  call set_offset
1877  end if
1878 
1879  if ( .not. init ) then
1880  call file_external_input_put_ref( 'DENS', dens_ref(:,:,:), error )
1881  call file_external_input_put_ref( 'VELX', velx_ref(:,:,:), error )
1882  call file_external_input_put_ref( 'VELY', vely_ref(:,:,:), error )
1883  call file_external_input_put_ref( 'PT', pott_ref(:,:,:), error )
1884  do iq = 1, qa
1885  iqb = bnd_iq(iq)
1886  if ( iqb > 0 ) then
1887  call file_external_input_put_ref( tracer_name(iq), qtrc_ref(:,:,:,iqb), error )
1888  end if
1889  end do
1890 
1891  if ( atmos_boundary_use_velz ) then
1892  call file_external_input_put_ref( 'VELZ', velz_ref(:,:,:), error )
1893  end if
1894  end if
1895 
1896  end if
1897 
1898  if ( .not. init ) call atmos_boundary_update( time )
1899 
1900 
1901  return
1902  end subroutine atmos_boundary_update_online_daughter
1903 
1904  !-----------------------------------------------------------------------------
1906  subroutine atmos_boundary_send( &
1907  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, QV, Qe )
1908  use scale_comm_cartesc_nest, only: &
1910  online_send_qa, &
1912  use scale_atmos_hydrometeor, only: &
1913  n_hyd
1914  use mod_atmos_phy_mp_vars, only: &
1915  qa_mp
1916  implicit none
1917 
1918  ! arguments
1919  real(rp), intent(in) :: dens(ka,ia,ja)
1920  real(rp), intent(in) :: momz(ka,ia,ja)
1921  real(rp), intent(in) :: momx(ka,ia,ja)
1922  real(rp), intent(in) :: momy(ka,ia,ja)
1923  real(rp), intent(in) :: rhot(ka,ia,ja)
1924  real(rp), intent(in), target :: qtrc(ka,ia,ja,qa_mp)
1925  real(rp), intent(in) :: qv (ka,ia,ja)
1926  real(rp), intent(in) :: qe (ka,ia,ja,n_hyd)
1927 
1928  ! works
1929  real(rp), pointer :: q(:,:,:,:)
1930 
1931  integer :: iq
1932  !---------------------------------------------------------------------------
1933 
1934  !$acc data copyin(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, QV, Qe)
1935 
1936  if ( online_send_diagqhyd ) then
1937  q => q_send_work
1938  !$acc kernels
1939  q(:,:,:,1) = qv(:,:,:)
1940  !$acc end kernels
1941  !$acc kernels
1942  do iq = 2, online_send_qa
1943  q(:,:,:,iq) = qe(:,:,:,iq-1)
1944  end do
1945  !$acc end kernels
1946  else
1947  q => qtrc
1948  end if
1949 
1950  call comm_cartesc_nest_nestdown_send( dens(:,:,:), & !(KA,IA,JA)
1951  momz(:,:,:), & !(KA,IA,JA)
1952  momx(:,:,:), & !(KA,IA,JA)
1953  momy(:,:,:), & !(KA,IA,JA)
1954  rhot(:,:,:), & !(KA,IA,JA)
1955  q(:,:,:,:) ) !(KA,IA,JA,NESTQA)
1956 
1957  !$acc end data
1958 
1959  return
1960  end subroutine atmos_boundary_send
1961 
1962  !-----------------------------------------------------------------------------
1964  subroutine atmos_boundary_recv
1965  use scale_prc, only: &
1966  prc_abort
1967  use scale_comm_cartesc_nest, only: &
1969  online_recv_qa, &
1971  use mod_atmos_phy_mp_driver, only: &
1973  implicit none
1974 
1975  ! works
1976  real(rp), pointer :: q(:,:,:,:)
1977  !---------------------------------------------------------------------------
1978 
1979  if ( online_recv_diagqhyd ) then
1980  q => q_recv_work
1981  else
1982  q => qtrc_ref(:,:,:,1:online_recv_qa)
1983  end if
1984 
1985  call comm_cartesc_nest_nestdown_recv( dens_ref(:,:,:), & !(KA,IA,JA)
1986  velz_ref(:,:,:), & !(KA,IA,JA)
1987  velx_ref(:,:,:), & !(KA,IA,JA)
1988  vely_ref(:,:,:), & !(KA,IA,JA)
1989  pott_ref(:,:,:), & !(KA,IA,JA)
1990  q(:,:,:,:) ) !(KA,IA,JA,NESTQA)
1991 
1992  if ( online_recv_diagqhyd ) then
1993  call atmos_phy_mp_driver_qhyd2qtrc( ka, ks, ke, ia, 1, ia, ja, 1, ja, &
1994  q(:,:,:,1), q(:,:,:,2:), & ! (in)
1995  qtrc_ref(:,:,:,:) ) ! (out)
1996  end if
1997 
1998  return
1999  end subroutine atmos_boundary_recv
2000 
2001  !-----------------------------------------------------------------------------
2002  subroutine set_boundary
2003  use scale_prc_cartesc, only: &
2004  prc_has_w, &
2005  prc_has_e, &
2006  prc_has_s, &
2007  prc_has_n
2008  use mod_atmos_vars, only: &
2009  dens, &
2010  momz, &
2011  momx, &
2012  momy, &
2013  rhot, &
2014  qtrc
2015  implicit none
2016 
2017  integer :: i, j, k, iq, iqb
2018 
2019  !$acc data copy(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC)
2020 
2021  ! fill HALO in western region
2022  if ( .NOT. prc_has_w ) then
2023  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2024  !$omp shared(JA,JS,IS,KS,KE,QA,DENS,MOMX,RHOT,QTRC) &
2025  !$omp shared(ATMOS_BOUNDARY_DENS,ATMOS_BOUNDARY_VELX) &
2026  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC) &
2027  !$omp shared(BND_QA,BND_IQ) &
2028  !$omp private(i,j,k,iq,iqb)
2029  !$acc kernels
2030  !$acc loop collapse(3) independent
2031  do j = 1, ja
2032  do i = 1, is-1
2033  do k = ks, ke
2034  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2035  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2036  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2037  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2038  !$acc loop seq
2039  do iq = 1, qa
2040  iqb = bnd_iq(iq)
2041  if ( iqb > 0 ) then
2042  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2043  else
2044  qtrc(k,i,j,iq) = qtrc(k,is,j,iq)
2045  end if
2046  end do
2047  end do
2048  end do
2049  end do
2050  !$acc end kernels
2051  !$omp parallel do
2052  !$acc kernels
2053  do j = 1, ja-1
2054  do i = 1, is-1
2055  do k = ks, ke
2056  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2057  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2058  end do
2059  end do
2060  end do
2061  !$acc end kernels
2062  !$acc kernels
2063  do i = 1, is-1
2064  do k = ks, ke
2065  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) &
2066  * atmos_boundary_dens(k,i,ja)
2067  end do
2068  end do
2069  !$acc end kernels
2070  if ( atmos_boundary_use_velz ) then
2071  !$omp parallel do
2072  !$acc kernels
2073  do j = 1, ja
2074  do i = 1, is-1
2075  do k = ks, ke-1
2076  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2077  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2078  end do
2079  end do
2080  end do
2081  !$acc end kernels
2082  else
2083  !$omp parallel do
2084  !$acc kernels
2085  do j = 1, ja
2086  !$acc loop independent
2087  do i = 1, is-1
2088  do k = ks, ke-1
2089  momz(k,i,j) = momz(k,is,j)
2090  end do
2091  end do
2092  end do
2093  !$acc end kernels
2094  end if
2095  end if
2096 
2097  ! fill HALO in eastern region
2098  if ( .NOT. prc_has_e ) then
2099  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2100  !$omp shared(JA,IE,IA,KS,KE,QA) &
2101  !$omp shared(DENS,RHOT,QTRC) &
2102  !$omp shared(ATMOS_BOUNDARY_DENS,ATMOS_BOUNDARY_VELX) &
2103  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC) &
2104  !$omp shared(BND_QA,BND_IQ) &
2105  !$omp private(i,j,k,iq,iqb)
2106  !$acc kernels
2107  !$acc loop collapse(3) independent
2108  do j = 1, ja
2109  do i = ie+1, ia
2110  do k = ks, ke
2111  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2112  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2113  !$acc loop seq
2114  do iq = 1, qa
2115  iqb = bnd_iq(iq)
2116  if ( iqb > 0 ) then
2117  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2118  else
2119  qtrc(k,i,j,iq) = qtrc(k,ie,j,iq)
2120  end if
2121  end do
2122  end do
2123  end do
2124  end do
2125  !$acc end kernels
2126  !$omp parallel do
2127  !$acc kernels
2128  do j = 1, ja
2129  do i = ie, ia-1
2130  do k = ks, ke
2131  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2132  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2133  end do
2134  end do
2135  end do
2136  !$acc end kernels
2137  !$omp parallel do
2138  !$acc kernels
2139  do j = 1, ja
2140  do k = ks, ke
2141  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) * atmos_boundary_dens(k,ia,j)
2142  end do
2143  end do
2144  !$acc end kernels
2145  !$omp parallel do
2146  !$acc kernels
2147  do j = 1, ja-1
2148  do i = ie+1, ia
2149  do k = ks, ke
2150  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2151  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2152  end do
2153  end do
2154  end do
2155  !$acc end kernels
2156  !$acc kernels
2157  do i = ie+1, ia
2158  do k = ks, ke
2159  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) &
2160  * atmos_boundary_dens(k,i,ja)
2161  end do
2162  end do
2163  !$acc end kernels
2164  if ( atmos_boundary_use_velz ) then
2165  !$omp parallel do
2166  !$acc kernels
2167  do j = 1, ja
2168  do i = ie+1, ia
2169  do k = ks, ke-1
2170  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2171  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2172  end do
2173  end do
2174  end do
2175  !$acc end kernels
2176  else
2177  !$omp parallel do
2178  !$acc kernels
2179  do j = 1, ja
2180  do i = ie+1, ia
2181  do k = ks, ke-1
2182  momz(k,i,j) = momz(k,ie,j)
2183  end do
2184  end do
2185  end do
2186  !$acc end kernels
2187  end if
2188  end if
2189 
2190  ! fill HALO in southern region
2191  if ( .NOT. prc_has_s ) then
2192  !$acc kernels
2193  !$acc loop collapse(3) independent
2194  do j = 1, js-1
2195  do i = 1, ia
2196  do k = ks, ke
2197  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2198  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2199  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2200  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2201  !$acc loop seq
2202  do iq = 1, qa
2203  iqb = bnd_iq(iq)
2204  if ( iqb > 0 ) then
2205  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2206  else
2207  qtrc(k,i,j,iq) = qtrc(k,i,js,iq)
2208  end if
2209  end do
2210  end do
2211  end do
2212  end do
2213  !$acc end kernels
2214  !$acc kernels
2215  do j = 1, js-1
2216  do i = 1, ia-1
2217  do k = ks, ke
2218  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2219  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2220  end do
2221  end do
2222  end do
2223  !$acc end kernels
2224  !$acc kernels
2225  do j = 1, js-1
2226  do k = ks, ke
2227  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) &
2228  * atmos_boundary_dens(k,ia,j)
2229  end do
2230  end do
2231  !$acc end kernels
2232  if ( atmos_boundary_use_velz ) then
2233  !$acc kernels
2234  do j = 1, js-1
2235  do i = 1, ia
2236  do k = ks, ke-1
2237  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2238  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2239  end do
2240  end do
2241  end do
2242  !$acc end kernels
2243  else
2244  !$acc kernels
2245  !$acc loop independent
2246  do j = 1, js-1
2247  do i = 1, ia
2248  do k = ks, ke-1
2249  momz(k,i,j) = momz(k,i,js)
2250  end do
2251  end do
2252  end do
2253  !$acc end kernels
2254  end if
2255  end if
2256 
2257  ! fill HALO in northern region
2258  if ( .NOT. prc_has_n ) then
2259  !$acc kernels
2260  !$acc loop collapse(3) independent
2261  do j = je+1, ja
2262  do i = 1, ia
2263  do k = ks, ke
2264  dens(k,i,j) = atmos_boundary_dens(k,i,j)
2265  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
2266  !$acc loop seq
2267  do iq = 1, qa
2268  iqb = bnd_iq(iq)
2269  if ( iqb > 0 ) then
2270  qtrc(k,i,j,iq) = atmos_boundary_qtrc(k,i,j,iqb)
2271  else
2272  qtrc(k,i,j,iq) = qtrc(k,i,je,iq)
2273  end if
2274  end do
2275  end do
2276  end do
2277  end do
2278  !$acc end kernels
2279  !$acc kernels
2280  do j = je, ja-1
2281  do i = 1, ia
2282  do k = ks, ke
2283  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
2284  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
2285  end do
2286  end do
2287  end do
2288  !$acc end kernels
2289  !$acc kernels
2290  do i = 1, ia
2291  do k = ks, ke
2292  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) * atmos_boundary_dens(k,i,ja)
2293  end do
2294  end do
2295  !$acc end kernels
2296  !$acc kernels
2297  do j = je+1, ja
2298  do i = 1, ia-1
2299  do k = ks, ke
2300  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
2301  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
2302  end do
2303  end do
2304  end do
2305  !$acc end kernels
2306  !$acc kernels
2307  do j = je+1, ja
2308  do k = ks, ke
2309  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) &
2310  * atmos_boundary_dens(k,ia,j)
2311  end do
2312  end do
2313  !$acc end kernels
2314  if ( atmos_boundary_use_velz ) then
2315  !$acc kernels
2316  do j = je+1, ja
2317  do i = 1, ia
2318  do k = ks, ke-1
2319  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
2320  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
2321  end do
2322  end do
2323  end do
2324  !$acc end kernels
2325  else
2326  !$acc kernels
2327  do j = je+1, ja
2328  do i = 1, ia
2329  do k = ks, ke-1
2330  momz(k,i,j) = momz(k,i,je)
2331  end do
2332  end do
2333  end do
2334  !$acc end kernels
2335  end if
2336  end if
2337 
2338  !$acc end data
2339 
2340  return
2341  end subroutine set_boundary
2342 
2343  subroutine history_bnd( &
2344  ATMOS_BOUNDARY_DENS, &
2345  ATMOS_BOUNDARY_VELZ, &
2346  ATMOS_BOUNDARY_VELX, &
2347  ATMOS_BOUNDARY_VELY, &
2348  ATMOS_BOUNDARY_POTT, &
2349  ATMOS_BOUNDARY_QTRC )
2350  use scale_file_history, only: &
2351  file_history_in
2352  use mod_atmos_phy_mp_vars, only: &
2353  qs_mp
2354  implicit none
2355  real(RP), intent(in) :: ATMOS_BOUNDARY_DENS(KA,IA,JA)
2356  real(RP), intent(in) :: ATMOS_BOUNDARY_VELZ(KA,IA,JA)
2357  real(RP), intent(in) :: ATMOS_BOUNDARY_VELX(KA,IA,JA)
2358  real(RP), intent(in) :: ATMOS_BOUNDARY_VELY(KA,IA,JA)
2359  real(RP), intent(in) :: ATMOS_BOUNDARY_POTT(KA,IA,JA)
2360  real(RP), intent(in) :: ATMOS_BOUNDARY_QTRC(KA,IA,JA,BND_QA)
2361 
2362  integer :: iq, iqb
2363 
2364  call file_history_in( atmos_boundary_dens(:,:,:), 'DENS_BND', 'Boundary Density', 'kg/m3' )
2365  call file_history_in( atmos_boundary_velz(:,:,:), 'VELZ_BND', 'Boundary velocity z-direction', 'm/s', dim_type='ZHXY' )
2366  call file_history_in( atmos_boundary_velx(:,:,:), 'VELX_BND', 'Boundary velocity x-direction', 'm/s', dim_type='ZXHY' )
2367  call file_history_in( atmos_boundary_vely(:,:,:), 'VELY_BND', 'Boundary velocity y-direction', 'm/s', dim_type='ZXYH' )
2368  call file_history_in( atmos_boundary_pott(:,:,:), 'PT_BND', 'Boundary potential temperature', 'K' )
2369  do iq = 1, qa
2370  iqb = bnd_iq(iq)
2371  if ( iqb > 0 ) then
2372  call file_history_in( atmos_boundary_qtrc(:,:,:,iqb), trim(tracer_name(iq))//'_BND', &
2373  trim(tracer_name(iq))//' in boundary', 'kg/kg' )
2374  end if
2375  enddo
2376 
2377  return
2378  end subroutine history_bnd
2379 
2380  subroutine calc_mass
2381  use scale_prc_cartesc, only: &
2382  prc_has_w, &
2383  prc_has_e, &
2384  prc_has_s, &
2385  prc_has_n
2386  use scale_statistics, only: &
2387  statistics_total
2388  use scale_atmos_grid_cartesc_real, only: &
2389  vol => atmos_grid_cartesc_real_vol, &
2390  totvol => atmos_grid_cartesc_real_totvol, &
2391  areazxv_y => atmos_grid_cartesc_real_areazxv_y, &
2392  totareazuy_x => atmos_grid_cartesc_real_totareazuy_x, &
2393  totareazxv_y => atmos_grid_cartesc_real_totareazxv_y
2394  use scale_atmos_refstate, only: &
2395  dens_ref => atmos_refstate_dens
2396  use mod_atmos_vars, only: &
2397  dens
2398  implicit none
2399 
2400  real(DP) :: masstot, masstot_current
2401  real(DP) :: massflx
2402  real(DP) :: offset_band, offset_bias
2403  real(DP) :: ref_tot
2404  real(DP) :: flx_w, flx_e, flx_s, flx_n
2405  real(DP) :: ref_w, ref_e, ref_s, ref_n
2406 
2407  real(RP), target :: work_x(KA,JA), work_y(KA,IA)
2408  real(RP), pointer :: ptr(:,:)
2409 
2410  integer :: k, i, j
2411 
2412  !$acc data create(work_x, work_y)
2413 
2414  ! total mass
2415  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
2416  dens_ref(:,:,:), & ! (in)
2417  "DENS_bnd", & ! (in)
2418  vol(:,:,:), totvol, & ! (in)
2419  log_suppress = .true., global = .true., & ! (in)
2420  sum = masstot ) ! (out)
2421 
2422 !!$ call STATISTICS_total( KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2423 !!$ DENS(:,:,:), "DENS_bnd_update", & ! (in)
2424 !!$ VOL(:,:,:), TOTVOL, & ! (in)
2425 !!$ log_suppress = .true., global = .true., & ! (in)
2426 !!$ sum = masstot_current ) ! (out)
2427 
2428 
2429  ! West
2430  if ( .NOT. prc_has_w ) then
2431  !$omp parallel do
2432  !$acc kernels
2433  do j = js, je
2434  do k = ks, ke
2435  work_x(k,j) = velx_ref(k,is-1,j) &
2436  * ( dens_ref(k,is-1,j) + dens_ref(k,is,j) ) * 0.5_rp
2437  end do
2438  end do
2439  !$acc end kernels
2440  ptr => work_x
2441  else
2442  ptr => zero_x
2443  end if
2444  call statistics_total( ka, ks, ke, ja, js, je, &
2445  ptr(:,:), "MFLUX_bnd_w", & ! (in)
2446  areazuy_w(:,:), totareazuy_x(is-1), & ! (in)
2447  log_suppress = .true., global = .true., & ! (in)
2448  sum = flx_w ) ! (out)
2449  if ( .NOT. prc_has_w ) then
2450  !$omp parallel do
2451  !$acc kernels
2452  do j = js, je
2453  do k = ks, ke
2454  work_x(k,j) = dens_ref(k,is,j)
2455  end do
2456  end do
2457  !$acc end kernels
2458  end if
2459  call statistics_total( ka, ks, ke, ja, js, je, &
2460  ptr(:,:), "DENS_ref_w", & ! (in)
2461  areazuy_w(:,:), totareazuy_x(is-1), & ! (in)
2462  log_suppress = .true., global = .true., & ! (in)
2463  sum = ref_w ) ! (out)
2464 
2465  ! East
2466  if ( .NOT. prc_has_e ) then
2467  !$omp parallel do
2468  !$acc kernels
2469  do j = js, je
2470  do k = ks, ke
2471  work_x(k,j) = velx_ref(k,ie,j) &
2472  * ( dens_ref(k,ie,j) + dens_ref(k,ie+1,j) ) * 0.5_rp
2473  end do
2474  end do
2475  !$acc end kernels
2476  ptr => work_x
2477  else
2478  ptr => zero_x
2479  end if
2480  call statistics_total( ka, ks, ke, ja, js, je, &
2481  ptr(:,:), "MFLUX_bnd_e", & ! (in)
2482  areazuy_e(:,:), totareazuy_x(ie), & ! (in)
2483  log_suppress = .true., global = .true., & ! (in)
2484  sum = flx_e ) ! (out)
2485  if ( .NOT. prc_has_e ) then
2486  !$omp parallel do
2487  !$acc kernels
2488  do j = js, je
2489  do k = ks, ke
2490  work_x(k,j) = dens_ref(k,ie,j)
2491  end do
2492  end do
2493  !$acc end kernels
2494  end if
2495  call statistics_total( ka, ks, ke, ja, js, je, &
2496  ptr(:,:), "DENS_ref_e", & ! (in)
2497  areazuy_e(:,:), totareazuy_x(ie), & ! (in)
2498  log_suppress = .true., global = .true., & ! (in)
2499  sum = ref_e ) ! (out)
2500 
2501  ! South
2502  if ( .NOT. prc_has_s ) then
2503  !$omp parallel do
2504  !$acc kernels
2505  do i = is, ie
2506  do k = ks, ke
2507  work_y(k,i) = vely_ref(k,i,js-1) &
2508  * ( dens_ref(k,i,js-1) + dens_ref(k,i,js) ) * 0.5_rp
2509  end do
2510  end do
2511  !$acc end kernels
2512  ptr => work_y
2513  else
2514  ptr => zero_y
2515  end if
2516  call statistics_total( ka, ks, ke, ia, is, ie, &
2517  ptr(:,:), "MFLUX_bnd_s", & ! (in)
2518  areazxv_y(:,:,js-1), totareazxv_y(js-1), & ! (in)
2519  log_suppress = .true., global = .true., & ! (in)
2520  sum = flx_s ) ! (out)
2521  if ( .NOT. prc_has_s ) then
2522  !$omp parallel do
2523  !$acc kernels
2524  do i = is, ie
2525  do k = ks, ke
2526  work_y(k,i) = dens_ref(k,i,js)
2527  end do
2528  end do
2529  !$acc end kernels
2530  end if
2531  call statistics_total( ka, ks, ke, ia, is, ie, &
2532  ptr(:,:), "DENS_ref_s", & ! (in)
2533  areazxv_y(:,:,js-1), totareazxv_y(js-1), & ! (in)
2534  log_suppress = .true., global = .true., & ! (in)
2535  sum = ref_s ) ! (out)
2536 
2537  ! North
2538  if ( .NOT. prc_has_n ) then
2539  !$omp parallel do
2540  !$acc kernels
2541  do i = is, ie
2542  do k = ks, ke
2543  work_y(k,i) = vely_ref(k,i,je) &
2544  * ( dens_ref(k,i,je) + dens_ref(k,i,je+1) ) * 0.5_rp
2545  end do
2546  end do
2547  !$acc end kernels
2548  ptr => work_y
2549  else
2550  ptr => zero_y
2551  end if
2552  call statistics_total( ka, ks, ke, ia, is, ie, &
2553  ptr(:,:), "MFLX_bnd_n", & ! (in)
2554  areazxv_y(:,:,je), totareazxv_y(je), & ! (in)
2555  log_suppress = .true., global = .true., & ! (in)
2556  sum = flx_n ) ! (out)
2557  if ( .NOT. prc_has_n ) then
2558  !$omp parallel do
2559  !$acc kernels
2560  do i = is, ie
2561  do k = ks, ke
2562  work_y(k,i) = dens_ref(k,i,je)
2563  end do
2564  end do
2565  !$acc end kernels
2566  end if
2567  call statistics_total( ka, ks, ke, ia, is, ie, &
2568  ptr(:,:), "DENS_ref_n", & ! (in)
2569  areazxv_y(:,:,je), totareazxv_y(je), & ! (in)
2570  log_suppress = .true., global = .true., & ! (in)
2571  sum = ref_n ) ! (out)
2572 
2573  massflx = flx_w - flx_e + flx_s - flx_n
2574 
2575  offset_band = ( masstot - masstot_now ) / update_dt &
2576  - ( massflx + massflx_now ) * 0.5_dp
2577 ! offset_bias = ( MASSTOT_now - masstot_current ) / UPDATE_DT
2578 
2579  log_info("ATMOS_BOUNDARY_calc_mass",*) "Offset_band is: ", offset_band, "(", masstot, masstot_now, massflx, massflx_now, ")"
2580 
2581  ref_tot = ref_w + ref_e + ref_s + ref_n
2582  offset_band = offset_band / ref_tot
2583 ! offset_bias = offset_bias / ref_tot
2584 
2585  log_info_cont(*) " per dens ", offset_band
2586 
2587  ! density of the reference state is used as weight
2588  if ( .not. prc_has_w ) then
2589  !$omp parallel do
2590  !$acc kernels
2591  do j = js, je
2592  do k = ks, ke
2593  mflux_offset_x(k,j,1,1) = offset_band * dens_ref(k,is,j)
2594 ! MFLUX_OFFSET_X(k,j,1,2) = offset_bias * DENS_ref(k,IS,j)
2595  end do
2596  end do
2597  !$acc end kernels
2598  end if
2599  if ( .not. prc_has_e ) then
2600  !$omp parallel do
2601  !$acc kernels
2602  do j = js, je
2603  do k = ks, ke
2604  mflux_offset_x(k,j,2,1) = - offset_band * dens_ref(k,ie,j)
2605 ! MFLUX_OFFSET_X(k,j,2,2) = - offset_bias * DENS_ref(k,IE,j)
2606  end do
2607  end do
2608  !$acc end kernels
2609  end if
2610  if ( .not. prc_has_s ) then
2611  !$omp parallel do
2612  !$acc kernels
2613  do i = is, ie
2614  do k = ks, ke
2615  mflux_offset_y(k,i,1,1) = offset_band * dens_ref(k,i,js)
2616 ! MFLUX_OFFSET_Y(k,i,1,2) = offset_bias * DENS_ref(k,i,JS)
2617  end do
2618  end do
2619  !$acc end kernels
2620  end if
2621  if ( .not. prc_has_n ) then
2622  !$omp parallel do
2623  !$acc kernels
2624  do i = is, ie
2625  do k = ks, ke
2626  mflux_offset_y(k,i,2,1) = - offset_band * dens_ref(k,i,je)
2627 ! MFLUX_OFFSET_Y(k,i,2,2) = - offset_bias * DENS_ref(k,i,JE)
2628  end do
2629  end do
2630  !$acc end kernels
2631  end if
2632 
2633  masstot_now = masstot
2634  massflx_now = massflx
2635 
2636  !$acc end data
2637 
2638  return
2639  end subroutine calc_mass
2640 
2641 
2642  subroutine set_offset
2643  integer :: k, i, j, n
2644 
2645  !$omp parallel do
2646  !$acc kernels
2647  do j = js, je
2648  do k = ks, ke
2649  do n = 1, 2
2650  atmos_boundary_mflux_offset_x(k,j,n) = mflux_offset_x(k,j,n,1) * offset_time_fact(now_step) !&
2651 ! + MFLUX_OFFSET_X(k,j,n,2)
2652  end do
2653  end do
2654  end do
2655  !$acc end kernels
2656 
2657  !$omp parallel do
2658  !$acc kernels
2659  do i = is, ie
2660  do k = ks, ke
2661  do n = 1, 2
2662  atmos_boundary_mflux_offset_y(k,i,n) = mflux_offset_y(k,i,n,1) * offset_time_fact(now_step) !&
2663 ! + MFLUX_OFFSET_Y(k,i,n,2)
2664  end do
2665  end do
2666  end do
2667  !$acc end kernels
2668 
2669  return
2670  end subroutine set_offset
2671 
2672 end module mod_atmos_bnd_driver
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:105
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_time::time_nowdaysec
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:72
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:350
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:78
scale_file_external_input::file_external_input_query
subroutine, public file_external_input_query(varname, time_current, do_readdata)
Check time to read.
Definition: scale_file_external_input.F90:2027
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
scale_tracer::tracer_unit
character(len=h_short), dimension(qa_max), public tracer_unit
Definition: scale_tracer.F90:41
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:51
scale_file_cartesc::file_cartesc_enddef
subroutine, public file_cartesc_enddef(fid)
Exit netCDF file define mode.
Definition: scale_file_cartesC.F90:964
mod_atmos_bnd_driver::atmos_boundary_set_online
subroutine atmos_boundary_set_online(time)
Set boundary value for real case experiment [online daughter].
Definition: mod_atmos_bnd_driver.F90:1447
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:3360
mod_atmos_bnd_driver::atmos_boundary_driver_finalize
subroutine, public atmos_boundary_driver_finalize
Finalize boundary value.
Definition: mod_atmos_bnd_driver.F90:1553
scale_atmos_refstate::atmos_refstate_dens
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_dens
refernce density [kg/m3]
Definition: scale_atmos_refstate.F90:39
scale_comm_cartesc_nest::comm_cartesc_nest_test_recv
subroutine, public comm_cartesc_nest_test_recv
[check communication status] Inter-communication (daughter side)
Definition: scale_comm_cartesC_nest.F90:3107
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:79
scale_comm_cartesc_nest
module Communication CartesianC nesting
Definition: scale_comm_cartesC_nest.F90:12
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:1678
mod_atmos_bnd_driver::atmos_boundary_driver_update
subroutine, public atmos_boundary_driver_update(time)
Update boundary value with a constant time boundary.
Definition: mod_atmos_bnd_driver.F90:1632
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_comm_cartesc_nest::online_recv_diagqhyd
logical, public online_recv_diagqhyd
Definition: scale_comm_cartesC_nest.F90:87
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
scale_comm_cartesc_nest::online_send_qa
integer, public online_send_qa
number of tracer sent to the daughter domain
Definition: scale_comm_cartesC_nest.F90:90
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
scale_comm_cartesc_nest::online_parent_nstep
integer, public online_parent_nstep
parent nsteps
Definition: scale_comm_cartesC_nest.F90:93
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:73
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:114
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:80
mod_atmos_bnd_driver::atmos_boundary_driver_setup
subroutine, public atmos_boundary_driver_setup
Setup.
Definition: mod_atmos_bnd_driver.F90:179
mod_atmos_phy_ch_vars::qs_ch
integer, public qs_ch
Definition: mod_atmos_phy_ch_vars.F90:62
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:2643
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:49
scale_comm_cartesc_nest::comm_cartesc_nest_nestdown_recv
subroutine, public comm_cartesc_nest_nestdown_recv(DENS_recv, VELZ_recv, VELX_recv, VELY_recv, POTT_recv, QTRC_recv)
Boundary data transfer from parent to daughter: nestdown (daughter side)
Definition: scale_comm_cartesC_nest.F90:2237
mod_atmos_phy_mp_driver
module atmosphere / physics / cloud microphysics
Definition: mod_atmos_phy_mp_driver.F90:12
scale_comm_cartesc_nest::comm_cartesc_nest_nestdown_send
subroutine, public comm_cartesc_nest_nestdown_send(DENS_send, MOMZ_send, MOMX_send, MOMY_send, RHOT_send, QTRC_send)
Boundary data transfer from parent to daughter: nestdown (parent side)
Definition: scale_comm_cartesC_nest.F90:2043
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:50
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:71
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_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:80
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
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::online_send_diagqhyd
logical, public online_send_diagqhyd
Definition: scale_comm_cartesC_nest.F90:88
mod_atmos_bnd_driver::atmos_boundary_set_file
subroutine atmos_boundary_set_file(time)
Set boundary value for real case experiment.
Definition: mod_atmos_bnd_driver.F90:1383
scale_comm_cartesc_nest::online_boundary_use_qhyd
logical, public online_boundary_use_qhyd
Definition: scale_comm_cartesC_nest.F90:85
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:47
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:45
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:84
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:74
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:77
scale_file_external_input::i_prev
integer, parameter, public i_prev
[index] previous
Definition: scale_file_external_input.F90:174
mod_atmos_bnd_driver::calc_mass
subroutine calc_mass
Definition: mod_atmos_bnd_driver.F90:2381
scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue_send
subroutine, public comm_cartesc_nest_recvwait_issue_send
Sub-command for data transfer from parent to daughter: nestdown (parent side)
Definition: scale_comm_cartesC_nest.F90:2513
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:67
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:1044
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:48
scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue_recv
subroutine, public comm_cartesc_nest_recvwait_issue_recv
Sub-command for data transfer from parent to daughter: nestdown (daughter side)
Definition: scale_comm_cartesC_nest.F90:2562
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:39
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:78
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:81
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_comm_cartesc_nest::online_recv_qa
integer, public online_recv_qa
number of tracer received from the parent domain
Definition: scale_comm_cartesC_nest.F90:89
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:79
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:97
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
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_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
mod_atmos_bnd_driver::atmos_boundary_velz
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_velz
Definition: mod_atmos_bnd_driver.F90:46
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_comm_cartesc_nest::comm_cartesc_nest_test_send
subroutine, public comm_cartesc_nest_test_send
[check communication status] Inter-communication (parent side)
Definition: scale_comm_cartesC_nest.F90:3075
scale_file::file_aggregate
logical, public file_aggregate
Definition: scale_file.F90:196
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:796
mod_atmos_bnd_driver::history_bnd
subroutine history_bnd(ATMOS_BOUNDARY_DENS, ATMOS_BOUNDARY_VELZ, ATMOS_BOUNDARY_VELX, ATMOS_BOUNDARY_VELY, ATMOS_BOUNDARY_POTT, ATMOS_BOUNDARY_QTRC)
Definition: mod_atmos_bnd_driver.F90:2350
scale_comm_cartesc_nest::online_use_velz
logical, public online_use_velz
Definition: scale_comm_cartesC_nest.F90:83
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:1018
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:79
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
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:81
mod_atmos_bnd_driver::real
logical, public real
Definition: mod_atmos_bnd_driver.F90:64
scale_file_cartesc::file_cartesc_open
subroutine, public file_cartesc_open(basename, fid, single, aggregate)
open a netCDF file for read
Definition: scale_file_cartesC.F90:760
mod_atmos_bnd_driver::set_boundary
subroutine set_boundary
Definition: mod_atmos_bnd_driver.F90:2003
mod_atmos_phy_mp_vars::qe_mp
integer, public qe_mp
Definition: mod_atmos_phy_mp_vars.F90:80
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:72
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:63
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:74
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:88
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_comm_cartesc_nest::online_parent_dtsec
real(dp), public online_parent_dtsec
parent DT [sec]
Definition: scale_comm_cartesC_nest.F90:92
scale_file_external_input
module file / external_input
Definition: scale_file_external_input.F90:12
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:43
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:68
mod_atmos_bnd_driver::atmos_boundary_driver_set
subroutine, public atmos_boundary_driver_set(time)
set
Definition: mod_atmos_bnd_driver.F90:548
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:48
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:1553
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
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:95
scale_comm_cartesc_nest::comm_cartesc_nest_recv_cancel_recv
subroutine, public comm_cartesc_nest_recv_cancel_recv
Sub-command for data transfer from parent to daughter: nestdown (daughter side)
Definition: scale_comm_cartesC_nest.F90:2670
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:82