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