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  !
36 
37  !-----------------------------------------------------------------------------
38  !
39  !++ Public parameters & variables
40  !
41  integer, public :: bnd_qa
42 
43  real(RP), public, allocatable :: atmos_boundary_dens(:,:,:)
44  real(RP), public, allocatable :: atmos_boundary_velz(:,:,:)
45  real(RP), public, allocatable :: atmos_boundary_velx(:,:,:)
46  real(RP), public, allocatable :: atmos_boundary_vely(:,:,:)
47  real(RP), public, allocatable :: atmos_boundary_pott(:,:,:)
48  real(RP), public, allocatable :: atmos_boundary_qtrc(:,:,:,:)
49 
50  real(RP), public, allocatable :: atmos_boundary_alpha_dens(:,:,:)
51  real(RP), public, allocatable :: atmos_boundary_alpha_velz(:,:,:)
52  real(RP), public, allocatable :: atmos_boundary_alpha_velx(:,:,:)
53  real(RP), public, allocatable :: atmos_boundary_alpha_vely(:,:,:)
54  real(RP), public, allocatable :: atmos_boundary_alpha_pott(:,:,:)
55  real(RP), public, allocatable :: atmos_boundary_alpha_qtrc(:,:,:,:)
56 
57 
58  real(RP), public :: atmos_boundary_smoother_fact = 0.2_rp ! fact for smoother to damping
59 
60  logical, public :: atmos_boundary_update_flag = .false.
61 
62  !-----------------------------------------------------------------------------
63  !
64  !++ Private procedure
65  !
66  private :: atmos_boundary_var_fillhalo
67  private :: atmos_boundary_alpha_fillhalo
68  private :: atmos_boundary_ref_fillhalo
69  private :: atmos_boundary_setalpha
70  private :: atmos_boundary_setinitval
71  private :: atmos_boundary_read
72  private :: atmos_boundary_write
73  private :: atmos_boundary_generate
74  private :: atmos_boundary_initialize_file
75  private :: atmos_boundary_initialize_online
76  private :: atmos_boundary_update_file
77  private :: atmos_boundary_update_online_parent
78  private :: atmos_boundary_update_online_daughter
79  private :: atmos_boundary_firstsend
80  private :: atmos_boundary_send
81  private :: atmos_boundary_recv
82 
83  abstract interface
84  subroutine getbnd( &
85  bnd_DENS, &
86  bnd_VELZ, &
87  bnd_VELX, &
88  bnd_VELY, &
89  bnd_POTT, &
90  bnd_QTRC, &
91  now_step, &
92  update_step )
93  use scale_precision
94  implicit none
95 
96  real(RP), intent(out) :: bnd_dens(:,:,:)
97  real(RP), intent(out) :: bnd_velz(:,:,:)
98  real(RP), intent(out) :: bnd_velx(:,:,:)
99  real(RP), intent(out) :: bnd_vely(:,:,:)
100  real(RP), intent(out) :: bnd_pott(:,:,:)
101  real(RP), intent(out) :: bnd_qtrc(:,:,:,:)
102  integer, intent(in) :: now_step
103  integer, intent(in) :: update_step
104  end subroutine getbnd
105  end interface
106 
107  procedure(getbnd), pointer :: get_boundary => null()
108  private :: get_boundary
109  private :: get_boundary_same_parent
110  private :: get_boundary_nearest_neighbor
111  private :: get_boundary_lerp_initpoint
112  private :: get_boundary_lerp_midpoint
113 
114  !
115  !-----------------------------------------------------------------------------
116  !
117  !++ Private parameters & variables
118  !
119  character(len=H_SHORT), private :: atmos_boundary_type = 'NONE'
120  character(len=H_LONG), private :: atmos_boundary_in_basename = ''
121  logical, private :: atmos_boundary_in_check_coordinates = .true.
122  character(len=H_LONG), private :: atmos_boundary_out_basename = ''
123  character(len=H_MID), private :: atmos_boundary_out_title = 'SCALE-RM BOUNDARY CONDITION'
124  character(len=H_SHORT), private :: atmos_boundary_out_dtype = 'DEFAULT'
125 
126  logical, private :: atmos_boundary_use_dens = .false. ! read from file?
127  logical, private :: atmos_boundary_use_velz = .false. ! read from file?
128  logical, private :: atmos_boundary_use_velx = .false. ! read from file?
129  logical, private :: atmos_boundary_use_vely = .false. ! read from file?
130  logical, private :: atmos_boundary_use_pott = .false. ! read from file?
131  logical, private :: atmos_boundary_use_qv = .false. ! read from file?
132  logical, private :: atmos_boundary_use_qhyd = .false. ! read from file?
133  logical, private :: atmos_boundary_use_chem = .false. ! read from file?
134 
135  real(RP), private :: atmos_boundary_value_velz = 0.0_rp ! velocity w at boundary, 0 [m/s]
136  real(RP), private :: atmos_boundary_value_velx = 0.0_rp ! velocity u at boundary, 0 [m/s]
137  real(RP), private :: atmos_boundary_value_vely = 0.0_rp ! velocity v at boundary, 0 [m/s]
138  real(RP), private :: atmos_boundary_value_pott = 300.0_rp ! potential temp. at boundary, 300 [K]
139  real(RP), private :: atmos_boundary_value_qtrc = 0.0_rp ! tracer at boundary, 0 [kg/kg]
140 
141  real(RP), private :: atmos_boundary_alphafact_dens = 1.0_rp ! alpha factor again default
142  real(RP), private :: atmos_boundary_alphafact_velz = 1.0_rp ! alpha factor again default
143  real(RP), private :: atmos_boundary_alphafact_velx = 1.0_rp ! alpha factor again default
144  real(RP), private :: atmos_boundary_alphafact_vely = 1.0_rp ! alpha factor again default
145  real(RP), private :: atmos_boundary_alphafact_pott = 1.0_rp ! alpha factor again default
146  real(RP), private :: atmos_boundary_alphafact_qtrc = 1.0_rp ! alpha factor again default
147 
148  real(RP), private :: atmos_boundary_fracz = 1.0_rp ! fraction of boundary region for dumping (z) (0-1)
149  real(RP), private :: atmos_boundary_fracx = 1.0_rp ! fraction of boundary region for dumping (x) (0-1)
150  real(RP), private :: atmos_boundary_fracy = 1.0_rp ! fraction of boundary region for dumping (y) (0-1)
151  real(RP), private :: atmos_boundary_tauz ! maximum value for damping tau (z) [s]
152  real(RP), private :: atmos_boundary_taux ! maximum value for damping tau (x) [s]
153  real(RP), private :: atmos_boundary_tauy ! maximum value for damping tau (y) [s]
154 
155  real(DP), private :: atmos_boundary_update_dt = 0.0_dp ! inteval time of boudary data update [s]
156  integer, private :: update_nstep
157 
158  real(RP), private, allocatable :: atmos_boundary_ref_dens(:,:,:,:) ! reference DENS (with HALO)
159  real(RP), private, allocatable :: atmos_boundary_ref_velz(:,:,:,:) ! reference VELZ (with HALO)
160  real(RP), private, allocatable :: atmos_boundary_ref_velx(:,:,:,:) ! reference VELX (with HALO)
161  real(RP), private, allocatable :: atmos_boundary_ref_vely(:,:,:,:) ! reference VELY (with HALO)
162  real(RP), private, allocatable :: atmos_boundary_ref_pott(:,:,:,:) ! reference POTT (with HALO)
163  real(RP), private, allocatable, target :: atmos_boundary_ref_qtrc(:,:,:,:,:) ! reference QTRC (with HALO)
164 
165  ! work
166  real(RP), private, allocatable, target :: q_work(:,:,:,:) ! QV + Qe
167 
168 
169  character(len=H_LONG), private :: atmos_boundary_interp_type = 'lerp_initpoint' ! type of boundary interporation
170 
171  integer, private :: atmos_boundary_start_date(6) = (/ -9999, 0, 0, 0, 0, 0 /) ! boundary initial date
172 
173  integer, private :: atmos_boundary_fid = -1
174 
175  integer, private :: now_step
176  integer, private :: boundary_timestep = 0
177  logical, private :: atmos_boundary_linear_v = .false. ! linear or non-linear profile of relax region
178  logical, private :: atmos_boundary_linear_h = .true. ! linear or non-linear profile of relax region
179  real(RP), private :: atmos_boundary_exp_h = 2.0_rp ! factor of non-linear profile of relax region
180  logical, private :: atmos_boundary_online = .false. ! boundary online update by communicate inter-domain
181  logical, private :: atmos_boundary_online_master = .false. ! master domain in communicate inter-domain
182 
183  logical, private :: do_parent_process = .false.
184  logical, private :: do_daughter_process = .false.
185  logical, private :: l_bnd = .false.
186 
187  real(DP), private :: boundary_time_initdaysec
188 
189  integer, private :: ref_size = 3
190  integer, private :: ref_old = 1
191  integer, private :: ref_now = 2
192  integer, private :: ref_new = 3
193 
194  !-----------------------------------------------------------------------------
195 contains
196  !-----------------------------------------------------------------------------
198  subroutine atmos_boundary_driver_setup
199  use scale_prc, only: &
200  prc_abort
201  use scale_const, only: &
203  use scale_time, only: &
204  dt => time_dtsec
205  use scale_comm_cartesc_nest, only: &
208  use_nesting, &
209  offline, &
212  nestqa => comm_cartesc_nest_bnd_qa
213  use mod_atmos_phy_mp_vars, only: &
214  qa_mp
215  use mod_atmos_phy_ch_vars, only: &
216  qa_ch
217  implicit none
218 
219  namelist / param_atmos_boundary / &
220  atmos_boundary_type, &
221  atmos_boundary_in_basename, &
222  atmos_boundary_in_check_coordinates, &
223  atmos_boundary_out_basename, &
224  atmos_boundary_out_title, &
225  atmos_boundary_out_dtype, &
226  atmos_boundary_use_velz, &
227  atmos_boundary_use_velx, &
228  atmos_boundary_use_vely, &
229  atmos_boundary_use_pott, &
230  atmos_boundary_use_dens, &
231  atmos_boundary_use_qv, &
232  atmos_boundary_use_qhyd, &
233  atmos_boundary_use_chem, &
234  atmos_boundary_value_velz, &
235  atmos_boundary_value_velx, &
236  atmos_boundary_value_vely, &
237  atmos_boundary_value_pott, &
238  atmos_boundary_value_qtrc, &
239  atmos_boundary_alphafact_dens, &
240  atmos_boundary_alphafact_velz, &
241  atmos_boundary_alphafact_velx, &
242  atmos_boundary_alphafact_vely, &
243  atmos_boundary_alphafact_pott, &
244  atmos_boundary_alphafact_qtrc, &
246  atmos_boundary_fracz, &
247  atmos_boundary_fracx, &
248  atmos_boundary_fracy, &
249  atmos_boundary_tauz, &
250  atmos_boundary_taux, &
251  atmos_boundary_tauy, &
252  atmos_boundary_update_dt, &
253  atmos_boundary_start_date, &
254  atmos_boundary_linear_v, &
255  atmos_boundary_linear_h, &
256  atmos_boundary_exp_h, &
257  atmos_boundary_interp_type
258 
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  !--- read namelist
271  rewind(io_fid_conf)
272  read(io_fid_conf,nml=param_atmos_boundary,iostat=ierr)
273  if( ierr < 0 ) then !--- missing
274  log_info("ATMOS_BOUNDARY_setup",*) 'Not found namelist. Default used.'
275  elseif( ierr > 0 ) then !--- fatal error
276  log_error("ATMOS_BOUNDARY_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_BOUNDARY. Check!'
277  call prc_abort
278  endif
279  log_nml(param_atmos_boundary)
280 
281  ! setting switches
282  if( .NOT. use_nesting ) then
283  atmos_boundary_online = .false.
284  else
285  if( offline ) then
286  atmos_boundary_online = .false.
287  else
288  atmos_boundary_online = .true.
289  endif
290  endif
291  do_parent_process = .false.
292  do_daughter_process = .false.
293  atmos_boundary_online_master = .false.
294  if ( atmos_boundary_online ) then
295  if ( online_iam_parent ) then
296  do_parent_process = .true.
297  if ( .NOT. online_iam_daughter ) then
298  atmos_boundary_online_master = .true.
299  endif
300  endif
301  if ( online_iam_daughter ) then
302  do_daughter_process = .true.
303  atmos_boundary_use_qhyd = online_boundary_use_qhyd
304  endif
305  endif
306 
307  if( atmos_boundary_use_qhyd ) then
308  bnd_qa = qa_mp
309  if( atmos_boundary_use_chem ) then
310  bnd_qa = bnd_qa + qa_ch
311  endif
312  else if ( qa_mp > 0 ) then
313  bnd_qa = 1
314  else
315  bnd_qa = 0
316  end if
317 
318  allocate( atmos_boundary_dens(ka,ia,ja) )
319  allocate( atmos_boundary_velz(ka,ia,ja) )
320  allocate( atmos_boundary_velx(ka,ia,ja) )
321  allocate( atmos_boundary_vely(ka,ia,ja) )
322  allocate( atmos_boundary_pott(ka,ia,ja) )
323  allocate( atmos_boundary_qtrc(ka,ia,ja,bnd_qa) )
330 
331  allocate( atmos_boundary_alpha_dens(ka,ia,ja) )
332  allocate( atmos_boundary_alpha_velz(ka,ia,ja) )
333  allocate( atmos_boundary_alpha_velx(ka,ia,ja) )
334  allocate( atmos_boundary_alpha_vely(ka,ia,ja) )
335  allocate( atmos_boundary_alpha_pott(ka,ia,ja) )
337  atmos_boundary_alpha_dens(:,:,:) = 0.0_rp
338  atmos_boundary_alpha_velz(:,:,:) = 0.0_rp
339  atmos_boundary_alpha_velx(:,:,:) = 0.0_rp
340  atmos_boundary_alpha_vely(:,:,:) = 0.0_rp
341  atmos_boundary_alpha_pott(:,:,:) = 0.0_rp
342  atmos_boundary_alpha_qtrc(:,:,:,:) = 0.0_rp
343 
344  if ( atmos_boundary_type == 'REAL' .OR. do_daughter_process ) then
345  l_bnd = .true.
346  else
347  l_bnd = .false.
348  end if
349 
350  if ( l_bnd ) then
351 
352  select case(atmos_boundary_interp_type)
353  case('same_parent')
354  get_boundary => get_boundary_same_parent
355  case('nearest_neighbor')
356  get_boundary => get_boundary_nearest_neighbor
357  case('lerp_initpoint')
358  get_boundary => get_boundary_lerp_initpoint
359  case('lerp_midpoint')
360  get_boundary => get_boundary_lerp_midpoint
361  case default
362  log_error("ATMOS_BOUNDARY_setup",*) 'Wrong parameter in ATMOS_BOUNDARY_interp_TYPE. Check!'
363  call prc_abort
364  end select
365 
366  allocate( atmos_boundary_ref_dens(ka,ia,ja,ref_size) )
367  allocate( atmos_boundary_ref_velz(ka,ia,ja,ref_size) )
368  allocate( atmos_boundary_ref_velx(ka,ia,ja,ref_size) )
369  allocate( atmos_boundary_ref_vely(ka,ia,ja,ref_size) )
370  allocate( atmos_boundary_ref_pott(ka,ia,ja,ref_size) )
371  allocate( atmos_boundary_ref_qtrc(ka,ia,ja,bnd_qa,ref_size) )
372  atmos_boundary_ref_dens(:,:,:,:) = const_undef
373  atmos_boundary_ref_velz(:,:,:,:) = const_undef
374  atmos_boundary_ref_velx(:,:,:,:) = const_undef
375  atmos_boundary_ref_vely(:,:,:,:) = const_undef
376  atmos_boundary_ref_pott(:,:,:,:) = const_undef
377  atmos_boundary_ref_qtrc(:,:,:,:,:) = const_undef
378 
379  ! initialize boundary value (reading file or waiting parent domain)
380  if ( do_daughter_process ) then
381  call atmos_boundary_initialize_online
382  else
383  if ( atmos_boundary_in_basename /= '' ) then
384  call atmos_boundary_initialize_file
385  else
386  log_error("ATMOS_BOUNDARY_setup",*) 'You need specify ATMOS_BOUNDARY_IN_BASENAME'
387  call prc_abort
388  endif
389  endif
390 
391  call atmos_boundary_setalpha
392 
394 
395  elseif ( atmos_boundary_type == 'NONE' ) then
396 
398 
399  elseif ( atmos_boundary_type == 'CONST' ) then
400 
401  call atmos_boundary_setalpha
402 
404 
405  elseif ( atmos_boundary_type == 'INIT' ) then
406 
407  call atmos_boundary_setalpha
408 
410 
411  elseif ( atmos_boundary_type == 'OFFLINE' ) then
412 
413  if ( atmos_boundary_in_basename /= '' ) then
414  call atmos_boundary_read
415  else
416  log_error("ATMOS_BOUNDARY_setup",*) 'You need specify ATMOS_BOUNDARY_IN_BASENAME'
417  call prc_abort
418  endif
419 
421 
422  else
423  log_error("ATMOS_BOUNDARY_setup",*) 'unsupported ATMOS_BOUNDARY_TYPE. Check!', trim(atmos_boundary_type)
424  call prc_abort
425  endif
426 
428 
429  !----- report data -----
430  log_newline
431  log_info("ATMOS_BOUNDARY_setup",*) 'Atmospheric boundary parameters '
432  log_info_cont(*) 'Atmospheric boundary type : ', atmos_boundary_type
433  log_newline
434  log_info_cont(*) 'Is VELZ used in atmospheric boundary? : ', atmos_boundary_use_velz
435  log_info_cont(*) 'Is VELX used in atmospheric boundary? : ', atmos_boundary_use_velx
436  log_info_cont(*) 'Is VELY used in atmospheric boundary? : ', atmos_boundary_use_vely
437  log_info_cont(*) 'Is POTT used in atmospheric boundary? : ', atmos_boundary_use_pott
438  log_info_cont(*) 'Is DENS used in atmospheric boundary? : ', atmos_boundary_use_dens
439  log_info_cont(*) 'Is QV used in atmospheric boundary? : ', atmos_boundary_use_qv
440  log_info_cont(*) 'Is QHYD used in atmospheric boundary? : ', atmos_boundary_use_qhyd
441  log_info_cont(*) 'Is CHEM used in atmospheric boundary? : ', atmos_boundary_use_chem
442  log_newline
443  log_info_cont(*) 'Atmospheric boundary VELZ values : ', atmos_boundary_value_velz
444  log_info_cont(*) 'Atmospheric boundary VELX values : ', atmos_boundary_value_velx
445  log_info_cont(*) 'Atmospheric boundary VELY values : ', atmos_boundary_value_vely
446  log_info_cont(*) 'Atmospheric boundary POTT values : ', atmos_boundary_value_pott
447  log_info_cont(*) 'Atmospheric boundary QTRC values : ', atmos_boundary_value_qtrc
448  log_newline
449  log_info_cont(*) 'Atmospheric boundary smoother factor : ', atmos_boundary_smoother_fact
450  log_info_cont(*) 'Atmospheric boundary z-fraction : ', atmos_boundary_fracz
451  log_info_cont(*) 'Atmospheric boundary x-fraction : ', atmos_boundary_fracx
452  log_info_cont(*) 'Atmospheric boundary y-fraction : ', atmos_boundary_fracy
453  log_info_cont(*) 'Atmospheric boundary z-relaxation time : ', atmos_boundary_tauz
454  log_info_cont(*) 'Atmospheric boundary x-relaxation time : ', atmos_boundary_taux
455  log_info_cont(*) 'Atmospheric boundary y-relaxation time : ', atmos_boundary_tauy
456  log_newline
457  log_info_cont(*) 'Atmospheric boundary update dt : ', atmos_boundary_update_dt
458  log_info_cont(*) 'Atmospheric boundary start date : ', atmos_boundary_start_date(:)
459  log_newline
460  log_info_cont(*) 'Linear profile in vertically relax region : ', atmos_boundary_linear_v
461  log_info_cont(*) 'Linear profile in horizontally relax region : ', atmos_boundary_linear_h
462  log_info_cont(*) 'Non-linear factor in horizontally relax region : ', atmos_boundary_exp_h
463  log_newline
464  log_info_cont(*) 'Online nesting for lateral boundary : ', atmos_boundary_online
465 
466  log_info_cont(*) 'Does lateral boundary exist in this domain? : ', l_bnd
467  if ( l_bnd ) then
468  log_info_cont(*) 'Lateral boundary interporation type : ', atmos_boundary_interp_type
469  endif
470 
471  if ( online_boundary_diagqhyd ) then
472  allocate( q_work(ka,ia,ja,nestqa) )
473  end if
474 
475  return
476  end subroutine atmos_boundary_driver_setup
477 
478  !-----------------------------------------------------------------------------
480  subroutine atmos_boundary_driver_set
481  use mod_atmos_vars, only: &
482  dens, &
483  momz, &
484  momx, &
485  momy, &
486  rhot, &
487  qtrc, &
488  qv, &
489  qe
490  use mod_atmos_phy_mp_vars, only: &
491  qs_mp, &
492  qe_mp
493  implicit none
494 
495  if ( do_parent_process ) then !online [parent]
496  call atmos_boundary_firstsend( &
497  dens, momz, momx, momy, rhot, qtrc(:,:,:,qs_mp:qe_mp), qv, qe )
498  end if
499 
500  if ( l_bnd ) then
501 
502  ! initialize boundary value (reading file or waiting parent domain)
503  if ( do_daughter_process ) then
505  else
506  if ( atmos_boundary_in_basename /= '' ) then
508  endif
509  endif
510 
511  elseif ( atmos_boundary_type == 'CONST' ) then
512 
513  call atmos_boundary_generate
514 
515  elseif ( atmos_boundary_type == 'INIT' ) then
516 
517  call atmos_boundary_setinitval( dens, & ! [IN]
518  momz, & ! [IN]
519  momx, & ! [IN]
520  momy, & ! [IN]
521  rhot, & ! [IN]
522  qtrc ) ! [IN]
523  endif
524 
525  if( atmos_boundary_out_basename /= '' ) then
526  call atmos_boundary_write
527  endif
528 
529  if ( atmos_boundary_update_flag ) then
530 
531  call history_bnd( &
538  end if
539 
540  return
541  end subroutine atmos_boundary_driver_set
542 
543  !-----------------------------------------------------------------------------
545  subroutine atmos_boundary_var_fillhalo
546  use scale_comm_cartesc, only: &
547  comm_vars8, &
548  comm_wait
549  implicit none
550 
551  integer :: i, j, iq
552  !---------------------------------------------------------------------------
553 
554  do j = 1, ja
555  do i = 1, ia
561 
567 
568  do iq = 1, bnd_qa
569  atmos_boundary_qtrc( 1:ks-1,i,j,iq) = atmos_boundary_qtrc(ks,i,j,iq)
570  atmos_boundary_qtrc(ke+1:ka, i,j,iq) = atmos_boundary_qtrc(ke,i,j,iq)
571  end do
572  end do
573  end do
574 
575  call comm_vars8( atmos_boundary_dens(:,:,:), 1 )
576  call comm_vars8( atmos_boundary_velz(:,:,:), 2 )
577  call comm_vars8( atmos_boundary_velx(:,:,:), 3 )
578  call comm_vars8( atmos_boundary_vely(:,:,:), 4 )
579  call comm_vars8( atmos_boundary_pott(:,:,:), 5 )
580  do iq = 1, bnd_qa
581  call comm_vars8( atmos_boundary_qtrc(:,:,:,iq), 5+iq )
582  end do
583 
584  call comm_wait ( atmos_boundary_dens(:,:,:), 1, .false. )
585  call comm_wait ( atmos_boundary_velz(:,:,:), 2, .false. )
586  call comm_wait ( atmos_boundary_velx(:,:,:), 3, .false. )
587  call comm_wait ( atmos_boundary_vely(:,:,:), 4, .false. )
588  call comm_wait ( atmos_boundary_pott(:,:,:), 5, .false. )
589  do iq = 1, bnd_qa
590  call comm_wait ( atmos_boundary_qtrc(:,:,:,iq), 5+iq, .false. )
591  end do
592 
593  return
594  end subroutine atmos_boundary_var_fillhalo
595 
596  !-----------------------------------------------------------------------------
598  subroutine atmos_boundary_alpha_fillhalo
599  use scale_comm_cartesc, only: &
600  comm_vars8, &
601  comm_wait
602  implicit none
603 
604  integer :: i, j, iq
605  !---------------------------------------------------------------------------
606 
607  do j = 1, ja
608  do i = 1, ia
614 
620 
621  do iq = 1, bnd_qa
624  end do
625  enddo
626  enddo
627 
628  call comm_vars8( atmos_boundary_alpha_dens(:,:,:), 1 )
629  call comm_vars8( atmos_boundary_alpha_velz(:,:,:), 2 )
630  call comm_vars8( atmos_boundary_alpha_velx(:,:,:), 3 )
631  call comm_vars8( atmos_boundary_alpha_vely(:,:,:), 4 )
632  call comm_vars8( atmos_boundary_alpha_pott(:,:,:), 5 )
633  do iq = 1, bnd_qa
634  call comm_vars8( atmos_boundary_alpha_qtrc(:,:,:,iq), 5+iq )
635  end do
636 
637  call comm_wait ( atmos_boundary_alpha_dens(:,:,:), 1, .false. )
638  call comm_wait ( atmos_boundary_alpha_velz(:,:,:), 2, .false. )
639  call comm_wait ( atmos_boundary_alpha_velx(:,:,:), 3, .false. )
640  call comm_wait ( atmos_boundary_alpha_vely(:,:,:), 4, .false. )
641  call comm_wait ( atmos_boundary_alpha_pott(:,:,:), 5, .false. )
642  do iq = 1, bnd_qa
643  call comm_wait ( atmos_boundary_alpha_qtrc(:,:,:,iq), 5+iq, .false. )
644  end do
645 
646  return
647  end subroutine atmos_boundary_alpha_fillhalo
648 
649  !-----------------------------------------------------------------------------
651  subroutine atmos_boundary_ref_fillhalo( &
652  ref_idx )
653  use scale_comm_cartesc, only: &
654  comm_vars8, &
655  comm_wait
656  implicit none
657 
658  ! arguments
659  integer, intent(in) :: ref_idx
660 
661  ! works
662  integer :: i, j, iq
663  !---------------------------------------------------------------------------
664 
665  !$omp parallel do default(none) private(i,j,iq) OMP_SCHEDULE_ collapse(2) &
666  !$omp shared(JSB,JEB,ISB,IEB,ATMOS_BOUNDARY_ref_DENS,ref_idx,ATMOS_BOUNDARY_ref_VELZ) &
667  !$omp shared(ATMOS_BOUNDARY_ref_VELX,ATMOS_BOUNDARY_ref_VELY,ATMOS_BOUNDARY_ref_POTT) &
668  !$omp shared(KA,KS,BND_QA,ATMOS_BOUNDARY_ref_QTRC,KE)
669  do j = jsb, jeb
670  do i = isb, ieb
671  atmos_boundary_ref_dens( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_dens(ks,i,j,ref_idx)
672  atmos_boundary_ref_velz( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_velz(ks,i,j,ref_idx)
673  atmos_boundary_ref_velx( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_velx(ks,i,j,ref_idx)
674  atmos_boundary_ref_vely( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_vely(ks,i,j,ref_idx)
675  atmos_boundary_ref_pott( 1:ks-1,i,j,ref_idx) = atmos_boundary_ref_pott(ks,i,j,ref_idx)
676 
677  atmos_boundary_ref_dens(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_dens(ke,i,j,ref_idx)
678  atmos_boundary_ref_velz(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_velz(ke,i,j,ref_idx)
679  atmos_boundary_ref_velx(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_velx(ke,i,j,ref_idx)
680  atmos_boundary_ref_vely(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_vely(ke,i,j,ref_idx)
681  atmos_boundary_ref_pott(ke+1:ka, i,j,ref_idx) = atmos_boundary_ref_pott(ke,i,j,ref_idx)
682 
683  do iq = 1, bnd_qa
684  atmos_boundary_ref_qtrc( 1:ks-1,i,j,iq,ref_idx) = atmos_boundary_ref_qtrc(ks,i,j,iq,ref_idx)
685  atmos_boundary_ref_qtrc(ke+1:ka, i,j,iq,ref_idx) = atmos_boundary_ref_qtrc(ke,i,j,iq,ref_idx)
686  end do
687  end do
688  end do
689 
690  call comm_vars8( atmos_boundary_ref_dens(:,:,:,ref_idx), 1 )
691  call comm_vars8( atmos_boundary_ref_velz(:,:,:,ref_idx), 2 )
692  call comm_vars8( atmos_boundary_ref_velx(:,:,:,ref_idx), 3 )
693  call comm_vars8( atmos_boundary_ref_vely(:,:,:,ref_idx), 4 )
694  call comm_vars8( atmos_boundary_ref_pott(:,:,:,ref_idx), 5 )
695 
696  do iq = 1, bnd_qa
697  call comm_vars8( atmos_boundary_ref_qtrc(:,:,:,iq,ref_idx), 5+iq )
698  end do
699 
700  call comm_wait ( atmos_boundary_ref_dens(:,:,:,ref_idx), 1, .false. )
701  call comm_wait ( atmos_boundary_ref_velz(:,:,:,ref_idx), 2, .false. )
702  call comm_wait ( atmos_boundary_ref_velx(:,:,:,ref_idx), 3, .false. )
703  call comm_wait ( atmos_boundary_ref_vely(:,:,:,ref_idx), 4, .false. )
704  call comm_wait ( atmos_boundary_ref_pott(:,:,:,ref_idx), 5, .false. )
705 
706  do iq = 1, bnd_qa
707  call comm_wait ( atmos_boundary_ref_qtrc(:,:,:,iq,ref_idx), 5+iq, .false. )
708  end do
709 
710  return
711  end subroutine atmos_boundary_ref_fillhalo
712 
713  !-----------------------------------------------------------------------------
715  subroutine atmos_boundary_setalpha
716  use scale_const, only: &
717  eps => const_eps, &
718  pi => const_pi
719  use scale_atmos_grid_cartesc, only: &
720  cbfz => atmos_grid_cartesc_cbfz, &
721  cbfx => atmos_grid_cartesc_cbfx, &
722  cbfy => atmos_grid_cartesc_cbfy, &
723  fbfz => atmos_grid_cartesc_fbfz, &
724  fbfx => atmos_grid_cartesc_fbfx, &
726  use scale_comm_cartesc_nest, only: &
728 
729  real(RP) :: coef_z, alpha_z1, alpha_z2
730  real(RP) :: coef_x, alpha_x1, alpha_x2
731  real(RP) :: coef_y, alpha_y1, alpha_y2
732  real(RP) :: ee1, ee2
733 
734  integer :: i, j, k, iq
735  !---------------------------------------------------------------------------
736 
737  ! check invalid fraction
738  atmos_boundary_fracz = max( min( atmos_boundary_fracz, 1.0_rp ), eps )
739  atmos_boundary_fracx = max( min( atmos_boundary_fracx, 1.0_rp ), eps )
740  atmos_boundary_fracy = max( min( atmos_boundary_fracy, 1.0_rp ), eps )
741 
742  if ( atmos_boundary_tauz <= 0.0_rp ) then ! invalid tau
743  coef_z = 0.0_rp
744  else
745  coef_z = 1.0_rp / atmos_boundary_tauz
746  endif
747 
748  if ( atmos_boundary_taux <= 0.0_rp ) then ! invalid tau
749  coef_x = 0.0_rp
750  else
751  coef_x = 1.0_rp / atmos_boundary_taux
752  endif
753 
754  if ( atmos_boundary_tauy <= 0.0_rp ) then ! invalid tau
755  coef_y = 0.0_rp
756  else
757  coef_y = 1.0_rp / atmos_boundary_tauy
758  endif
759 
760  !$omp parallel do default(none) &
761  !$omp shared(JA,IA,KA,CBFZ,ATMOS_BOUNDARY_FRACZ,FBFZ,ATMOS_BOUNDARY_LINEAR_V,coef_z,CBFX) &
762  !$omp shared(ATMOS_BOUNDARY_FRACX,PI,FBFX,ATMOS_BOUNDARY_LINEAR_H,coef_x) &
763  !$omp shared(ATMOS_BOUNDARY_EXP_H,CBFY,ATMOS_BOUNDARY_FRACY,FBFY,coef_y,l_bnd) &
764  !$omp shared(do_daughter_process,ONLINE_USE_VELZ,ATMOS_BOUNDARY_USE_VELZ,ATMOS_BOUNDARY_alpha_VELZ) &
765  !$omp shared(ATMOS_BOUNDARY_ALPHAFACT_VELZ,ATMOS_BOUNDARY_USE_DENS,ATMOS_BOUNDARY_alpha_DENS) &
766  !$omp shared(ATMOS_BOUNDARY_ALPHAFACT_DENS,ATMOS_BOUNDARY_USE_VELX,ATMOS_BOUNDARY_alpha_VELX) &
767  !$omp shared(ATMOS_BOUNDARY_USE_VELY,ATMOS_BOUNDARY_alpha_VELY,ATMOS_BOUNDARY_ALPHAFACT_VELY) &
768  !$omp shared(ATMOS_BOUNDARY_USE_POTT,ATMOS_BOUNDARY_alpha_POTT,ATMOS_BOUNDARY_ALPHAFACT_POTT) &
769  !$omp shared(ATMOS_BOUNDARY_USE_QV,ATMOS_BOUNDARY_alpha_QTRC,ATMOS_BOUNDARY_ALPHAFACT_QTRC) &
770  !$omp shared(ATMOS_BOUNDARY_USE_QHYD,BND_QA,ATMOS_BOUNDARY_ALPHAFACT_VELX) &
771  !$omp private(i,j,k,ee1,ee2,alpha_z1,alpha_z2,alpha_x1,alpha_x2,alpha_y1,alpha_y2,iq) OMP_SCHEDULE_ collapse(2)
772  do j = 1, ja
773  do i = 1, ia
774  do k = 1, ka
775  ee1 = cbfz(k)
776  if ( ee1 <= 1.0_rp - atmos_boundary_fracz ) then
777  ee1 = 0.0_rp
778  else
779  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracz ) / atmos_boundary_fracz
780  endif
781 
782  ee2 = fbfz(k)
783  if ( ee2 <= 1.0_rp - atmos_boundary_fracz ) then
784  ee2 = 0.0_rp
785  else
786  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracz ) / atmos_boundary_fracz
787  endif
788 
789  alpha_z1 = 0.0_rp
790  alpha_z2 = 0.0_rp
791  if ( atmos_boundary_linear_v ) then
792  alpha_z1 = coef_z * ee1
793  alpha_z2 = coef_z * ee2
794  else
795  if ( ee1 > 0.0_rp .AND. ee1 <= 0.5_rp ) then
796  alpha_z1 = coef_z * 0.5_rp * ( 1.0_rp - cos( ee1*pi ) )
797  elseif( ee1 > 0.5_rp .AND. ee1 <= 1.0_rp ) then
798  alpha_z1 = coef_z * 0.5_rp * ( 1.0_rp + sin( (ee1-0.5_rp)*pi ) )
799  endif
800  if ( ee2 > 0.0_rp .AND. ee2 <= 0.5_rp ) then
801  alpha_z2 = coef_z * 0.5_rp * ( 1.0_rp - cos( ee2*pi ) )
802  elseif( ee2 > 0.5_rp .AND. ee2 <= 1.0_rp ) then
803  alpha_z2 = coef_z * 0.5_rp * ( 1.0_rp + sin( (ee2-0.5_rp)*pi ) )
804  endif
805  endif
806 
807  ee1 = cbfx(i)
808  if ( ee1 <= 1.0_rp - atmos_boundary_fracx ) then
809  ee1 = 0.0_rp
810  else
811  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracx ) / atmos_boundary_fracx
812  endif
813 
814  ee2 = fbfx(i)
815  if ( ee2 <= 1.0_rp - atmos_boundary_fracx ) then
816  ee2 = 0.0_rp
817  else
818  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracx ) / atmos_boundary_fracx
819  endif
820 
821  if ( atmos_boundary_linear_h ) then
822  alpha_x1 = coef_x * ee1
823  alpha_x2 = coef_x * ee2
824  else
825  alpha_x1 = coef_x * ee1 * exp( -(1.0_rp-ee1) * atmos_boundary_exp_h )
826  alpha_x2 = coef_x * ee2 * exp( -(1.0_rp-ee2) * atmos_boundary_exp_h )
827  end if
828 
829  ee1 = cbfy(j)
830  if ( ee1 <= 1.0_rp - atmos_boundary_fracy ) then
831  ee1 = 0.0_rp
832  else
833  ee1 = ( ee1 - 1.0_rp + atmos_boundary_fracy ) / atmos_boundary_fracy
834  endif
835 
836  ee2 = fbfy(j)
837  if ( ee2 <= 1.0_rp - atmos_boundary_fracy ) then
838  ee2 = 0.0_rp
839  else
840  ee2 = ( ee2 - 1.0_rp + atmos_boundary_fracy ) / atmos_boundary_fracy
841  endif
842 
843  if ( atmos_boundary_linear_h ) then
844  alpha_y1 = coef_y * ee1
845  alpha_y2 = coef_y * ee2
846  else
847  alpha_y1 = coef_y * ee1 * exp( -(1.0_rp-ee1) * atmos_boundary_exp_h )
848  alpha_y2 = coef_y * ee2 * exp( -(1.0_rp-ee2) * atmos_boundary_exp_h )
849  end if
850 
851 
852  if ( l_bnd ) then
853  if ( do_daughter_process ) then ! online
854  if ( online_use_velz ) then
855  if( atmos_boundary_use_velz ) then
856  atmos_boundary_alpha_velz(k,i,j) = max( alpha_z2, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
857  else
858  atmos_boundary_alpha_velz(k,i,j) = max( alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
859  endif
860  else
861  if ( atmos_boundary_use_velz ) then
862  atmos_boundary_alpha_velz(k,i,j) = alpha_z2 * atmos_boundary_alphafact_velz
863  end if
864  end if
865  else ! offline
866  if ( atmos_boundary_use_velz ) then
867  atmos_boundary_alpha_velz(k,i,j) = max( alpha_z2, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
868  endif
869  end if
870  if ( atmos_boundary_use_dens ) then
871  atmos_boundary_alpha_dens(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_dens
872  else
873  atmos_boundary_alpha_dens(k,i,j) = 0.0_rp
874 ! ATMOS_BOUNDARY_alpha_DENS(k,i,j) = max( alpha_x1, alpha_y1 ) * ATMOS_BOUNDARY_ALPHAFACT_DENS
875  endif
876  if ( atmos_boundary_use_velx ) then
877  atmos_boundary_alpha_velx(k,i,j) = max( alpha_z1, alpha_x2, alpha_y1 ) * atmos_boundary_alphafact_velx
878  else
879  atmos_boundary_alpha_velx(k,i,j) = max( alpha_x2, alpha_y1 ) * atmos_boundary_alphafact_velx
880  endif
881  if ( atmos_boundary_use_vely ) then
882  atmos_boundary_alpha_vely(k,i,j) = max( alpha_z1, alpha_x1, alpha_y2 ) * atmos_boundary_alphafact_vely
883  else
884  atmos_boundary_alpha_vely(k,i,j) = max( alpha_x1, alpha_y2 ) * atmos_boundary_alphafact_vely
885  endif
886  if ( atmos_boundary_use_pott ) then
887  atmos_boundary_alpha_pott(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_pott
888  else
889  atmos_boundary_alpha_pott(k,i,j) = max( alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_pott
890  endif
891  if ( atmos_boundary_use_qv ) then
892  atmos_boundary_alpha_qtrc(k,i,j,1) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
893  else
894  atmos_boundary_alpha_qtrc(k,i,j,1) = max( alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
895  endif
896  if ( atmos_boundary_use_qhyd ) then
897  do iq = 2, bnd_qa
898  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
899  end do
900  else
901  do iq = 2, bnd_qa
902  atmos_boundary_alpha_qtrc(k,i,j,iq) = max( alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_qtrc
903  end do
904  endif
905  else
906  atmos_boundary_alpha_dens(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_dens
907  atmos_boundary_alpha_velz(k,i,j) = max( alpha_z2, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_velz
908  atmos_boundary_alpha_velx(k,i,j) = max( alpha_z1, alpha_x2, alpha_y1 ) * atmos_boundary_alphafact_velx
909  atmos_boundary_alpha_vely(k,i,j) = max( alpha_z1, alpha_x1, alpha_y2 ) * atmos_boundary_alphafact_vely
910  atmos_boundary_alpha_pott(k,i,j) = max( alpha_z1, alpha_x1, alpha_y1 ) * atmos_boundary_alphafact_pott
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  end if
915  enddo
916  enddo
917  enddo
918 
919  if ( l_bnd ) then
920  if ( .NOT. online_use_velz .AND. .NOT. atmos_boundary_use_velz ) then
921  atmos_boundary_alpha_velz(:,:,:) = 0.0_rp
922  end if
923  else
924  if ( .NOT. atmos_boundary_use_dens ) then
925  atmos_boundary_alpha_dens(:,:,:) = 0.0_rp
926  end if
927  if ( .NOT. atmos_boundary_use_velz ) then
928  atmos_boundary_alpha_velz(:,:,:) = 0.0_rp
929  end if
930  if ( .NOT. atmos_boundary_use_velx ) then
931  atmos_boundary_alpha_velx(:,:,:) = 0.0_rp
932  end if
933  if ( .NOT. atmos_boundary_use_vely ) then
934  atmos_boundary_alpha_vely(:,:,:) = 0.0_rp
935  end if
936  if ( .NOT. atmos_boundary_use_pott ) then
937  atmos_boundary_alpha_pott(:,:,:) = 0.0_rp
938  end if
939  if ( .NOT. atmos_boundary_use_qv ) then
940  if ( bnd_qa > 0 ) atmos_boundary_alpha_qtrc(:,:,:,1) = 0.0_rp
941  end if
942  if ( .NOT. atmos_boundary_use_qhyd ) then
943  do iq = 2, bnd_qa
944  atmos_boundary_alpha_qtrc(:,:,:,iq) = 0.0_rp
945  end do
946  end if
947  end if
948 
949 
950  call atmos_boundary_alpha_fillhalo
951 
952  return
953  end subroutine atmos_boundary_setalpha
954 
955  !-----------------------------------------------------------------------------
957  subroutine atmos_boundary_setinitval( &
958  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC )
959  implicit none
960 
961  real(RP), intent(in) :: DENS(ka,ia,ja)
962  real(RP), intent(in) :: MOMZ(ka,ia,ja)
963  real(RP), intent(in) :: MOMX(ka,ia,ja)
964  real(RP), intent(in) :: MOMY(ka,ia,ja)
965  real(RP), intent(in) :: RHOT(ka,ia,ja)
966  real(RP), intent(in) :: QTRC(ka,ia,ja,bnd_qa)
967 
968  integer :: i, j, k, iq
969  !---------------------------------------------------------------------------
970 
971  do j = 1, ja
972  do i = 1, ia
973  do k = ks, ke
974  atmos_boundary_dens(k,i,j) = dens(k,i,j)
975  atmos_boundary_velz(k,i,j) = momz(k,i,j) / ( dens(k,i,j)+dens(k+1,i, j ) ) * 2.0_rp
976  atmos_boundary_pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
977  do iq = 1, bnd_qa
978  atmos_boundary_qtrc(k,i,j,iq) = qtrc(k,i,j,iq)
979  end do
980  enddo
981  enddo
982  enddo
983 
984  do j = 1, ja
985  do i = 1, ia-1
986  do k = ks, ke
987  atmos_boundary_velx(k,i,j) = momx(k,i,j) / ( dens(k,i,j)+dens(k, i+1,j ) ) * 2.0_rp
988  enddo
989  enddo
990  enddo
991  do j = 1, ja
992  do k = ks, ke
993  atmos_boundary_velx(k,ia,j) = momx(k,ia,j) / dens(k,ia,j)
994  enddo
995  enddo
996 
997  do j = 1, ja-1
998  do i = 1, ia
999  do k = ks, ke
1000  atmos_boundary_vely(k,i,j) = momy(k,i,j) / ( dens(k,i,j)+dens(k, i, j+1) ) * 2.0_rp
1001  enddo
1002  enddo
1003  enddo
1004  do i = 1, ia
1005  do k = ks, ke
1006  atmos_boundary_vely(k,i,ja) = momy(k,i,ja) / dens(k,i,ja)
1007  enddo
1008  enddo
1009 
1010  call atmos_boundary_var_fillhalo
1011 
1012  return
1013  end subroutine atmos_boundary_setinitval
1014 
1015  !-----------------------------------------------------------------------------
1017  subroutine atmos_boundary_read
1018  use scale_prc, only: &
1019  prc_abort
1020  use scale_file_cartesc, only: &
1022  file_cartesc_check_coordinates, &
1023  file_cartesc_read, &
1025  implicit none
1026 
1027  integer :: fid
1028  integer :: iq
1029  !---------------------------------------------------------------------------
1030 
1031  call file_cartesc_open( atmos_boundary_in_basename, fid )
1032 
1033  if ( atmos_boundary_in_check_coordinates ) then
1034  call file_cartesc_check_coordinates( fid, atmos=.true. )
1035  end if
1036 
1037  if ( atmos_boundary_use_dens &
1038  .OR. atmos_boundary_use_velz &
1039  .OR. atmos_boundary_use_velx &
1040  .OR. atmos_boundary_use_vely &
1041  .OR. atmos_boundary_use_pott &
1042  ) then
1043  call file_cartesc_read( fid, 'DENS', 'ZXY', atmos_boundary_dens(:,:,:) )
1044  end if
1045  if ( atmos_boundary_use_dens ) then
1046  call file_cartesc_read( fid, 'ALPHA_DENS', 'ZXY', atmos_boundary_alpha_dens(:,:,:) )
1047  endif
1048 
1049  if ( atmos_boundary_use_velz ) then
1050  call file_cartesc_read( fid, 'VELZ', 'ZHXY', atmos_boundary_velz(:,:,:) )
1051  call file_cartesc_read( fid, 'ALPHA_VELZ', 'ZHXY', atmos_boundary_alpha_velz(:,:,:) )
1052  endif
1053 
1054  if ( atmos_boundary_use_velx ) then
1055  call file_cartesc_read( fid, 'VELX', 'ZXHY', atmos_boundary_velx(:,:,:) )
1056  call file_cartesc_read( fid, 'ALPHA_VELX', 'ZXHY', atmos_boundary_alpha_velx(:,:,:) )
1057  endif
1058 
1059  if ( atmos_boundary_use_vely ) then
1060  call file_cartesc_read( fid, 'VELY', 'ZXYH', atmos_boundary_vely(:,:,:) )
1061  call file_cartesc_read( fid, 'ALPHA_VELY', 'ZXYH', atmos_boundary_alpha_vely(:,:,:) )
1062  endif
1063 
1064  if ( atmos_boundary_use_pott ) then
1065  call file_cartesc_read( fid, 'POTT', 'ZXY', atmos_boundary_pott(:,:,:) )
1066  call file_cartesc_read( fid, 'ALPHA_POTT', 'ZXY', atmos_boundary_alpha_pott(:,:,:) )
1067  endif
1068 
1069  if ( atmos_boundary_use_qv ) then
1070  call file_cartesc_read( fid, 'QV', 'ZXY', atmos_boundary_qtrc(:,:,:,1) )
1071  call file_cartesc_read( fid, 'ALPHA_QV', 'ZXY', atmos_boundary_alpha_qtrc(:,:,:,1) )
1072  endif
1073 
1074  if ( atmos_boundary_use_qhyd ) then
1075  do iq = 2, bnd_qa
1076  call file_cartesc_read( fid, tracer_name(iq), 'ZXY', atmos_boundary_qtrc(:,:,:,iq) )
1077  call file_cartesc_read( fid, 'ALPHA_'//trim(tracer_name(iq)), 'ZXY', atmos_boundary_alpha_qtrc(:,:,:,iq) )
1078  end do
1079  endif
1080 
1081  call file_cartesc_close( fid )
1082 
1083 
1084  call atmos_boundary_var_fillhalo
1085  call atmos_boundary_alpha_fillhalo
1086 
1087  return
1088  end subroutine atmos_boundary_read
1089 
1090  !-----------------------------------------------------------------------------
1092  subroutine atmos_boundary_write
1093  use scale_file_cartesc, only: &
1094  file_cartesc_write
1095  use scale_comm_cartesc_nest, only: &
1097  implicit none
1098 
1099  integer :: iq
1100  !---------------------------------------------------------------------------
1101 
1102  if ( atmos_boundary_use_dens &
1103  .OR. atmos_boundary_use_velz &
1104  .OR. atmos_boundary_use_velx &
1105  .OR. atmos_boundary_use_vely &
1106  .OR. atmos_boundary_use_pott &
1107  ) then
1108  call file_cartesc_write( atmos_boundary_dens(:,:,:), &
1109  atmos_boundary_out_basename, atmos_boundary_out_title, &
1110  'DENS', 'Reference Density', 'kg/m3', 'ZXY', &
1111  atmos_boundary_out_dtype )
1112  end if
1113  if ( atmos_boundary_use_dens .OR. l_bnd ) then
1114  call file_cartesc_write( atmos_boundary_alpha_dens(:,:,:), &
1115  atmos_boundary_out_basename, atmos_boundary_out_title, &
1116  'ALPHA_DENS', 'Alpha for DENS', '1', 'ZXY', &
1117  atmos_boundary_out_dtype )
1118  endif
1119 
1120  if ( atmos_boundary_use_velz .OR. (l_bnd .AND. online_use_velz) ) then
1121  call file_cartesc_write( atmos_boundary_velz(:,:,:), &
1122  atmos_boundary_out_basename, atmos_boundary_out_title, &
1123  'VELZ', 'Reference Velocity w', 'm/s', 'ZHXY', &
1124  atmos_boundary_out_dtype )
1125  call file_cartesc_write( atmos_boundary_alpha_velz(:,:,:), &
1126  atmos_boundary_out_basename, atmos_boundary_out_title, &
1127  'ALPHA_VELZ', 'Alpha for VELZ', '1', 'ZHXY', &
1128  atmos_boundary_out_dtype )
1129  endif
1130 
1131  if ( atmos_boundary_use_velx .OR. l_bnd ) then
1132  call file_cartesc_write( atmos_boundary_velx(:,:,:), &
1133  atmos_boundary_out_basename, atmos_boundary_out_title, &
1134  'VELX', 'Reference Velocity u', 'm/s', 'ZXHY', &
1135  atmos_boundary_out_dtype )
1136  call file_cartesc_write( atmos_boundary_alpha_velx(:,:,:), &
1137  atmos_boundary_out_basename, atmos_boundary_out_title, &
1138  'ALPHA_VELX', 'Alpha for VELX', '1', 'ZXHY', &
1139  atmos_boundary_out_dtype )
1140  endif
1141 
1142  if ( atmos_boundary_use_vely .OR. l_bnd ) then
1143  call file_cartesc_write( atmos_boundary_vely(:,:,:), &
1144  atmos_boundary_out_basename, atmos_boundary_out_title, &
1145  'VELY', 'Reference Velocity y', 'm/s', 'ZXYH', &
1146  atmos_boundary_out_dtype )
1147  call file_cartesc_write( atmos_boundary_alpha_vely(:,:,:), &
1148  atmos_boundary_out_basename, atmos_boundary_out_title, &
1149  'ALPHA_VELY', 'Alpha for VELY', '1', 'ZXYH', &
1150  atmos_boundary_out_dtype )
1151  endif
1152 
1153  if ( atmos_boundary_use_pott .OR. l_bnd ) then
1154  call file_cartesc_write( atmos_boundary_pott(:,:,:), &
1155  atmos_boundary_out_basename, atmos_boundary_out_title, &
1156  'POTT', 'Reference POTT', 'K', 'ZXY', &
1157  atmos_boundary_out_dtype )
1158  call file_cartesc_write( atmos_boundary_alpha_pott(:,:,:), &
1159  atmos_boundary_out_basename, atmos_boundary_out_title, &
1160  'ALPHA_POTT', 'Alpha for POTT', '1', 'ZXY', &
1161  atmos_boundary_out_dtype )
1162  endif
1163 
1164  if ( atmos_boundary_use_qv .OR. l_bnd ) then
1165  call file_cartesc_write( atmos_boundary_qtrc(:,:,:,1), &
1166  atmos_boundary_out_basename, atmos_boundary_out_title, &
1167  'QV', 'Reference QV', 'kg/kg', 'ZXY', &
1168  atmos_boundary_out_dtype )
1169  call file_cartesc_write( atmos_boundary_alpha_qtrc(:,:,:,1), &
1170  atmos_boundary_out_basename, atmos_boundary_out_title, &
1171  'ALPHA_QV', 'Alpha for QV', '1', 'ZXY', &
1172  atmos_boundary_out_dtype )
1173  endif
1174 
1175  if ( atmos_boundary_use_qhyd ) then
1176  do iq = 2, bnd_qa
1177  call file_cartesc_write( atmos_boundary_qtrc(:,:,:,iq), &
1178  atmos_boundary_out_basename, atmos_boundary_out_title, &
1179  tracer_name(iq), 'Reference '//trim(tracer_name(iq)), &
1180  tracer_unit(iq), 'ZXY', &
1181  atmos_boundary_out_dtype )
1182  call file_cartesc_write( atmos_boundary_alpha_qtrc(:,:,:,iq), &
1183  atmos_boundary_out_basename, atmos_boundary_out_title, &
1184  'ALPHA_'//trim(tracer_name(iq)), 'Alpha for '//trim(tracer_name(iq)), &
1185  '1', 'ZXY', &
1186  atmos_boundary_out_dtype )
1187  end do
1188  endif
1189 
1190  return
1191  end subroutine atmos_boundary_write
1192 
1193  !-----------------------------------------------------------------------------
1195  subroutine atmos_boundary_generate
1196  use scale_atmos_refstate, only: &
1198  implicit none
1199 
1200  integer :: i, j, k, iq
1201  !---------------------------------------------------------------------------
1202 
1203  do j = 1, ja
1204  do i = 1, ia
1205  do k = 1, ka
1207  atmos_boundary_velz(k,i,j) = atmos_boundary_value_velz
1208  atmos_boundary_velx(k,i,j) = atmos_boundary_value_velx
1209  atmos_boundary_vely(k,i,j) = atmos_boundary_value_vely
1210  atmos_boundary_pott(k,i,j) = atmos_boundary_value_pott
1211  do iq = 1, bnd_qa
1212  atmos_boundary_qtrc(k,i,j,iq) = atmos_boundary_value_qtrc
1213  end do
1214  enddo
1215  enddo
1216  enddo
1217 
1218  call atmos_boundary_var_fillhalo
1219 
1220  return
1221  end subroutine atmos_boundary_generate
1222 
1223  !-----------------------------------------------------------------------------
1225  subroutine atmos_boundary_initialize_file
1226  use scale_calendar, only: &
1230  use scale_time, only: &
1231  time_nowdate
1232  use scale_file_cartesc, only: &
1234  file_cartesc_check_coordinates
1235  implicit none
1236 
1237  integer :: boundary_time_startday
1238  real(DP) :: boundary_time_startsec
1239  real(DP) :: boundary_time_startms
1240  integer :: boundary_time_offset_year
1241  character(len=27) :: boundary_chardate
1242 
1243  if ( atmos_boundary_start_date(1) == -9999 ) then
1244  atmos_boundary_start_date = time_nowdate
1245  end if
1246 
1247  !--- calculate time of the initial step in boundary file [no offset]
1248  boundary_time_startms = 0.0_dp
1249  boundary_time_offset_year = 0
1250  call calendar_date2char( boundary_chardate, & ! [OUT]
1251  atmos_boundary_start_date(:), & ! [IN]
1252  boundary_time_startms ) ! [IN]
1253 
1254  call calendar_date2daysec( boundary_time_startday, & ! [OUT]
1255  boundary_time_startsec, & ! [OUT]
1256  atmos_boundary_start_date(:), & ! [IN]
1257  boundary_time_startms, & ! [IN]
1258  boundary_time_offset_year ) ! [IN]
1259 
1260  boundary_time_initdaysec = calendar_combine_daysec( boundary_time_startday, boundary_time_startsec )
1261 
1262  log_info("ATMOS_BOUNDARY_initialize_file",'(1x,A,A)') 'BOUNDARY START Date : ', boundary_chardate
1263 
1264  call file_cartesc_open( atmos_boundary_in_basename, atmos_boundary_fid )
1265 
1266  if ( atmos_boundary_in_check_coordinates ) then
1267  call file_cartesc_check_coordinates( atmos_boundary_fid, atmos=.true. )
1268  end if
1269 
1270  return
1271  end subroutine atmos_boundary_initialize_file
1272 
1273  !-----------------------------------------------------------------------------
1275  subroutine atmos_boundary_set_file
1276  use scale_prc, only: &
1277  prc_abort
1278  use scale_time, only: &
1279  time_nowdate, &
1280  time_dtsec
1281  use scale_calendar, only: &
1284  implicit none
1285 
1286  real(RP) :: bnd_DENS(ka,ia,ja) ! damping coefficient for DENS (0-1)
1287  real(RP) :: bnd_VELZ(ka,ia,ja) ! damping coefficient for VELZ (0-1)
1288  real(RP) :: bnd_VELX(ka,ia,ja) ! damping coefficient for VELX (0-1)
1289  real(RP) :: bnd_VELY(ka,ia,ja) ! damping coefficient for VELY (0-1)
1290  real(RP) :: bnd_POTT(ka,ia,ja) ! damping coefficient for POTT (0-1)
1291  real(RP) :: bnd_QTRC(ka,ia,ja,bnd_qa) ! damping coefficient for QTRC (0-1)
1292 
1293  integer :: run_time_startdate(6)
1294  integer :: run_time_startday
1295  real(DP) :: run_time_startsec
1296  real(DP) :: run_time_startms
1297  integer :: run_time_offset_year
1298  real(DP) :: run_time_nowdaysec
1299 
1300  real(DP) :: boundary_diff_daysec
1301  real(RP) :: boundary_inc_offset
1302  integer :: fillgaps_steps
1303 
1304  integer :: i, j, k, iq
1305  !---------------------------------------------------------------------------
1306 
1307  if ( atmos_boundary_update_dt <= 0.0_dp ) then
1308  log_error("ATMOS_BOUNDARY_set_file",*) 'You need specify ATMOS_BOUNDARY_UPDATE_DT as larger than 0.0'
1309  call prc_abort
1310  endif
1311  update_nstep = nint( atmos_boundary_update_dt / time_dtsec )
1312  if ( abs(update_nstep * time_dtsec - atmos_boundary_update_dt) > 1e-10_dp ) then
1313  log_error("ATMOS_BOUNDARY_set_file",*) 'ATMOS_BOUNDARY_UPDATE_DT is not multiple of DT'
1314  call prc_abort
1315  end if
1316 
1317  !--- recalculate time of the run [no offset]
1318  run_time_startdate(:) = time_nowdate(:)
1319  run_time_startms = 0.0_dp
1320  run_time_offset_year = 0
1321 
1322  call calendar_date2daysec( run_time_startday, & ! [OUT]
1323  run_time_startsec, & ! [OUT]
1324  run_time_startdate(:), & ! [IN]
1325  run_time_startms, & ! [IN]
1326  run_time_offset_year ) ! [IN]
1327 
1328  run_time_nowdaysec = calendar_combine_daysec( run_time_startday, run_time_startsec )
1329 
1330  boundary_diff_daysec = run_time_nowdaysec - boundary_time_initdaysec
1331  boundary_timestep = 1 + int( boundary_diff_daysec / atmos_boundary_update_dt )
1332  boundary_inc_offset = mod( boundary_diff_daysec, atmos_boundary_update_dt )
1333  fillgaps_steps = int( boundary_inc_offset / time_dtsec )
1334 
1335  log_info("ATMOS_BOUNDARY_set_file",*) 'BOUNDARY TIMESTEP NUMBER FOR INIT:', boundary_timestep
1336  log_info("ATMOS_BOUNDARY_set_file",*) 'BOUNDARY OFFSET:', boundary_inc_offset
1337  log_info("ATMOS_BOUNDARY_set_file",*) 'BOUNDARY FILLGAPS STEPS:', fillgaps_steps
1338 
1339  ! read boundary data from input file
1340  call atmos_boundary_update_file( ref_now )
1341 
1342  boundary_timestep = boundary_timestep + 1
1343  call atmos_boundary_update_file( ref_new )
1344 
1345  ! copy now to old
1346  !$omp parallel do default(none) private(i,j,k,iq) OMP_SCHEDULE_ collapse(2) &
1347  !$omp shared(JA,IA,KA,ATMOS_BOUNDARY_ref_DENS,ref_old,ref_now,ATMOS_BOUNDARY_ref_VELX) &
1348  !$omp shared(ATMOS_BOUNDARY_ref_VELY,ATMOS_BOUNDARY_ref_POTT,BND_QA,ATMOS_BOUNDARY_ref_QTRC)
1349  do j = 1, ja
1350  do i = 1, ia
1351  do k = 1, ka
1352  atmos_boundary_ref_dens(k,i,j,ref_old) = atmos_boundary_ref_dens(k,i,j,ref_now)
1353  atmos_boundary_ref_velx(k,i,j,ref_old) = atmos_boundary_ref_velx(k,i,j,ref_now)
1354  atmos_boundary_ref_vely(k,i,j,ref_old) = atmos_boundary_ref_vely(k,i,j,ref_now)
1355  atmos_boundary_ref_pott(k,i,j,ref_old) = atmos_boundary_ref_pott(k,i,j,ref_now)
1356  do iq = 1, bnd_qa
1357  atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now)
1358  end do
1359  end do
1360  end do
1361  end do
1362 
1363  ! set boundary data
1364  do j = 1, ja
1365  do i = 1, ia
1366  do k = 1, ka
1367  atmos_boundary_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_now)
1368  atmos_boundary_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_now)
1369  atmos_boundary_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_now)
1370  atmos_boundary_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_now)
1371  do iq = 1, bnd_qa
1372  atmos_boundary_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now)
1373  end do
1374  end do
1375  end do
1376  end do
1377 
1378  if ( atmos_boundary_use_velz ) then
1379  do j = 1, ja
1380  do i = 1, ia
1381  do k = 1, ka
1382  atmos_boundary_velz(k,i,j) = atmos_boundary_value_velz
1383  end do
1384  end do
1385  end do
1386  end if
1387 
1388  now_step = fillgaps_steps
1389 
1390  ! get time boundary
1391  call get_boundary( bnd_dens(:,:,:), & ! [OUT]
1392  bnd_velz(:,:,:), & ! [OUT]
1393  bnd_velx(:,:,:), & ! [OUT]
1394  bnd_vely(:,:,:), & ! [OUT]
1395  bnd_pott(:,:,:), & ! [OUT]
1396  bnd_qtrc(:,:,:,:), & ! [OUT]
1397  now_step, & ! [IN]
1398  update_nstep ) ! [IN]
1399 
1400  ! fill in gaps of the offset
1401  do j = 1, ja
1402  do i = 1, ia
1403  do k = 1, ka
1404  atmos_boundary_dens(k,i,j) = bnd_dens(k,i,j)
1405  atmos_boundary_velx(k,i,j) = bnd_velx(k,i,j)
1406  atmos_boundary_vely(k,i,j) = bnd_vely(k,i,j)
1407  atmos_boundary_pott(k,i,j) = bnd_pott(k,i,j)
1408  do iq = 1, bnd_qa
1409  atmos_boundary_qtrc(k,i,j,iq) = bnd_qtrc(k,i,j,iq)
1410  end do
1411  end do
1412  end do
1413  end do
1414 
1415  return
1416  end subroutine atmos_boundary_set_file
1417 
1418  !-----------------------------------------------------------------------------
1420  subroutine atmos_boundary_initialize_online
1421  use scale_prc, only: &
1422  prc_abort
1423  use scale_comm_cartesc_nest, only: &
1425  online_use_velz, &
1426  parent_dtsec, &
1427  nestqa => comm_cartesc_nest_bnd_qa
1428  implicit none
1429 
1430  ! parameters
1431  integer, parameter :: handle = 2
1432 
1433  atmos_boundary_update_dt = parent_dtsec(handle)
1434 
1435  if ( nestqa > bnd_qa ) then
1436  log_error("ATMOS_BOUNDARY_initialize_online",*) 'NEST_BND_QA exceeds BND_QA'
1437  log_error_cont(*) 'This must not be occur.'
1438  log_error_cont(*) 'Please send your configuration file to SCALE develop member.'
1439  call prc_abort
1440  end if
1441 
1442  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
1443 
1444  return
1445  end subroutine atmos_boundary_initialize_online
1446 
1447  !-----------------------------------------------------------------------------
1449  subroutine atmos_boundary_set_online
1450  use scale_prc, only: &
1451  prc_abort
1452  use scale_time, only: &
1453  time_dtsec, &
1454  time_nstep
1455  use scale_comm_cartesc_nest, only: &
1456  online_use_velz, &
1457  parent_nstep
1458  implicit none
1459 
1460  ! parameters
1461  integer, parameter :: handle = 2
1462 
1463  ! works
1464  integer :: i, j, k, iq
1465  !---------------------------------------------------------------------------
1466 
1467  ! import data from parent domain
1468  boundary_timestep = 1
1469  log_info("ATMOS_BOUNDARY_set_online",*) 'BOUNDARY TIMESTEP NUMBER FOR INIT:', boundary_timestep
1470 
1471  call atmos_boundary_update_online_daughter( ref_now )
1472 
1473  boundary_timestep = boundary_timestep + 1
1474  log_info("ATMOS_BOUNDARY_set_online",*) 'BOUNDARY TIMESTEP NUMBER FOR INIT:', boundary_timestep
1475 
1476  call atmos_boundary_update_online_daughter( ref_new )
1477 
1478  ! copy now to old
1479  do j = 1, ja
1480  do i = 1, ia
1481  do k = 1, ka
1482  atmos_boundary_ref_dens(k,i,j,ref_old) = atmos_boundary_ref_dens(k,i,j,ref_now)
1483  atmos_boundary_ref_velx(k,i,j,ref_old) = atmos_boundary_ref_velx(k,i,j,ref_now)
1484  atmos_boundary_ref_vely(k,i,j,ref_old) = atmos_boundary_ref_vely(k,i,j,ref_now)
1485  atmos_boundary_ref_pott(k,i,j,ref_old) = atmos_boundary_ref_pott(k,i,j,ref_now)
1486  do iq = 1, bnd_qa
1487  atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now)
1488  end do
1489  end do
1490  end do
1491  end do
1492 
1493  ! set boundary data
1494  do j = 1, ja
1495  do i = 1, ia
1496  do k = 1, ka
1497  atmos_boundary_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_now)
1498  atmos_boundary_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_now)
1499  atmos_boundary_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_now)
1500  atmos_boundary_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_now)
1501  do iq = 1, bnd_qa
1502  atmos_boundary_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now)
1503  end do
1504  end do
1505  end do
1506  end do
1507 
1508  if ( online_use_velz ) then
1509  do j = 1, ja
1510  do i = 1, ia
1511  do k = 1, ka
1512  atmos_boundary_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_now)
1513  end do
1514  end do
1515  end do
1516  else if ( atmos_boundary_use_velz ) then
1517  do j = 1, ja
1518  do i = 1, ia
1519  do k = 1, ka
1520  atmos_boundary_velz(k,i,j) = atmos_boundary_value_velz
1521  end do
1522  end do
1523  end do
1524  end if
1525 
1526  update_nstep = nint( atmos_boundary_update_dt / time_dtsec )
1527  if ( update_nstep * time_dtsec /= atmos_boundary_update_dt ) then
1528  log_error("ATMOS_BOUNDARY_set_online",*) 'DT of the parent is not multiple of the DT'
1529  call prc_abort
1530  end if
1531  if ( update_nstep * parent_nstep(handle) /= time_nstep ) then
1532  log_error("ATMOS_BOUNDARY_set_online",*) 'DURATION must be the same as that of the parent'
1533  call prc_abort
1534  end if
1535 
1536  now_step = 0 ! should be set as zero in initialize process
1537 
1538  return
1539  end subroutine atmos_boundary_set_online
1540 
1541  !-----------------------------------------------------------------------------
1543  subroutine atmos_boundary_firstsend( &
1544  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, QV, Qe )
1545  use scale_comm_cartesc_nest, only: &
1546  nestqa => comm_cartesc_nest_bnd_qa
1547  use scale_atmos_hydrometeor, only: &
1548  n_hyd
1549  use mod_atmos_phy_mp_vars, only: &
1550  qa_mp
1551  implicit none
1552 
1553  ! arguments
1554  real(RP), intent(in) :: DENS(ka,ia,ja)
1555  real(RP), intent(in) :: MOMZ(ka,ia,ja)
1556  real(RP), intent(in) :: MOMX(ka,ia,ja)
1557  real(RP), intent(in) :: MOMY(ka,ia,ja)
1558  real(RP), intent(in) :: RHOT(ka,ia,ja)
1559  real(RP), intent(in) :: QTRC(ka,ia,ja,qa_mp)
1560  real(RP), intent(in) :: QV (ka,ia,ja)
1561  real(RP), intent(in) :: Qe (ka,ia,ja,n_hyd)
1562  !---------------------------------------------------------------------------
1563 
1564  ! send data at the first time
1565  if ( do_parent_process ) then !online [parent]
1566  ! issue send
1567  call atmos_boundary_send( dens, momz, momx, momy, rhot, qtrc, qv, qe )
1568  endif
1569 
1570  return
1571  end subroutine atmos_boundary_firstsend
1572 
1573  !-----------------------------------------------------------------------------
1579  nestqa => comm_cartesc_nest_bnd_qa
1580  use scale_file_cartesc, only: &
1582  implicit none
1583 
1584  ! works
1585  integer :: handle
1586  !---------------------------------------------------------------------------
1587 
1588  if ( do_parent_process ) then !online [parent]
1589  handle = 1
1590  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
1591  endif
1592 
1593  if ( do_daughter_process ) then !online [daughter]
1594  handle = 2
1595  call comm_cartesc_nest_recv_cancel( handle )
1596  endif
1597 
1598  if ( atmos_boundary_fid > 0 ) then
1599  call file_cartesc_close( atmos_boundary_fid )
1600  atmos_boundary_fid = -1
1601  end if
1602 
1603  return
1604  end subroutine atmos_boundary_driver_finalize
1605 
1606  !-----------------------------------------------------------------------------
1608  subroutine atmos_boundary_driver_update
1609  use scale_prc, only: &
1610  prc_abort
1611  use scale_prc_cartesc, only: &
1612  prc_has_w, &
1613  prc_has_e, &
1614  prc_has_s, &
1615  prc_has_n
1616  use scale_comm_cartesc_nest, only: &
1617  online_use_velz, &
1619  use mod_atmos_vars, only: &
1620  dens, &
1621  momz, &
1622  momx, &
1623  momy, &
1624  rhot, &
1625  qtrc, &
1626  qv, &
1627  qe
1628  use mod_atmos_phy_mp_vars, only: &
1629  qs_mp, &
1630  qe_mp
1631  implicit none
1632 
1633  real(RP) :: bnd_DENS(ka,ia,ja) ! damping coefficient for DENS (0-1)
1634  real(RP) :: bnd_VELZ(ka,ia,ja) ! damping coefficient for VELZ (0-1)
1635  real(RP) :: bnd_VELX(ka,ia,ja) ! damping coefficient for VELX (0-1)
1636  real(RP) :: bnd_VELY(ka,ia,ja) ! damping coefficient for VELY (0-1)
1637  real(RP) :: bnd_POTT(ka,ia,ja) ! damping coefficient for POTT (0-1)
1638  real(RP) :: bnd_QTRC(ka,ia,ja,bnd_qa) ! damping coefficient for QTRC (0-1)
1639 
1640  integer :: handle
1641  integer :: i, j, k, iq, iqa
1642  !---------------------------------------------------------------------------
1643 
1644  if ( do_parent_process ) then !online [parent]
1645  ! should be called every time step
1646  call atmos_boundary_update_online_parent( dens,momz,momx,momy,rhot,qtrc(:,:,:,qs_mp:qe_mp), qv, qe )
1647  endif
1648 
1649  if ( l_bnd ) then
1650  ! update referce vars
1651  if ( now_step >= update_nstep ) then
1652  now_step = 0
1653  boundary_timestep = boundary_timestep + 1
1654 
1655  call update_ref_index
1656 
1657  if ( do_daughter_process ) then !online [daughter]
1658  call atmos_boundary_update_online_daughter( ref_new )
1659  else
1660  call atmos_boundary_update_file( ref_new )
1661  end if
1662  end if
1663 
1664  ! step boundary
1665  now_step = now_step + 1
1666 
1667  ! get boundaryal coefficients
1668  call get_boundary( bnd_dens(:,:,:), & ! [OUT]
1669  bnd_velz(:,:,:), & ! [OUT]
1670  bnd_velx(:,:,:), & ! [OUT]
1671  bnd_vely(:,:,:), & ! [OUT]
1672  bnd_pott(:,:,:), & ! [OUT]
1673  bnd_qtrc(:,:,:,:), & ! [OUT]
1674  now_step, & ! [IN]
1675  update_nstep ) ! [IN]
1676 
1677  ! update boundary vars
1678  do j = 1, ja
1679  do i = 1, ia
1680  do k = 1, ka
1681  atmos_boundary_dens(k,i,j) = bnd_dens(k,i,j)
1682  atmos_boundary_velx(k,i,j) = bnd_velx(k,i,j)
1683  atmos_boundary_vely(k,i,j) = bnd_vely(k,i,j)
1684  atmos_boundary_pott(k,i,j) = bnd_pott(k,i,j)
1685  do iq = 1, bnd_qa
1686  atmos_boundary_qtrc(k,i,j,iq) = bnd_qtrc(k,i,j,iq)
1687  end do
1688  end do
1689  end do
1690  end do
1691  if ( online_use_velz ) then
1692  do j = 1, ja
1693  do i = 1, ia
1694  do k = 1, ka
1695  atmos_boundary_velz(k,i,j) = bnd_velz(k,i,j)
1696  end do
1697  end do
1698  end do
1699  end if
1700 
1701  ! fill HALO in western region
1702  if ( .NOT. prc_has_w ) then
1703  !$omp parallel do default(none) &
1704  !$omp shared(JS,IS,KA,DENS,ATMOS_BOUNDARY_DENS,MOMX,ATMOS_BOUNDARY_VELX,RHOT) &
1705  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC,BND_QA,QS_MP,QTRC,QA,JA) &
1706  !$omp private(i,j,k,iq,iqa) OMP_SCHEDULE_ collapse(2)
1707  do j = 1, ja
1708  do i = 1, is-1
1709  do k = 1, ka
1710  dens(k,i,j) = atmos_boundary_dens(k,i,j)
1711  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
1712  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
1713  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
1714  do iq = 1, bnd_qa
1715  iqa = iq + qs_mp - 1
1716  qtrc(k,i,j,iqa) = atmos_boundary_qtrc(k,i,j,iq)
1717  end do
1718  do iq = 1, qa
1719  if ( iq < qs_mp .or. iq >= bnd_qa + qs_mp ) then
1720  qtrc(k,i,j,iq) = qtrc(k,is,j,iq) &
1721  * ( 0.5_rp - sign(0.5_rp, atmos_boundary_velx(k,is-1,j)) )
1722  end if
1723  end do
1724  end do
1725  end do
1726  end do
1727  do j = 1, ja-1
1728  do i = 1, is-1
1729  do k = 1, ka
1730  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
1731  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
1732  end do
1733  end do
1734  end do
1735  do i = 1, is-1
1736  do k = 1, ka
1737  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) &
1738  * atmos_boundary_dens(k,i,ja)
1739  end do
1740  end do
1741  if ( online_use_velz ) then
1742  do j = 1, ja
1743  do i = 1, is-1
1744  do k = ks, ke-1
1745  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
1746  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
1747  end do
1748  end do
1749  end do
1750  else
1751  do j = 1, ja
1752  do i = 1, is-1
1753  do k = ks, ke-1
1754  momz(k,i,j) = momz(k,is,j)
1755  end do
1756  end do
1757  end do
1758  end if
1759  end if
1760 
1761  ! fill HALO in eastern region
1762  if ( .NOT. prc_has_e ) then
1763  !$omp parallel do default(none) &
1764  !$omp shared(JA,IE,IA,KA,DENS,ATMOS_BOUNDARY_DENS,ATMOS_BOUNDARY_VELX,RHOT) &
1765  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC,BND_QA,QS_MP,QTRC,QA) &
1766  !$omp private(i,j,k,iq,iqa) OMP_SCHEDULE_ collapse(2)
1767  do j = 1, ja
1768  do i = ie+1, ia
1769  do k = 1, ka
1770  dens(k,i,j) = atmos_boundary_dens(k,i,j)
1771  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
1772  do iq = 1, bnd_qa
1773  iqa = iq + qs_mp - 1
1774  qtrc(k,i,j,iqa) = atmos_boundary_qtrc(k,i,j,iq)
1775  end do
1776  do iq = 1, qa
1777  if ( iq < qs_mp .or. iq >= bnd_qa + qs_mp ) then
1778  qtrc(k,i,j,iq) = qtrc(k,ie,j,iq) &
1779  * ( 0.5_rp + sign(0.5_rp, atmos_boundary_velx(k,ie,j)) )
1780  end if
1781  end do
1782  end do
1783  end do
1784  end do
1785  do j = 1, ja
1786  do i = ie, ia-1
1787  do k = 1, ka
1788  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
1789  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
1790  end do
1791  end do
1792  end do
1793  do j = 1, ja
1794  do k = 1, ka
1795  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) * atmos_boundary_dens(k,ia,j)
1796  end do
1797  end do
1798  do j = 1, ja-1
1799  do i = ie+1, ia
1800  do k = 1, ka
1801  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
1802  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
1803  end do
1804  end do
1805  end do
1806  do i = ie+1, ia
1807  do k = 1, ka
1808  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) &
1809  * atmos_boundary_dens(k,i,ja)
1810  end do
1811  end do
1812  if ( online_use_velz ) then
1813  do j = 1, ja
1814  do i = ie+1, ia
1815  do k = ks, ke-1
1816  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
1817  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
1818  end do
1819  end do
1820  end do
1821  else
1822  do j = 1, ja
1823  do i = ie+1, ia
1824  do k = ks, ke-1
1825  momz(k,i,j) = momz(k,ie,j)
1826  end do
1827  end do
1828  end do
1829  end if
1830  end if
1831 
1832  ! fill HALO in southern region
1833  if ( .NOT. prc_has_s ) then
1834  !$omp parallel do default(none) &
1835  !$omp shared(JS,IA,KA,DENS,ATMOS_BOUNDARY_DENS,MOMY,ATMOS_BOUNDARY_VELY,RHOT) &
1836  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC,BND_QA,QS_MP,QTRC,QA) &
1837  !$omp private(i,j,k,iq,iqa) OMP_SCHEDULE_ collapse(2)
1838  do j = 1, js-1
1839  do i = 1, ia
1840  do k = 1, ka
1841  dens(k,i,j) = atmos_boundary_dens(k,i,j)
1842  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
1843  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
1844  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
1845  do iq = 1, bnd_qa
1846  iqa = iq + qs_mp - 1
1847  qtrc(k,i,j,iqa) = atmos_boundary_qtrc(k,i,j,iq)
1848  end do
1849  do iq = 1, qa
1850  if ( iq < qs_mp .or. iq >= bnd_qa + qs_mp ) then
1851  qtrc(k,i,j,iq) = qtrc(k,i,js,iq) &
1852  * ( 0.5_rp - sign(0.5_rp, atmos_boundary_vely(k,i,js-1)) )
1853  end if
1854  end do
1855  end do
1856  end do
1857  end do
1858  do j = 1, js-1
1859  do i = 1, ia-1
1860  do k = 1, ka
1861  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
1862  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
1863  end do
1864  end do
1865  end do
1866  do j = 1, js-1
1867  do k = 1, ka
1868  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) &
1869  * atmos_boundary_dens(k,ia,j)
1870  end do
1871  end do
1872  if ( online_use_velz ) then
1873  do j = 1, js-1
1874  do i = 1, ia
1875  do k = ks, ke-1
1876  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
1877  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
1878  end do
1879  end do
1880  end do
1881  else
1882  do j = 1, js-1
1883  do i = 1, ia
1884  do k = ks, ke-1
1885  momz(k,i,j) = momz(k,i,js)
1886  end do
1887  end do
1888  end do
1889  end if
1890  end if
1891 
1892  ! fill HALO in northern region
1893  if ( .NOT. prc_has_n ) then
1894  !$omp parallel do default(none) &
1895  !$omp shared(JE,JA,IA,KA,DENS,ATMOS_BOUNDARY_DENS,ATMOS_BOUNDARY_VELY,RHOT) &
1896  !$omp shared(ATMOS_BOUNDARY_POTT,ATMOS_BOUNDARY_QTRC,BND_QA,QS_MP,QTRC,QA) &
1897  !$omp private(i,j,k,iq,iqa) OMP_SCHEDULE_ collapse(2)
1898  do j = je+1, ja
1899  do i = 1, ia
1900  do k = 1, ka
1901  dens(k,i,j) = atmos_boundary_dens(k,i,j)
1902  rhot(k,i,j) = atmos_boundary_pott(k,i,j) * atmos_boundary_dens(k,i,j)
1903  do iq = 1, bnd_qa
1904  iqa = iq + qs_mp - 1
1905  qtrc(k,i,j,iqa) = atmos_boundary_qtrc(k,i,j,iq)
1906  end do
1907  do iq = bnd_qa+1, qa
1908  if ( iq < qs_mp .or. iq >= bnd_qa + qs_mp ) then
1909  qtrc(k,i,j,iq) = qtrc(k,i,je,iq) &
1910  * ( 0.5_rp + sign(0.5_rp, atmos_boundary_vely(k,i,je)) )
1911  end if
1912  end do
1913  end do
1914  end do
1915  end do
1916  do j = je, ja-1
1917  do i = 1, ia
1918  do k = 1, ka
1919  momy(k,i,j) = atmos_boundary_vely(k,i,j) &
1920  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i,j+1) ) * 0.5_rp
1921  end do
1922  end do
1923  end do
1924  do i = 1, ia
1925  do k = 1, ka
1926  momy(k,i,ja) = atmos_boundary_vely(k,i,ja) * atmos_boundary_dens(k,i,ja)
1927  end do
1928  end do
1929  do j = je+1, ja
1930  do i = 1, ia-1
1931  do k = 1, ka
1932  momx(k,i,j) = atmos_boundary_velx(k,i,j) &
1933  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k,i+1,j) ) * 0.5_rp
1934  end do
1935  end do
1936  end do
1937  do j = je+1, ja
1938  do k = 1, ka
1939  momx(k,ia,j) = atmos_boundary_velx(k,ia,j) &
1940  * atmos_boundary_dens(k,ia,j)
1941  end do
1942  end do
1943  if ( online_use_velz ) then
1944  do j = je+1, ja
1945  do i = 1, ia
1946  do k = ks, ke-1
1947  momz(k,i,j) = atmos_boundary_velz(k,i,j) &
1948  * ( atmos_boundary_dens(k,i,j) + atmos_boundary_dens(k+1,i,j) ) * 0.5_rp
1949  end do
1950  end do
1951  end do
1952  else
1953  do j = je+1, ja
1954  do i = 1, ia
1955  do k = ks, ke-1
1956  momz(k,i,j) = momz(k,i,je)
1957  end do
1958  end do
1959  end do
1960  end if
1961  end if
1962 
1963  elseif ( do_parent_process ) then
1964  ! do nothing
1965  else
1966  log_error("ATMOS_BOUNDARY_update",*) '[BUG] invalid path'
1967  call prc_abort
1968  end if
1969 
1970  call history_bnd( atmos_boundary_dens, &
1976 
1977  ! To be enable to do asynchronous communicaton
1978  if ( do_parent_process ) then !online [parent]
1979  handle = 1
1980  call comm_cartesc_nest_test( handle )
1981  endif
1982  if ( do_daughter_process ) then !online [daughter]
1983  handle = 2
1984  call comm_cartesc_nest_test( handle )
1985  endif
1986 
1987  return
1988  end subroutine atmos_boundary_driver_update
1989 
1990  !-----------------------------------------------------------------------------
1992  subroutine atmos_boundary_update_file( ref )
1993  use scale_file_cartesc, only: &
1994  file_cartesc_read, &
1996  use mod_atmos_phy_mp_vars, only: &
1997  qs_mp
1998  implicit none
1999 
2000  integer, intent(in) :: ref
2001 
2002  integer :: fid, iq
2003  !---------------------------------------------------------------------------
2004 
2005  log_info("ATMOS_BOUNDARY_update_file",*) "Atmos Boundary: read from boundary file(timestep=", boundary_timestep, ")"
2006 
2007  fid = atmos_boundary_fid
2008 
2009  call file_cartesc_read( fid, 'DENS', 'ZXY', atmos_boundary_ref_dens(:,:,:,ref), step=boundary_timestep )
2010  call file_cartesc_read( fid, 'VELX', 'ZXHY', atmos_boundary_ref_velx(:,:,:,ref), step=boundary_timestep )
2011  call file_cartesc_read( fid, 'VELY', 'ZXYH', atmos_boundary_ref_vely(:,:,:,ref), step=boundary_timestep )
2012  call file_cartesc_read( fid, 'POTT', 'ZXY', atmos_boundary_ref_pott(:,:,:,ref), step=boundary_timestep )
2013  do iq = 1, bnd_qa
2014  call file_cartesc_read( fid, tracer_name(qs_mp+iq-1), 'ZXY', atmos_boundary_ref_qtrc(:,:,:,iq,ref), step=boundary_timestep )
2015  end do
2016 
2017  call file_cartesc_flush( fid )
2018 
2019  ! fill HALO in reference
2020  call atmos_boundary_ref_fillhalo( ref )
2021 
2022  return
2023  end subroutine atmos_boundary_update_file
2024 
2025  !-----------------------------------------------------------------------------
2027  subroutine atmos_boundary_update_online_parent( &
2028  DENS, & ! [in]
2029  MOMZ, & ! [in]
2030  MOMX, & ! [in]
2031  MOMY, & ! [in]
2032  RHOT, & ! [in]
2033  QTRC, & ! [in]
2034  QV, & ! [in]
2035  Qe ) ! [in]
2036  use scale_comm_cartesc_nest, only: &
2038  nestqa => comm_cartesc_nest_bnd_qa
2039  use scale_atmos_hydrometeor, only: &
2040  n_hyd
2041  use mod_atmos_phy_mp_vars, only: &
2042  qa_mp
2043  implicit none
2044 
2045  ! arguments
2046  real(RP), intent(in) :: DENS(ka,ia,ja)
2047  real(RP), intent(in) :: MOMZ(ka,ia,ja)
2048  real(RP), intent(in) :: MOMX(ka,ia,ja)
2049  real(RP), intent(in) :: MOMY(ka,ia,ja)
2050  real(RP), intent(in) :: RHOT(ka,ia,ja)
2051  real(RP), intent(in) :: QTRC(ka,ia,ja,qa_mp)
2052  real(RP), intent(in) :: QV (ka,ia,ja)
2053  real(RP), intent(in) :: Qe (ka,ia,ja,n_hyd)
2054 
2055  integer, parameter :: handle = 1
2056  !---------------------------------------------------------------------------
2057 
2058  log_info("ATMOS_BOUNDARY_update_online_parent",*)"ATMOS BOUNDARY update online: PARENT"
2059 
2060  ! issue wait
2061  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
2062 
2063  ! issue send
2064  call atmos_boundary_send( dens, momz, momx, momy, rhot, qtrc, qv, qe )
2065 
2066  return
2067  end subroutine atmos_boundary_update_online_parent
2068 
2069  !-----------------------------------------------------------------------------
2071  subroutine atmos_boundary_update_online_daughter( &
2072  ref ) ! [in]
2073  use scale_comm_cartesc_nest, only: &
2075  nestqa => comm_cartesc_nest_bnd_qa
2076  implicit none
2077 
2078  ! arguments
2079  integer, intent(in) :: ref
2080 
2081  integer, parameter :: handle = 2
2082  !---------------------------------------------------------------------------
2083 
2084  log_info("ATMOS_BOUNDARY_update_online_daughter",'(1X,A,I5)') 'ATMOS BOUNDARY update online: DAUGHTER', boundary_timestep
2085 
2086  ! issue wait
2087  call atmos_boundary_recv( ref )
2088 
2089  ! fill HALO in reference
2090  call atmos_boundary_ref_fillhalo( ref )
2091 
2092  ! issue receive
2093  call comm_cartesc_nest_recvwait_issue( handle, nestqa )
2094 
2095  return
2096  end subroutine atmos_boundary_update_online_daughter
2097 
2098  !-----------------------------------------------------------------------------
2100  subroutine atmos_boundary_send( &
2101  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, QV, Qe )
2102  use scale_comm_cartesc_nest, only: &
2105  daughter_ka, &
2106  daughter_ia, &
2107  daughter_ja, &
2108  nestqa => comm_cartesc_nest_bnd_qa
2109  use scale_atmos_hydrometeor, only: &
2110  n_hyd
2111  use mod_atmos_phy_mp_vars, only: &
2112  qa_mp
2113  implicit none
2114 
2115  ! parameters
2116  integer, parameter :: handle = 1
2117 
2118  ! arguments
2119  real(RP), intent(in) :: DENS(ka,ia,ja)
2120  real(RP), intent(in) :: MOMZ(ka,ia,ja)
2121  real(RP), intent(in) :: MOMX(ka,ia,ja)
2122  real(RP), intent(in) :: MOMY(ka,ia,ja)
2123  real(RP), intent(in) :: RHOT(ka,ia,ja)
2124  real(RP), intent(in), target :: QTRC(ka,ia,ja,qa_mp)
2125  real(RP), intent(in) :: QV (ka,ia,ja)
2126  real(RP), intent(in) :: Qe (ka,ia,ja,n_hyd)
2127 
2128  ! works
2129  real(RP), pointer :: Q(:,:,:,:)
2130  real(RP) :: dummy_d( daughter_ka(handle), daughter_ia(handle), daughter_ja(handle), nestqa )
2131 
2132  integer :: iq
2133  !---------------------------------------------------------------------------
2134 
2135  if ( online_boundary_diagqhyd ) then
2136  q => q_work
2137  q(:,:,:,1) = qv(:,:,:)
2138  do iq = 2, nestqa
2139  q(:,:,:,iq) = qe(:,:,:,iq-1)
2140  end do
2141  else
2142  q => qtrc
2143  end if
2144 
2145  call comm_cartesc_nest_nestdown( handle, &
2146  nestqa, &
2147  dens(:,:,:), & !(KA,IA,JA)
2148  momz(:,:,:), & !(KA,IA,JA)
2149  momx(:,:,:), & !(KA,IA,JA)
2150  momy(:,:,:), & !(KA,IA,JA)
2151  rhot(:,:,:), & !(KA,IA,JA)
2152  q(:,:,:,:), & !(KA,IA,JA,NESTQA)
2153  dummy_d(:,:,:,1), & !(KA,IA,JA)
2154  dummy_d(:,:,:,1), & !(KA,IA,JA)
2155  dummy_d(:,:,:,1), & !(KA,IA,JA)
2156  dummy_d(:,:,:,1), & !(KA,IA,JA)
2157  dummy_d(:,:,:,1), & !(KA,IA,JA)
2158  dummy_d(:,:,:,:) ) !(KA,IA,JA,NESTQA)
2159 
2160  return
2161  end subroutine atmos_boundary_send
2162 
2163  !-----------------------------------------------------------------------------
2165  subroutine atmos_boundary_recv( &
2166  ref_idx )
2167  use scale_prc, only: &
2168  prc_abort
2169  use scale_comm_cartesc_nest, only: &
2172  parent_ka, &
2173  parent_ia, &
2174  parent_ja, &
2175  nestqa => comm_cartesc_nest_bnd_qa
2176  use mod_atmos_phy_mp_driver, only: &
2178  implicit none
2179 
2180  ! parameters
2181  integer, parameter :: handle = 2
2182 
2183  ! arguments
2184  integer, intent(in) :: ref_idx
2185 
2186  ! works
2187  real(RP), pointer :: Q(:,:,:,:)
2188  real(RP) :: dummy_p( parent_ka(handle), parent_ia(handle), parent_ja(handle), nestqa )
2189  !---------------------------------------------------------------------------
2190 
2191 !OCL XFILL
2192  dummy_p(:,:,:,:) = 0.0_rp
2193 
2194  if ( online_boundary_diagqhyd ) then
2195  q => q_work
2196  else
2197  q => atmos_boundary_ref_qtrc(:,:,:,1:nestqa,ref_idx)
2198  end if
2199 
2200  call comm_cartesc_nest_nestdown( handle, &
2201  nestqa, &
2202  dummy_p(:,:,:,1), & !(KA,IA,JA)
2203  dummy_p(:,:,:,1), & !(KA,IA,JA)
2204  dummy_p(:,:,:,1), & !(KA,IA,JA)
2205  dummy_p(:,:,:,1), & !(KA,IA,JA)
2206  dummy_p(:,:,:,1), & !(KA,IA,JA)
2207  dummy_p(:,:,:,1:nestqa), & !(KA,IA,JA,NESTQA)
2208  atmos_boundary_ref_dens(:,:,:,ref_idx), & !(KA,IA,JA)
2209  atmos_boundary_ref_velz(:,:,:,ref_idx), & !(KA,IA,JA)
2210  atmos_boundary_ref_velx(:,:,:,ref_idx), & !(KA,IA,JA)
2211  atmos_boundary_ref_vely(:,:,:,ref_idx), & !(KA,IA,JA)
2212  atmos_boundary_ref_pott(:,:,:,ref_idx), & !(KA,IA,JA)
2213  q(:,:,:,:) ) !(KA,IA,JA,NESTQA)
2214 
2215 
2216  if ( online_boundary_diagqhyd ) then
2217  call atmos_phy_mp_driver_qhyd2qtrc( ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2218  q(:,:,:,1), q(:,:,:,2:), & ! (in)
2219  atmos_boundary_ref_qtrc(:,:,:,:,ref_idx) ) ! (out)
2220  else if ( nestqa > 0 ) then
2221  atmos_boundary_ref_qtrc(:,:,:,1,ref_idx) = q(:,:,:,1)
2222  end if
2223 
2224  return
2225  end subroutine atmos_boundary_recv
2226 
2227  !-----------------------------------------------------------------------------
2229  subroutine get_boundary_same_parent( &
2230  bnd_DENS, &
2231  bnd_VELZ, &
2232  bnd_VELX, &
2233  bnd_VELY, &
2234  bnd_POTT, &
2235  bnd_QTRC, &
2236  now_step, &
2237  update_step )
2238  implicit none
2239 
2240  ! arguments
2241  real(RP), intent(out) :: bnd_DENS(:,:,:)
2242  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2243  real(RP), intent(out) :: bnd_VELX(:,:,:)
2244  real(RP), intent(out) :: bnd_VELY(:,:,:)
2245  real(RP), intent(out) :: bnd_POTT(:,:,:)
2246  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2247  integer, intent(in) :: now_step
2248  integer, intent(in) :: update_step
2249 
2250  ! works
2251  integer :: i, j, k, iq
2252  integer :: ref
2253  !---------------------------------------------------------------------------
2254 
2255  if ( now_step == update_step ) then
2256  ref = ref_new
2257  else
2258  ref = ref_now
2259  end if
2260 
2261  do j = 1, ja
2262  do i = 1, ia
2263  do k = 1, ka
2264  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref)
2265  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref)
2266  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref)
2267  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref)
2268  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref)
2269  do iq = 1, bnd_qa
2270  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref)
2271  end do
2272  end do
2273  end do
2274  end do
2275 
2276  return
2277  end subroutine get_boundary_same_parent
2278 
2279  !-----------------------------------------------------------------------------
2281  subroutine get_boundary_nearest_neighbor( &
2282  bnd_DENS, &
2283  bnd_VELZ, &
2284  bnd_VELX, &
2285  bnd_VELY, &
2286  bnd_POTT, &
2287  bnd_QTRC, &
2288  now_step, &
2289  update_step )
2290  implicit none
2291 
2292  ! parameters
2293  real(RP) :: EPS = 1.0e-4_rp
2294 
2295  ! arguments
2296  real(RP), intent(out) :: bnd_DENS(:,:,:)
2297  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2298  real(RP), intent(out) :: bnd_VELX(:,:,:)
2299  real(RP), intent(out) :: bnd_VELY(:,:,:)
2300  real(RP), intent(out) :: bnd_POTT(:,:,:)
2301  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2302  integer, intent(in) :: now_step
2303  integer, intent(in) :: update_step
2304 
2305  ! works
2306  integer :: i, j, k, iq
2307  integer :: ref_idx
2308 
2309  real(RP) :: real_nstep
2310  real(RP) :: half_nstep
2311  !---------------------------------------------------------------------------
2312 
2313  real_nstep = real( now_step, kind=rp )
2314  half_nstep = real( UPDATE_NSTEP, kind=RP ) * 0.5_rp
2315 
2316  ! this step before half of the parent step
2317  if( ( real_nstep - eps ) < half_nstep ) then
2318  ref_idx = ref_now
2319 
2320  ! this step after half of the parent step
2321  else if( ( real_nstep - 1.0_rp + eps ) > half_nstep ) then
2322  ref_idx = ref_new
2323 
2324  ! this step across half of the parent step
2325  else
2326  ref_idx = ref_now
2327 
2328  end if
2329 
2330  do j = 1, ja
2331  do i = 1, ia
2332  do k = 1, ka
2333  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_idx)
2334  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_idx)
2335  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_idx)
2336  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_idx)
2337  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_idx)
2338  do iq = 1, bnd_qa
2339  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_idx)
2340  end do
2341  end do
2342  end do
2343  end do
2344 
2345  return
2346  end subroutine get_boundary_nearest_neighbor
2347 
2348  !-----------------------------------------------------------------------------
2350  subroutine get_boundary_lerp_initpoint( &
2351  bnd_DENS, &
2352  bnd_VELZ, &
2353  bnd_VELX, &
2354  bnd_VELY, &
2355  bnd_POTT, &
2356  bnd_QTRC, &
2357  now_step, &
2358  update_step )
2359  implicit none
2360 
2361  ! arguments
2362  real(RP), intent(out) :: bnd_DENS(:,:,:)
2363  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2364  real(RP), intent(out) :: bnd_VELX(:,:,:)
2365  real(RP), intent(out) :: bnd_VELY(:,:,:)
2366  real(RP), intent(out) :: bnd_POTT(:,:,:)
2367  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2368  integer, intent(in) :: now_step
2369  integer, intent(in) :: update_step
2370 
2371  ! works
2372  integer :: i, j, k, iq
2373 
2374  real(RP) :: fact
2375  !---------------------------------------------------------------------------
2376 
2377  fact = REAL(now_step, kind=RP) / update_step
2378 
2379  !$omp parallel do default(none) private(i,j,k,iq) OMP_SCHEDULE_ collapse(2) &
2380  !$omp shared(JA,IA,KA,bnd_DENS,ATMOS_BOUNDARY_ref_DENS,ref_now,fact,ref_new,bnd_VELZ) &
2381  !$omp shared(ATMOS_BOUNDARY_ref_VELZ,bnd_VELX,ATMOS_BOUNDARY_ref_VELX,bnd_VELY) &
2382  !$omp shared(ATMOS_BOUNDARY_ref_VELY,bnd_POTT,ATMOS_BOUNDARY_ref_POTT,BND_QA,bnd_QTRC,ATMOS_BOUNDARY_ref_QTRC)
2383  do j = 1, ja
2384  do i = 1, ia
2385  do k = 1, ka
2386  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_now) * ( 1.0_rp-fact ) &
2387  + atmos_boundary_ref_dens(k,i,j,ref_new) * fact
2388  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_now) * ( 1.0_rp-fact ) &
2389  + atmos_boundary_ref_velz(k,i,j,ref_new) * fact
2390  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_now) * ( 1.0_rp-fact) &
2391  + atmos_boundary_ref_velx(k,i,j,ref_new) * fact
2392  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_now) * ( 1.0_rp-fact) &
2393  + atmos_boundary_ref_vely(k,i,j,ref_new) * fact
2394  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_now) * ( 1.0_rp-fact ) &
2395  + atmos_boundary_ref_pott(k,i,j,ref_new) * fact
2396  do iq = 1, bnd_qa
2397  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * ( 1.0_rp-fact ) &
2398  + atmos_boundary_ref_qtrc(k,i,j,iq,ref_new) * fact
2399  end do
2400  end do
2401  end do
2402  end do
2403 
2404  return
2405  end subroutine get_boundary_lerp_initpoint
2406 
2407  !-----------------------------------------------------------------------------
2409  subroutine get_boundary_lerp_midpoint( &
2410  bnd_DENS, &
2411  bnd_VELZ, &
2412  bnd_VELX, &
2413  bnd_VELY, &
2414  bnd_POTT, &
2415  bnd_QTRC, &
2416  now_step, &
2417  update_step )
2418  use scale_time, only: &
2419  time_dtsec
2420  implicit none
2421 
2422  ! parameters
2423  real(RP) :: EPS = 1.0e-4_rp
2424 
2425  ! arguments
2426  real(RP), intent(out) :: bnd_DENS(:,:,:)
2427  real(RP), intent(out) :: bnd_VELZ(:,:,:)
2428  real(RP), intent(out) :: bnd_VELX(:,:,:)
2429  real(RP), intent(out) :: bnd_VELY(:,:,:)
2430  real(RP), intent(out) :: bnd_POTT(:,:,:)
2431  real(RP), intent(out) :: bnd_QTRC(:,:,:,:)
2432  integer, intent(in) :: now_step
2433  integer, intent(in) :: update_step
2434 
2435  ! works
2436  integer :: i, j, k, iq
2437 
2438  real(RP) :: real_nstep
2439  real(RP) :: half_nstep
2440  real(RP) :: t1, t2
2441  !---------------------------------------------------------------------------
2442 
2443  real_nstep = real( now_step, kind=rp )
2444  half_nstep = real( UPDATE_NSTEP, kind=RP ) * 0.5_rp
2445 
2446  ! this step before half of the parent step
2447  if( ( real_nstep - eps ) < half_nstep ) then
2448 
2449  t1 = time_dtsec / atmos_boundary_update_dt * ( real_nstep + half_nstep - 0.5_rp )
2450 
2451  do j = 1, ja
2452  do i = 1, ia
2453  do k = 1, ka
2454  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_now) * t1 &
2455  - atmos_boundary_ref_dens(k,i,j,ref_old) * ( t1 - 1.0_rp )
2456  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_now) * t1 &
2457  - atmos_boundary_ref_velz(k,i,j,ref_old) * ( t1 - 1.0_rp )
2458  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_now) * t1 &
2459  - atmos_boundary_ref_velx(k,i,j,ref_old) * ( t1 - 1.0_rp )
2460  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_now) * t1 &
2461  - atmos_boundary_ref_vely(k,i,j,ref_old) * ( t1 - 1.0_rp )
2462  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_now) * t1 &
2463  - atmos_boundary_ref_pott(k,i,j,ref_old) * ( t1 - 1.0_rp )
2464  do iq = 1, bnd_qa
2465  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * t1 &
2466  - atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) * ( t1 - 1.0_rp )
2467  end do
2468  end do
2469  end do
2470  end do
2471 
2472  ! this step after half of the parent step
2473  else if( ( real_nstep - 1.0_rp + eps ) > half_nstep ) then
2474 
2475  t1 = time_dtsec / atmos_boundary_update_dt * ( real_nstep - half_nstep - 0.5_rp )
2476 
2477  do j = 1, ja
2478  do i = 1, ia
2479  do k = 1, ka
2480  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_new) * t1 &
2481  - atmos_boundary_ref_dens(k,i,j,ref_now) * ( t1 - 1.0_rp )
2482  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_new) * t1 &
2483  - atmos_boundary_ref_velz(k,i,j,ref_now) * ( t1 - 1.0_rp )
2484  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_new) * t1 &
2485  - atmos_boundary_ref_velx(k,i,j,ref_now) * ( t1 - 1.0_rp )
2486  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_new) * t1 &
2487  - atmos_boundary_ref_vely(k,i,j,ref_now) * ( t1 - 1.0_rp )
2488  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_new) * t1 &
2489  - atmos_boundary_ref_pott(k,i,j,ref_now) * ( t1 - 1.0_rp )
2490  do iq = 1, bnd_qa
2491  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_new) * t1 &
2492  - atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * ( t1 - 1.0_rp )
2493  end do
2494  end do
2495  end do
2496  end do
2497 
2498  ! this step across half of the parent step
2499  else
2500 
2501  t1 = time_dtsec / atmos_boundary_update_dt * ( real_nstep + half_nstep - 1.0_rp )
2502  t2 = time_dtsec / atmos_boundary_update_dt * ( real_nstep - half_nstep )
2503 
2504  do j = 1, ja
2505  do i = 1, ia
2506  do k = 1, ka
2507  bnd_dens(k,i,j) = atmos_boundary_ref_dens(k,i,j,ref_new) * t2 * 0.25_rp &
2508  + atmos_boundary_ref_dens(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2509  - atmos_boundary_ref_dens(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2510  bnd_velz(k,i,j) = atmos_boundary_ref_velz(k,i,j,ref_new) * t2 * 0.25_rp &
2511  + atmos_boundary_ref_velz(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2512  - atmos_boundary_ref_velz(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2513  bnd_velx(k,i,j) = atmos_boundary_ref_velx(k,i,j,ref_new) * t2 * 0.25_rp &
2514  + atmos_boundary_ref_velx(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2515  - atmos_boundary_ref_velx(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2516  bnd_vely(k,i,j) = atmos_boundary_ref_vely(k,i,j,ref_new) * t2 * 0.25_rp &
2517  + atmos_boundary_ref_vely(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2518  - atmos_boundary_ref_vely(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2519  bnd_pott(k,i,j) = atmos_boundary_ref_pott(k,i,j,ref_new) * t2 * 0.25_rp &
2520  + atmos_boundary_ref_pott(k,i,j,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2521  - atmos_boundary_ref_pott(k,i,j,ref_old) * ( t1 - 1.0_rp ) * 0.25_rp
2522  do iq = 1, bnd_qa
2523  bnd_qtrc(k,i,j,iq) = atmos_boundary_ref_qtrc(k,i,j,iq,ref_new) * t2 * 0.25_rp &
2524  + atmos_boundary_ref_qtrc(k,i,j,iq,ref_now) * ( t1 - t2 + 3.0_rp ) * 0.25_rp &
2525  - atmos_boundary_ref_qtrc(k,i,j,iq,ref_old) * ( t1 - 1.0_rp )
2526  end do
2527  end do
2528  end do
2529  end do
2530 
2531  end if
2532 
2533  return
2534  end subroutine get_boundary_lerp_midpoint
2535 
2536  !-----------------------------------------------------------------------------
2538  subroutine update_ref_index
2539  implicit none
2540 
2541  ! works
2542  integer :: ref_tmp
2543  !---------------------------------------------------------------------------
2544 
2545  ref_tmp = ref_old
2546  ref_old = ref_now
2547  ref_now = ref_new
2548  ref_new = ref_tmp
2549 
2550  return
2551  end subroutine update_ref_index
2552 
2553  subroutine history_bnd( &
2554  ATMOS_BOUNDARY_DENS, &
2555  ATMOS_BOUNDARY_VELZ, &
2556  ATMOS_BOUNDARY_VELX, &
2557  ATMOS_BOUNDARY_VELY, &
2558  ATMOS_BOUNDARY_POTT, &
2559  ATMOS_BOUNDARY_QTRC )
2560  use scale_file_history, only: &
2561  file_history_in
2562  use mod_atmos_phy_mp_vars, only: &
2563  qs_mp
2564  implicit none
2565  real(RP), intent(in) :: ATMOS_BOUNDARY_DENS(ka,ia,ja)
2566  real(RP), intent(in) :: ATMOS_BOUNDARY_VELZ(ka,ia,ja)
2567  real(RP), intent(in) :: ATMOS_BOUNDARY_VELX(ka,ia,ja)
2568  real(RP), intent(in) :: ATMOS_BOUNDARY_VELY(ka,ia,ja)
2569  real(RP), intent(in) :: ATMOS_BOUNDARY_POTT(ka,ia,ja)
2570  real(RP), intent(in) :: ATMOS_BOUNDARY_QTRC(ka,ia,ja,bnd_qa)
2571 
2572  integer :: iq, iqa
2573 
2574  call file_history_in( atmos_boundary_dens(:,:,:), 'DENS_BND', 'Boundary Density', 'kg/m3' )
2575  call file_history_in( atmos_boundary_velz(:,:,:), 'VELZ_BND', 'Boundary velocity z-direction', 'm/s', dim_type='ZHXY' )
2576  call file_history_in( atmos_boundary_velx(:,:,:), 'VELX_BND', 'Boundary velocity x-direction', 'm/s', dim_type='ZXHY' )
2577  call file_history_in( atmos_boundary_vely(:,:,:), 'VELY_BND', 'Boundary velocity y-direction', 'm/s', dim_type='ZXYH' )
2578  call file_history_in( atmos_boundary_pott(:,:,:), 'POTT_BND', 'Boundary potential temperature', 'K' )
2579  do iq = 1, bnd_qa
2580  iqa = qs_mp + iq - 1
2581  call file_history_in( atmos_boundary_qtrc(:,:,:,iq), trim(tracer_name(iqa))//'_BND', &
2582  trim(tracer_name(iqa))//' in boundary', 'kg/kg' )
2583  enddo
2584 
2585  return
2586  end subroutine history_bnd
2587 
2588 end module mod_atmos_bnd_driver
subroutine, public comm_cartesc_nest_nestdown(HANDLE, BND_QA, DENS_send, MOMZ_send, MOMX_send, MOMY_send, RHOT_send, QTRC_send, DENS_recv, VELZ_recv, VELX_recv, VELY_recv, POTT_recv, QTRC_recv)
Boundary data transfer from parent to daughter: nestdown.
real(rp), dimension(:,:,:), allocatable, target, public momz
integer, dimension(2), public daughter_ia
daughter max number in x-direction (with halo)
real(dp) function, public calendar_combine_daysec(absday, abssec)
Combine day and second.
integer, dimension(2), public daughter_ja
daughter max number in y-direction (with halo)
subroutine, public atmos_boundary_driver_update
Update boundary value with a constant time boundary.
logical, public online_iam_parent
a flag to say "I am a parent"
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_vely
subroutine, public calendar_date2char(chardate, ymdhms, subsec)
Convert from gregorian date to absolute day/second.
subroutine, public atmos_boundary_driver_set
set
integer, public time_nstep
total steps [number]
Definition: scale_time.F90:75
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_dens
integer, dimension(2), public daughter_ka
daughter max number in z-direction (with halo)
integer, public ia
of whole cells: x, local, with HALO
module Atmosphere / Physics Cloud Microphysics
logical, public atmos_boundary_update_flag
subroutine, public atmos_boundary_driver_setup
Setup.
module atmosphere / reference state
integer, dimension(2), public parent_nstep
parent step [number]
module ATMOSPHERIC Variables
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfz
center buffer factor (0-1): z
real(rp), dimension(:,:,:), allocatable, target, public momx
integer, public qa
subroutine, public atmos_boundary_driver_finalize
Finalize boundary value.
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_pott
integer, public ja
of whole cells: y, local, with HALO
module process / cartesC
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
logical, public prc_has_s
subroutine update_ref_index
Update indices of array of boundary references.
character(len=h_short), dimension(qa_max), public tracer_name
logical, public prc_has_n
logical, public prc_has_e
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_velz
real(rp), public const_undef
Definition: scale_const.F90:41
character(len=h_short), dimension(qa_max), public tracer_unit
module COMMUNICATION
integer, public comm_cartesc_nest_bnd_qa
number of tracer treated in nesting system
real(rp), dimension(0-1), allocatable, public dens
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, dimension(2), public parent_ia
parent max number in x-direction (with halo)
module TRACER
module Index
Definition: scale_index.F90:11
subroutine, public comm_cartesc_nest_recv_cancel(HANDLE)
Sub-command for data transfer from parent to daughter: nestdown.
real(rp), public atmos_boundary_smoother_fact
module atmosphere / hydrometeor
subroutine, public atmos_phy_mp_driver_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:38
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
logical, public online_iam_daughter
a flag to say "I am a daughter"
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_dens
refernce density [kg/m3]
subroutine, public comm_cartesc_nest_recvwait_issue(HANDLE, BND_QA)
Sub-command for data transfer from parent to daughter: nestdown.
module ATMOSPHERE / Boundary treatment
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
real(rp), dimension(:,:,:,:), allocatable, public atmos_boundary_alpha_qtrc
module TIME
Definition: scale_time.F90:16
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_velx
module atmosphere / grid / cartesC
integer, public ks
start point of inner domain: z, local
integer, dimension(2), public parent_ja
parent max number in y-direction (with halo)
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfx
center buffer factor (0-1): x
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
subroutine, public comm_cartesc_nest_test(HANDLE)
[check communication status] Inter-communication
module CONSTANT
Definition: scale_const.F90:11
real(rp), dimension(0-1), allocatable, public qtrc
module Communication CartesianC nesting
integer, public js
start point of inner domain: y, local
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfx
face buffer factor (0-1): x
real(rp), dimension(:,:,:), allocatable, target, public momy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfy
center buffer factor (0-1): y
module Atmosphere / Physics Chemistry
subroutine atmos_boundary_set_file
Set boundary value for real case experiment.
module profiler
Definition: scale_prof.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:33
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfz
face buffer factor (0-1): z
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfy
face buffer factor (0-1): y
module PRECISION
real(dp), dimension(2), public parent_dtsec
parent DT [sec]
module file / cartesianC
integer, dimension(2), public parent_ka
parent max number in z-direction (with halo)
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_pi
pi
Definition: scale_const.F90:31
module CALENDAR
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:69
real(rp), dimension(:,:,:,:), allocatable, target, public qe
module STDIO
Definition: scale_io.F90:10
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_velz
subroutine, public file_cartesc_flush(fid)
Flush all pending requests to a netCDF file (PnetCDF only)
integer, parameter, public n_hyd
module atmosphere / physics / cloud microphysics
integer, parameter, public rp
subroutine, public file_cartesc_open(basename, fid, aggregate)
open a netCDF file for read
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_pott
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_dens
subroutine, public calendar_date2daysec(absday, abssec, ymdhms, subsec, offset_year)
Convert from gregorian date to absolute day/second.
subroutine, public file_cartesc_close(fid)
Close a netCDF file.
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_vely
real(rp), dimension(:,:,:,:), allocatable, public atmos_boundary_qtrc
subroutine atmos_boundary_set_online
Set boundary value for real case experiment [online daughter].
module file_history
logical, public prc_has_w
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_velx