SCALE-RM
mod_realinput.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_stdio
23  use scale_grid_real, only: &
24  lon => real_lon, &
25  lat => real_lat, &
26  lonx => real_lonx, &
27  laty => real_laty, &
28  cz => real_cz, &
29  fz => real_fz
30  use scale_grid_nest, only: &
32  use scale_index
33  use scale_tracer
34  use scale_process, only: &
35  prc_ismaster, &
37  use scale_external_io, only: &
38  iscale, &
39  iwrfarw, &
40  inicam, &
41  igrads
42  use scale_comm, only: &
43  comm_bcast
44  !-----------------------------------------------------------------------------
45  implicit none
46  private
47  !-----------------------------------------------------------------------------
48  !
49  !++ Public procedure
50  !
51  public :: realinput_atmos
52  public :: realinput_surface
53 
54  !-----------------------------------------------------------------------------
55  !
56  !++ Public parameters & variables
57  !
58  !-----------------------------------------------------------------------------
59  !
60  !++ Private procedure
61  !
62  private :: parentatominput
63  private :: parentatomboundary
64  private :: parentsurfacesetup
65  private :: parentsurfaceinput
66  private :: parentoceanboundary
67  private :: interp_oceanland_data
68 
69  !-----------------------------------------------------------------------------
70  !
71  !++ Private parameters & variables
72  !
73  integer, parameter :: I_intrp_off = 0
74  integer, parameter :: I_intrp_mask = 1
75  integer, parameter :: I_intrp_fill = 2
76 
77  integer, parameter :: cosin = 1
78  integer, parameter :: sine = 2
79 
80  real(RP), allocatable :: lon_org (:,:)
81  real(RP), allocatable :: lat_org (:,:)
82  real(RP), allocatable :: cz_org(:,:,:)
83 
84  real(RP), allocatable :: dens_org(:,:,:)
85  real(RP), allocatable :: qtrc_org(:,:,:,:)
86 
87  real(RP), allocatable :: velz_org(:,:,:)
88  real(RP), allocatable :: velx_org(:,:,:)
89  real(RP), allocatable :: vely_org(:,:,:)
90  real(RP), allocatable :: pott_org(:,:,:)
91  real(RP), allocatable :: temp_org(:,:,:)
92  real(RP), allocatable :: pres_org(:,:,:)
93 
94  real(RP), allocatable :: hfact(:,:,:)
95  real(RP), allocatable :: vfact(:,:,:,:,:)
96  integer, allocatable :: igrd (:,:,:)
97  integer, allocatable :: jgrd (:,:,:)
98  integer, allocatable :: kgrd (:,:,:,:,:)
99  integer, allocatable :: ncopy(:,:,:)
100 
101  real(RP), allocatable :: tw_org(:,:)
102  real(RP), allocatable :: sst_org(:,:)
103  real(RP), allocatable :: albw_org(:,:,:)
104  real(RP), allocatable :: olon_org(:,:)
105  real(RP), allocatable :: olat_org(:,:)
106  real(RP), allocatable :: omask_org(:,:)
107 
108  integer, private :: itp_nh = 4
109  integer, private :: itp_nv = 2
110 
111  integer, private :: io_fid_grads_nml = -1
112  integer, private :: io_fid_grads_data = -1
113 
114  logical, private :: do_read_atom
115  logical, private :: do_read_land
116  logical, private :: do_read_ocean
117  logical, private :: rotate
118  logical, private :: use_waterratio
119  logical, private :: update_coord
120  logical, private :: use_temp
121  logical, private :: serial
122  logical, private :: serial_land
123  logical, private :: serial_ocean
124  logical, private :: first = .true.
125 
126  integer, private :: i_intrp_land_temp
127  integer, private :: i_intrp_land_water
128  integer, private :: i_intrp_land_sfc_temp
129  integer, private :: i_intrp_ocean_temp
130  integer, private :: i_intrp_ocean_sfc_temp
131 
132 
133  ! replace missing value
134  real(RP), private, parameter :: maskval_tg = 298.0_rp
135  real(RP), private, parameter :: maskval_strg = 0.02_rp
136  ! default value 0.02: set as value of forest at 40% of evapolation rate.
137  ! forest is considered as a typical landuse over Japan area.
138 
139 
140  ! for namelist
141  integer :: NUMBER_OF_FILES = 1
142  integer :: NUMBER_OF_TSTEPS = 1 ! num of time steps in one file
143  integer :: NUMBER_OF_SKIP_TSTEPS = 0 ! num of skipped first several data
144 
145  character(len=H_LONG) :: FILETYPE_ORG = ''
146  character(len=H_LONG) :: BASENAME_ORG = ''
147  logical :: BASENAME_ADD_NUM = .false.
148  character(len=H_LONG) :: BASENAME_BOUNDARY = 'boundary_atmos'
149  character(len=H_LONG) :: BOUNDARY_TITLE = 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
150  real(RP) :: BOUNDARY_UPDATE_DT = 0.0_rp ! inteval time of boudary data update [s]
151 
152  integer :: PARENT_MP_TYPE = 6 ! microphysics type of the parent model (number of classes)
153  ! 0: dry, 3:3-class, 5:5-class, 6:6-class, >6:double moment
154  logical :: SERIAL_PROC_READ = .true. ! read by one MPI process and broadcast
155  ! only for SCALE boundary
156  logical :: USE_FILE_DENSITY = .false. ! use density data from files
157  !-----------------------------------------------------------------------------
158 contains
159  !-----------------------------------------------------------------------------
160  subroutine realinput_atmos( &
161  flg_intrp )
162  use mod_atmos_vars, only: &
163  dens, &
164  momz, &
165  momx, &
166  momy, &
167  rhot, &
168  qtrc
169  use mod_atmos_admin, only: &
171  implicit none
172 
173  logical, intent(in) :: flg_intrp ! flag for interpolation of SBM(S10) from outer bulk-MP model
174 
175 
176  namelist / param_mkinit_real_atmos / &
177  number_of_files, &
178  number_of_tsteps, &
179  number_of_skip_tsteps, &
180  filetype_org, &
181  basename_org, &
182  basename_add_num, &
183  basename_boundary, &
184  boundary_title, &
185  boundary_update_dt, &
186  parent_mp_type, &
187  serial_proc_read, &
188  use_file_density
189 
190  character(len=H_LONG) :: basename_atmos = ''
191  character(len=H_LONG) :: basename_withnum = ''
192  character(len=5) :: num = ''
193 
194  ! atmos
195  real(RP), allocatable :: dens_org(:,:,:,:)
196  real(RP), allocatable :: momz_org(:,:,:,:)
197  real(RP), allocatable :: momx_org(:,:,:,:)
198  real(RP), allocatable :: momy_org(:,:,:,:)
199  real(RP), allocatable :: rhot_org(:,:,:,:)
200  real(RP), allocatable :: qtrc_org(:,:,:,:,:)
201 
202  integer :: mdlid
203  integer :: dims(6) ! dims 1-3: normal, 4-6: staggerd
204 
205  integer :: totaltimesteps = 1
206  integer :: timelen
207  integer :: skip_steps
208  integer :: ierr
209  logical :: flg_bin ! flag for SBM(S10) is used or not 0-> not used, 1-> used
210 
211  integer :: k, i, j, iq, n, ns, ne
212  !---------------------------------------------------------------------------
213 
214  if( io_l ) write(io_fid_log,*)
215  if( io_l ) write(io_fid_log,*) '+++ Module[RealCaseAtmos]/Categ[MKINIT]'
216 
217  !--- read namelist
218  rewind(io_fid_conf)
219  read(io_fid_conf,nml=param_mkinit_real_atmos,iostat=ierr)
220  if( ierr < 0 ) then !--- missing
221  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
222  elseif( ierr > 0 ) then !--- fatal error
223  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_ATMOS. Check!'
224  call prc_mpistop
225  endif
226  if( io_nml ) write(io_fid_nml,nml=param_mkinit_real_atmos)
227 
228  if( atmos_phy_mp_type == 'SUZUKI10' ) then
229  flg_bin = .true.
230  else
231  flg_bin = .false.
232  endif
233 
234  if ( filetype_org == "GrADS" ) then
235  basename_withnum = basename_org ! namelist file name
236  basename_atmos = ""
237  else
238  if ( number_of_files > 1 .or. basename_add_num ) then
239  basename_withnum = trim(basename_org)//"_00000"
240  else
241  basename_withnum = trim(basename_org)
242  end if
243  basename_atmos = basename_org
244  end if
245 
246  call parentatomsetup( dims(:), timelen, mdlid, & ![OUT]
247  basename_withnum, & ![IN]
248  filetype_org, & ![IN]
249  use_file_density, & ![IN]
250  serial_proc_read ) ![IN]
251 
252  if ( boundary_update_dt <= 0.0_rp ) then
253  write(*,*) 'xxx BOUNDARY_UPDATE_DT is necessary in real case preprocess'
254  call prc_mpistop
255  endif
256 
257  if ( timelen > 0 ) then
258  number_of_tsteps = timelen ! read from file
259  endif
260 
261  totaltimesteps = number_of_files * number_of_tsteps
262 
263  allocate( dens_org(ka,ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
264  allocate( momz_org(ka,ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
265  allocate( momx_org(ka,ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
266  allocate( momy_org(ka,ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
267  allocate( rhot_org(ka,ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
268  allocate( qtrc_org(ka,ia,ja,qa,1+number_of_skip_tsteps:totaltimesteps) )
269 
270  !--- read external file
271  do n = 1, number_of_files
272 
273  if ( number_of_files > 1 .or. basename_add_num ) then
274  write(num,'(I5.5)') n-1
275  basename_withnum = trim(basename_atmos)//"_"//num
276  else
277  basename_withnum = trim(basename_atmos)
278  end if
279 
280  if( io_l ) write(io_fid_log,*) ' '
281  if( io_l ) write(io_fid_log,*) '+++ Target File Name: ',trim(basename_withnum)
282  if( io_l ) write(io_fid_log,*) ' Time Steps in One File: ', number_of_tsteps
283 
284  ns = number_of_tsteps * (n - 1) + 1
285  ne = ns + (number_of_tsteps - 1)
286 
287  if ( ne <= number_of_skip_tsteps ) then
288  if( io_l ) write(io_fid_log,*) ' SKIP'
289  cycle
290  end if
291 
292  skip_steps = max(number_of_skip_tsteps - ns + 1, 0)
293  ns = max(ns, number_of_skip_tsteps+1)
294 
295  ! read all prepared data
296  call parentatominput( dens_org(:,:,:,ns:ne), &
297  momz_org(:,:,:,ns:ne), &
298  momx_org(:,:,:,ns:ne), &
299  momy_org(:,:,:,ns:ne), &
300  rhot_org(:,:,:,ns:ne), &
301  qtrc_org(:,:,:,:,ns:ne), &
302  basename_withnum, &
303  dims(:), &
304  mdlid, &
305  flg_bin, flg_intrp, &
306  parent_mp_type, &
307  number_of_tsteps, &
308  skip_steps )
309  enddo
310 
311  !--- input initial data
312  ns = number_of_skip_tsteps + 1 ! skip first several data
313  do j = 1, ja
314  do i = 1, ia
315  do k = ks, ke
316  dens(k,i,j) = dens_org(k,i,j,ns)
317  momz(k,i,j) = momz_org(k,i,j,ns)
318  momx(k,i,j) = momx_org(k,i,j,ns)
319  momy(k,i,j) = momy_org(k,i,j,ns)
320  rhot(k,i,j) = rhot_org(k,i,j,ns)
321 
322  do iq = 1, qa
323  qtrc(k,i,j,iq) = qtrc_org(k,i,j,iq,ns)
324  enddo
325  enddo
326  enddo
327  enddo
328 
329  !--- output boundary data
330  totaltimesteps = totaltimesteps - number_of_skip_tsteps ! skip first several data
331  call parentatomboundary( dens_org(:,:,:,ns:ne), &
332  momz_org(:,:,:,ns:ne), &
333  momx_org(:,:,:,ns:ne), &
334  momy_org(:,:,:,ns:ne), &
335  rhot_org(:,:,:,ns:ne), &
336  qtrc_org(:,:,:,:,ns:ne), &
337  totaltimesteps, &
338  boundary_update_dt, &
339  basename_boundary, &
340  boundary_title )
341 
342  deallocate( dens_org )
343  deallocate( momz_org )
344  deallocate( momx_org )
345  deallocate( momy_org )
346  deallocate( rhot_org )
347  deallocate( qtrc_org )
348 
349  return
350  end subroutine realinput_atmos
351 
352  !-----------------------------------------------------------------------------
353  subroutine realinput_surface
354  use scale_const, only: &
355  i_sw => const_i_sw, &
356  i_lw => const_i_lw
357  use mod_atmos_vars, only: &
358  dens, &
359  momz, &
360  momx, &
361  momy, &
362  rhot, &
363  qtrc
364  use mod_land_vars, only: &
365  land_temp, &
366  land_water, &
367  land_sfc_temp, &
369  use mod_urban_vars, only: &
370  urban_sfc_temp, &
372  urban_tc, &
373  urban_qc, &
374  urban_uc, &
375  urban_tr, &
376  urban_tb, &
377  urban_tg, &
378  urban_trl, &
379  urban_tbl, &
380  urban_tgl, &
381  urban_rainr, &
382  urban_rainb, &
383  urban_raing, &
384  urban_roff
385  use mod_ocean_vars, only: &
386  ocean_temp, &
387  ocean_sfc_temp, &
389  ocean_sfc_z0m, &
390  ocean_sfc_z0h, &
392  use mod_atmos_phy_sf_vars, only: &
398  use scale_landuse, only: &
399  fact_ocean => landuse_fact_ocean, &
400  fact_land => landuse_fact_land, &
401  fact_urban => landuse_fact_urban
402  implicit none
403 
404  logical :: use_file_landwater = .true. ! use land water data from files
405  real(RP) :: init_landwater_ratio = 0.5_rp ! Ratio of land water to storage is constant, if USE_FILE_LANDWATER is ".false."
406  real(RP) :: init_ocean_alb_lw = 0.04_rp ! initial LW albedo on the ocean
407  real(RP) :: init_ocean_alb_sw = 0.10_rp ! initial SW albedo on the ocean
408  real(RP) :: init_ocean_z0w = 1.0e-3_rp ! initial surface roughness on the ocean
409  character(len=H_SHORT) :: intrp_land_temp = 'off'
410  character(len=H_SHORT) :: intrp_land_water = 'off'
411  character(len=H_SHORT) :: intrp_land_sfc_temp = 'off'
412  character(len=H_SHORT) :: intrp_ocean_temp = 'off'
413  character(len=H_SHORT) :: intrp_ocean_sfc_temp = 'off'
414  integer :: intrp_iter_max = 100
415  character(len=H_SHORT) :: soilwater_ds2vc = 'limit'
416  logical :: soilwater_ds2vc_flag ! true: 'critical', false: 'limit'
417  logical :: elevation_collection = .true.
418 
419  namelist / param_mkinit_real_land / &
420  number_of_files, &
421  number_of_tsteps, &
422  number_of_skip_tsteps, &
423  filetype_org, &
424  basename_org, &
425  basename_add_num, &
426  use_file_landwater, &
427  init_landwater_ratio, &
428  intrp_land_temp, &
429  intrp_land_water, &
430  intrp_land_sfc_temp, &
431  intrp_iter_max, &
432  soilwater_ds2vc, &
433  elevation_collection, &
434  serial_proc_read
435 
436  namelist / param_mkinit_real_ocean / &
437  number_of_files, &
438  number_of_tsteps, &
439  number_of_skip_tsteps, &
440  filetype_org, &
441  basename_org, &
442  basename_add_num, &
443  basename_boundary, &
444  boundary_title, &
445  boundary_update_dt, &
446  init_ocean_alb_lw, &
447  init_ocean_alb_sw, &
448  init_ocean_z0w, &
449  intrp_ocean_temp, &
450  intrp_ocean_sfc_temp, &
451  intrp_iter_max, &
452  serial_proc_read
453 
454  character(len=H_LONG) :: filetype_land
455  character(len=H_LONG) :: filetype_ocean
456  character(len=H_LONG) :: basename_land
457  character(len=H_LONG) :: basename_ocean
458  character(len=H_LONG) :: basename_withnum = ''
459  character(len=5) :: num = ''
460  logical :: serial_proc_read_land
461  logical :: serial_proc_read_ocean
462 
463  ! land
464  real(RP) :: land_temp_org(lkmax,ia,ja)
465  real(RP) :: land_water_org(lkmax,ia,ja)
466  real(RP) :: land_sfc_temp_org(ia,ja)
467  real(RP) :: land_sfc_albedo_org(ia,ja,2)
468 
469  ! urban
470  real(RP) :: urban_tc_org(ia,ja)
471  real(RP) :: urban_qc_org(ia,ja)
472  real(RP) :: urban_uc_org(ia,ja)
473  real(RP) :: urban_sfc_temp_org(ia,ja)
474  real(RP) :: urban_sfc_albedo_org(ia,ja,2)
475 
476  ! ocean
477  real(RP), allocatable :: ocean_temp_org(:,:,:)
478  real(RP), allocatable :: ocean_sfc_temp_org(:,:,:)
479  real(RP), allocatable :: ocean_sfc_albedo_org(:,:,:,:)
480  real(RP), allocatable :: ocean_sfc_z0_org(:,:,:)
481 
482  integer :: mdlid_land, mdlid_ocean
483  integer :: ldims(3), odims(2)
484 
485  integer :: totaltimesteps = 1
486  integer :: timelen
487  integer :: skip_steps
488  integer :: lit
489  integer :: lfn
490  integer :: ierr
491 
492  integer :: k, i, j, n, ns, ne
493  !---------------------------------------------------------------------------
494 
495  if( io_l ) write(io_fid_log,*)
496  if( io_l ) write(io_fid_log,*) '+++ Module[RealCaseSurface]/Categ[MKINIT]'
497 
498 
499  ! LAND/URBAN
500 
501  !--- read namelist
502  rewind(io_fid_conf)
503  read(io_fid_conf,nml=param_mkinit_real_land,iostat=ierr)
504  if( ierr < 0 ) then !--- missing
505  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
506  elseif( ierr > 0 ) then !--- fatal error
507  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_LAND. Check!'
508  call prc_mpistop
509  endif
510  if( io_nml ) write(io_fid_nml,nml=param_mkinit_real_land)
511 
512  filetype_land = filetype_org
513 
514  lfn = number_of_skip_tsteps / number_of_tsteps
515  if ( filetype_land .ne. "GrADS" .and. ( number_of_files > 1 .or. basename_add_num ) ) then
516  write(num,'(I5.5)') lfn
517  basename_land = trim(basename_org)//"_"//num
518  else
519  basename_land = trim(basename_org)
520  end if
521 
522  serial_land = serial_proc_read
523 
524  lit = mod(number_of_skip_tsteps,number_of_tsteps)+1
525 
526  !--- read external file
527  if( io_l ) write(io_fid_log,*) ' '
528  if( io_l ) write(io_fid_log,*) '+++ Target File Name (Land): ',trim(basename_land)
529  if( io_l ) write(io_fid_log,*) ' Time Steps in One File: ', number_of_tsteps
530  if( io_l ) write(io_fid_log,*) ' Time Step to read: ', lit
531 
532 
533 
534  ! OCEAN
535 
536  !--- read namelist
537  rewind(io_fid_conf)
538  read(io_fid_conf,nml=param_mkinit_real_ocean,iostat=ierr)
539  if( ierr < 0 ) then !--- missing
540  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
541  elseif( ierr > 0 ) then !--- fatal error
542  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_OCEAN. Check!'
543  call prc_mpistop
544  endif
545  if( io_nml ) write(io_fid_nml,nml=param_mkinit_real_ocean)
546 
547  filetype_ocean = filetype_org
548 
549  if ( filetype_ocean .ne. "GrADS" .and. ( number_of_files > 1 .or. basename_add_num ) ) then
550  basename_ocean = trim(basename_org)//"_00000"
551  else
552  basename_ocean = trim(basename_org)
553  end if
554 
555  select case( soilwater_ds2vc )
556  case( 'critical' )
557  soilwater_ds2vc_flag = .true.
558  case('limit' )
559  soilwater_ds2vc_flag = .false.
560  case default
561  write(*,*) 'xxx Unsupported SOILWATER_DS2CV TYPE:', trim(soilwater_ds2vc)
562  call prc_mpistop
563  end select
564 
565  serial_ocean = serial_proc_read
566 
567 
568  call parentsurfacesetup( ldims, odims, & ![OUT]
569  mdlid_land, & ![OUT]
570  mdlid_ocean, & ![OUT]
571  timelen, & ![OUT]
572  basename_land, & ![IN]
573  basename_ocean, & ![IN]
574  filetype_land, & ![IN]
575  filetype_ocean, & ![IN]
576  use_file_landwater, & ![IN]
577  intrp_land_temp, & ![IN]
578  intrp_land_water, & ![IN]
579  intrp_land_sfc_temp, & ![IN]
580  intrp_ocean_temp, & ![IN]
581  intrp_ocean_sfc_temp ) ![IN]
582 
583  if ( timelen > 0 ) then
584  number_of_tsteps = timelen ! read from file
585  endif
586 
587  totaltimesteps = number_of_files * number_of_tsteps
588 
589  allocate( ocean_temp_org(ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
590  allocate( ocean_sfc_temp_org(ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
591  allocate( ocean_sfc_albedo_org(ia,ja,2,1+number_of_skip_tsteps:totaltimesteps ) )
592  allocate( ocean_sfc_z0_org(ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
593 
594  if ( mdlid_land == igrads ) then
595  if ( number_of_files > 1 .or. basename_add_num ) then
596  write(num,'(I5.5)') lfn
597  basename_land = "_"//num
598  else
599  basename_land = ""
600  end if
601  end if
602 
603  if ( mdlid_ocean == igrads ) then
604  basename_org = ""
605  end if
606 
607  !--- read external file
608  do n = 1, number_of_files
609 
610  if ( number_of_files > 1 .or. basename_add_num ) then
611  write(num,'(I5.5)') n-1
612  basename_ocean = trim(basename_org)//"_"//num
613  else
614  basename_ocean = trim(basename_org)
615  end if
616 
617  if( io_l ) write(io_fid_log,*) ' '
618  if( io_l ) write(io_fid_log,*) '+++ Target File Name (Ocean): ', trim(basename_ocean)
619  if( io_l ) write(io_fid_log,*) ' Time Steps in One File: ', number_of_tsteps
620 
621  ns = number_of_tsteps * (n - 1) + 1
622  ne = ns + (number_of_tsteps - 1)
623 
624  if ( ne <= number_of_skip_tsteps ) then
625  if( io_l ) write(io_fid_log,*) ' SKIP'
626  cycle
627  end if
628 
629  skip_steps = max(number_of_skip_tsteps - ns + 1, 0)
630  ns = max(ns, number_of_skip_tsteps+1)
631 
632  ! read all prepared data
633  call parentsurfaceinput( land_temp_org, &
634  land_water_org, &
635  land_sfc_temp_org, &
636  land_sfc_albedo_org, &
637  urban_tc_org, &
638  urban_qc_org, &
639  urban_uc_org, &
640  urban_sfc_temp_org, &
641  urban_sfc_albedo_org, &
642  ocean_temp_org(:,:,ns:ne), &
643  ocean_sfc_temp_org(:,:,ns:ne), &
644  ocean_sfc_albedo_org(:,:,:,ns:ne), &
645  ocean_sfc_z0_org(:,:,ns:ne), &
646  dens, &
647  momz, &
648  momx, &
649  momy, &
650  rhot, &
651  qtrc, &
652  basename_land, &
653  basename_ocean, &
654  mdlid_land, mdlid_ocean, &
655  ldims, odims, &
656  use_file_landwater, &
657  init_landwater_ratio, &
658  init_ocean_alb_lw, &
659  init_ocean_alb_sw, &
660  init_ocean_z0w, &
661  intrp_iter_max, &
662  soilwater_ds2vc_flag, &
663  elevation_collection, &
664  number_of_tsteps, &
665  skip_steps, lit )
666 
667  enddo
668 
669 
670  !--- input initial data
671  do j = 1, ja
672  do i = 1, ia
673  land_sfc_temp(i,j) = land_sfc_temp_org(i,j)
674  land_sfc_albedo(i,j,i_lw) = land_sfc_albedo_org(i,j,i_lw)
675  land_sfc_albedo(i,j,i_sw) = land_sfc_albedo_org(i,j,i_sw)
676  do k = 1, lkmax
677  land_temp(k,i,j) = land_temp_org(k,i,j)
678  land_water(k,i,j) = land_water_org(k,i,j)
679  end do
680 
681  urban_sfc_temp(i,j) = urban_sfc_temp_org(i,j)
682  urban_sfc_albedo(i,j,i_lw) = urban_sfc_albedo_org(i,j,i_lw)
683  urban_sfc_albedo(i,j,i_sw) = urban_sfc_albedo_org(i,j,i_sw)
684  do k = uks, uke
685  urban_trl(k,i,j) = urban_sfc_temp_org(i,j)
686  urban_tbl(k,i,j) = urban_sfc_temp_org(i,j)
687  urban_tgl(k,i,j) = urban_sfc_temp_org(i,j)
688  end do
689  urban_tc(i,j) = urban_tc_org(i,j)
690  urban_qc(i,j) = urban_qc_org(i,j)
691  urban_uc(i,j) = urban_uc_org(i,j)
692  urban_tr(i,j) = urban_sfc_temp_org(i,j)
693  urban_tb(i,j) = urban_sfc_temp_org(i,j)
694  urban_tg(i,j) = urban_sfc_temp_org(i,j)
695  urban_rainr(i,j) = 0.0_rp
696  urban_rainb(i,j) = 0.0_rp
697  urban_raing(i,j) = 0.0_rp
698  urban_roff(i,j) = 0.0_rp
699  enddo
700  enddo
701 
702 
703  ns = number_of_skip_tsteps + 1 ! skip first several data
704  do j = 1, ja
705  do i = 1, ia
706  ocean_temp(i,j) = ocean_temp_org(i,j,ns)
707  ocean_sfc_temp(i,j) = ocean_sfc_temp_org(i,j,ns)
708  ocean_sfc_albedo(i,j,i_lw) = ocean_sfc_albedo_org(i,j,i_lw,ns)
709  ocean_sfc_albedo(i,j,i_sw) = ocean_sfc_albedo_org(i,j,i_sw,ns)
710  ocean_sfc_z0m(i,j) = ocean_sfc_z0_org(i,j,ns)
711  ocean_sfc_z0h(i,j) = ocean_sfc_z0_org(i,j,ns)
712  ocean_sfc_z0e(i,j) = ocean_sfc_z0_org(i,j,ns)
713  enddo
714  enddo
715 
716  do j = 1, ja
717  do i = 1, ia
718  atmos_phy_sf_sfc_temp(i,j) = fact_ocean(i,j) * ocean_sfc_temp(i,j) &
719  + fact_land(i,j) * land_sfc_temp(i,j) &
720  + fact_urban(i,j) * urban_sfc_temp(i,j)
721  atmos_phy_sf_sfc_albedo(i,j,i_lw) = fact_ocean(i,j) * ocean_sfc_albedo(i,j,i_lw) &
722  + fact_land(i,j) * land_sfc_albedo(i,j,i_lw) &
723  + fact_urban(i,j) * urban_sfc_albedo(i,j,i_lw)
724  atmos_phy_sf_sfc_albedo(i,j,i_sw) = fact_ocean(i,j) * ocean_sfc_albedo(i,j,i_sw) &
725  + fact_land(i,j) * land_sfc_albedo(i,j,i_sw) &
726  + fact_urban(i,j) * urban_sfc_albedo(i,j,i_sw)
730  end do
731  end do
732 
733 
734  !--- output boundary data
735  totaltimesteps = totaltimesteps - number_of_skip_tsteps ! skip first several data
736  if ( totaltimesteps > 1 ) then
737  if ( boundary_update_dt <= 0.0_rp ) then
738  write(*,*) 'xxx BOUNDARY_UPDATE_DT is necessary in real case preprocess'
739  call prc_mpistop
740  endif
741 
742  call parentoceanboundary( ocean_temp_org(:,:,ns:ne), &
743  ocean_sfc_temp_org(:,:,ns:ne), &
744  ocean_sfc_albedo_org(:,:,:,ns:ne), &
745  ocean_sfc_z0_org(:,:,ns:ne), &
746  totaltimesteps, &
747  boundary_update_dt, &
748  basename_boundary, &
749  boundary_title )
750 
751  end if
752 
753  deallocate( ocean_temp_org )
754  deallocate( ocean_sfc_temp_org )
755  deallocate( ocean_sfc_albedo_org )
756  deallocate( ocean_sfc_z0_org )
757 
758  return
759  end subroutine realinput_surface
760 
761 
762  !-----------------------------------------------------------------------------
764  subroutine parentatomsetup( &
765  dims, &
766  timelen, &
767  mdlid, &
768  basename_org, &
769  filetype, &
770  use_file_density_in, &
771  serial_in )
772  use scale_external_io, only: &
773  iscale, &
774  iwrfarw, &
775  inicam, &
776  igrads
777  use mod_realinput_scale, only: &
779  use mod_realinput_wrfarw, only: &
781  use mod_realinput_nicam, only: &
783  use mod_realinput_grads, only: &
785  implicit none
786 
787  integer, intent(out) :: dims(6)
788  integer, intent(out) :: timelen
789  integer, intent(out) :: mdlid
790  character(len=*), intent(in) :: basename_org
791  character(len=*), intent(in) :: filetype
792  logical, intent(in) :: serial_in ! read by a serial process
793  logical, intent(in) :: use_file_density_in ! use density data from files
794 
795  !---------------------------------------------------------------------------
796 
797  if( io_l ) write(io_fid_log,*)
798  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[RealinputAtmos]/Categ[Setup]'
799 
800  serial = serial_in
801  if( serial ) then
802  if( prc_ismaster ) then
803  do_read_atom = .true.
804  else
805  do_read_atom = .false.
806  endif
807  else
808  do_read_atom = .true.
809  endif
810 
811  select case(trim(filetype))
812  case('SCALE-RM')
813 
814  mdlid = iscale
815  do_read_atom = .true.
816  serial = .false.
817  call parentatomsetupscale( dims ) ! (out)
818  update_coord = .false.
819  use_file_density = use_file_density_in
820  use_temp = .false.
821  rotate = .false.
822  timelen = -1
823 
824  case('WRF-ARW')
825 
826  mdlid = iwrfarw
827  if ( do_read_atom ) call parentatomsetupwrfarw( dims, timelen, & ! (out)
828  basename_org ) ! (in)
829  update_coord = .true.
830  use_file_density = .false.
831  use_temp = .true.
832  rotate = .true.
833 
834  case('NICAM-NETCDF')
835 
836  mdlid = inicam
837  if ( do_read_atom ) call parentatomsetupnicam( dims, timelen, & ! (out)
838  basename_org ) ! (in)
839  update_coord = .false.
840  use_file_density = .false.
841  use_temp = .true.
842  rotate = .true.
843 
844  case('GrADS')
845 
846  mdlid = igrads
847  if ( do_read_atom ) call parentatomsetupgrads( dims, & ! (out)
848  basename_org ) ! (in)
849  update_coord = .true.
850  use_file_density = use_file_density_in
851  use_temp = .true.
852  rotate = .true.
853  timelen = -1
854 
855  case default
856 
857  write(*,*) 'xxx Unsupported FILE TYPE:', trim(filetype)
858  call prc_mpistop
859 
860  endselect
861 
862  if( serial ) then
863  call comm_bcast( dims(:), 6 )
864  call comm_bcast( timelen )
865  endif
866 
867  if( io_l ) write(io_fid_log,*) '+++ Horizontal Interpolation Level:', &
869  itp_nh = int( nest_interp_level )
870  itp_nv = 2
871 
872  allocate( hfact( ia, ja, itp_nh ) )
873  allocate( vfact( ka, ia, ja, itp_nh, itp_nv ) )
874  allocate( igrd( ia, ja, itp_nh ) )
875  allocate( jgrd( ia, ja, itp_nh ) )
876  allocate( kgrd( ka, ia, ja, itp_nh, itp_nv ) )
877  allocate( ncopy( ia, ja, itp_nh ) )
878 
879  allocate( lon_org( dims(2), dims(3) ) )
880  allocate( lat_org( dims(2), dims(3) ) )
881  allocate( cz_org( dims(1)+2, dims(2), dims(3) ) )
882 
883  allocate( velz_org( dims(1)+2, dims(2), dims(3) ) )
884  allocate( velx_org( dims(1)+2, dims(2), dims(3) ) )
885  allocate( vely_org( dims(1)+2, dims(2), dims(3) ) )
886  allocate( pott_org( dims(1)+2, dims(2), dims(3) ) )
887  allocate( temp_org( dims(1)+2, dims(2), dims(3) ) )
888  allocate( pres_org( dims(1)+2, dims(2), dims(3) ) )
889  allocate( qtrc_org( dims(1)+2, dims(2), dims(3), qa ) )
890  allocate( dens_org( dims(1)+2, dims(2), dims(3) ) )
891 
892  return
893  end subroutine parentatomsetup
894 
895  !-----------------------------------------------------------------------------
897  subroutine parentatominput( &
898  dens, &
899  momz, &
900  momx, &
901  momy, &
902  rhot, &
903  qtrc, &
904  basename_org, &
905  dims, &
906  mdlid, &
907  flg_bin, &
908  flg_intrp, &
909  mptype_parent, &
910  timelen, &
911  skiplen )
912  use scale_comm, only: &
913  comm_vars8, &
914  comm_wait
915  use scale_gridtrans, only: &
916  rotc => gtrans_rotc
917  use scale_atmos_thermodyn, only: &
918  thermodyn_pott => atmos_thermodyn_pott
919  use scale_atmos_hydrometeor, only: &
920  hydrometeor_diagnose_number_concentration => atmos_hydrometeor_diagnose_number_concentration, &
921  i_qv, &
922  i_qc, &
923  qls, &
924  qle
925  use scale_atmos_hydrostatic, only: &
926  hydrostatic_buildrho_real => atmos_hydrostatic_buildrho_real
927  use scale_interpolation_nest, only: &
931  use mod_atmos_admin, only: &
933  use mod_realinput_scale, only: &
936  use mod_realinput_wrfarw, only: &
939  use mod_realinput_nicam, only: &
942  use mod_realinput_grads, only: &
945  implicit none
946 
947  real(RP), intent(out) :: dens(:,:,:,:)
948  real(RP), intent(out) :: momz(:,:,:,:)
949  real(RP), intent(out) :: momx(:,:,:,:)
950  real(RP), intent(out) :: momy(:,:,:,:)
951  real(RP), intent(out) :: rhot(:,:,:,:)
952  real(RP), intent(out) :: qtrc(:,:,:,:,:)
953  character(len=*), intent(in) :: basename_org
954  integer, intent(in) :: dims(6)
955  integer, intent(in) :: mdlid ! model type id
956  logical, intent(in) :: flg_bin ! flag for SBM(S10) is used or not 0-> not used, 1-> used
957  logical, intent(in) :: flg_intrp ! flag for interpolation of SBM(S10) from outer bulk-MP model
958  integer, intent(in) :: mptype_parent ! microphysics type of the parent model (number of classes)
959  integer, intent(in) :: timelen ! time steps in one file
960  integer, intent(in) :: skiplen ! skip steps
961 
962  real(RP) :: velz (KA,IA,JA)
963  real(RP) :: velx (KA,IA,JA)
964  real(RP) :: vely (KA,IA,JA)
965  real(RP) :: llvelx(KA,IA,JA)
966  real(RP) :: llvely(KA,IA,JA)
967  real(RP) :: work (KA,IA,JA)
968  real(RP) :: pott (KA,IA,JA)
969  real(RP) :: temp (KA,IA,JA)
970  real(RP) :: pres (KA,IA,JA)
971 
972  real(RP) :: qc(KA,IA,JA)
973 
974  integer :: k, i, j, iq
975  integer :: n, nn
976  character(len=H_SHORT) :: mptype_run
977  !---------------------------------------------------------------------------
978 
979  if( io_l ) write(io_fid_log,*)
980  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[RealinputAtmos]/Categ[Input]'
981 
982  select case(atmos_phy_mp_type)
983  case("DRY", "NONE")
984  mptype_run = 'dry'
985  case("KESSLER")
986  mptype_run = 'single'
987  case("TOMITA08")
988  mptype_run = 'single'
989  case("SN14")
990  mptype_run = 'double'
991  case("SUZUKI10")
992  mptype_run = 'single-bin'
993  case default
994  write(*,*) 'xxx Unsupported ATMOS_PHY_MP_TYPE (', trim(atmos_phy_mp_type), '). Check!'
995  call prc_mpistop
996  end select
997 
998 
999  if ( do_read_atom ) then
1000 
1001  select case( mdlid )
1002  case( iscale ) ! TYPE: SCALE-RM
1003 
1004  call parentatomopenscale( lon_org, lat_org, & ! (out)
1005  cz_org, & ! (out)
1006  basename_org, & ! (in)
1007  dims ) ! (in)
1008 
1009  case( iwrfarw ) ! TYPE: WRF-ARW
1010 
1012 
1013  case( inicam ) ! TYPE: NICAM-NETCDF
1014 
1015  call parentatomopennicam( lon_org, lat_org, & ! (out)
1016  cz_org, & ! (out)
1017  basename_org, & ! (in)
1018  dims ) ! (in)
1019 
1020  case( igrads ) ! TYPE: GrADS format
1021 
1022  call parentatomopengrads
1023 
1024  end select
1025 
1026  end if
1027 
1028  do n = skiplen+1, timelen
1029  nn = n - skiplen
1030 
1031  if ( do_read_atom ) then
1032 
1033  select case( mdlid )
1034  case( iscale ) ! TYPE: SCALE-RM
1035 
1036  call parentatominputscale( velz_org, velx_org, vely_org, & ! (out)
1037  pres_org, dens_org, pott_org, & ! (out)
1038  qtrc_org, & ! (out)
1039  cz_org, & ! (in)
1040  flg_bin, flg_intrp, & ! (in)
1041  basename_org, mptype_parent, & ! (in)
1042  dims, n ) ! (in)
1043 
1044  case( iwrfarw ) ! TYPE: WRF-ARW
1045 
1046  call parentatominputwrfarw( velz_org, velx_org, vely_org, & ! (out)
1047  pres_org, temp_org, qtrc_org, & ! (out)
1048  lon_org, lat_org, cz_org, & ! (out)
1049  basename_org, mptype_parent, & ! (in)
1050  dims, n ) ! (in)
1051 
1052  case( inicam ) ! TYPE: NICAM-NETCDF
1053 
1054  call parentatominputnicam( velz_org, velx_org, vely_org, & ! (out)
1055  pres_org, temp_org, qtrc_org, & ! (out)
1056  basename_org, dims, n ) ! (in)
1057 
1058  case( igrads ) ! TYPE: GrADS format
1059 
1060  call parentatominputgrads( velz_org, velx_org, vely_org, & ! (out)
1061  pres_org, dens_org, temp_org, & ! (out)
1062  qtrc_org, & ! (out)
1063  lon_org, lat_org, cz_org, & ! (out)
1064  basename_org, dims, n ) ! (in)
1065 
1066  end select
1067 
1068  if ( use_temp ) then
1069  do j = 1, dims(3)
1070  do i = 1, dims(2)
1071  do k = 1, dims(1)+2
1072  call thermodyn_pott( pott_org(k,i,j), & ! [OUT]
1073  temp_org(k,i,j), & ! [IN]
1074  pres_org(k,i,j), & ! [IN]
1075  qtrc_org(k,i,j,:), & ! [IN]
1076  tracer_cv(:), & ! [IN]
1077  tracer_r(:), & ! [IN]
1078  tracer_mass(:) ) ! [IN]
1079  end do
1080  end do
1081  end do
1082  end if
1083 
1084  end if
1085 
1086  if ( serial ) then
1087  call comm_bcast( velz_org, dims(1)+2, dims(2), dims(3) )
1088  call comm_bcast( velx_org, dims(1)+2, dims(2), dims(3) )
1089  call comm_bcast( vely_org, dims(1)+2, dims(2), dims(3) )
1090  call comm_bcast( pott_org, dims(1)+2, dims(2), dims(3) )
1091  call comm_bcast( qtrc_org, dims(1)+2, dims(2), dims(3), qa )
1092  if ( use_file_density ) then
1093  call comm_bcast( dens_org, dims(1)+2, dims(2), dims(3) )
1094  else
1095  call comm_bcast( pres_org, dims(1)+2, dims(2), dims(3) )
1096  end if
1097 
1098  if ( first .or. update_coord ) then
1099 
1100  call comm_bcast( lon_org, dims(2), dims(3) )
1101  call comm_bcast( lat_org, dims(2), dims(3) )
1102  call comm_bcast( cz_org, dims(1)+2, dims(2), dims(3) )
1103 
1104  end if
1105 
1106  end if
1107 
1108  if ( first .or. update_coord ) then
1109 
1110  call intrpnest_domain_compatibility( lon_org(:,:), lat_org(:,:), cz_org(:,:,:), &
1111  lon(:,:), lat(:,:), cz(ks:ke,:,:) )
1112 
1113  ! full level
1114  call intrpnest_interp_fact_llz( hfact, vfact, & ! [OUT]
1115  kgrd, igrd, jgrd, & ! [OUT]
1116  ncopy, & ! [OUT]
1117  cz, lat, lon, & ! [IN]
1118  ks, ke, ia, ja, & ! [IN]
1119  cz_org, lat_org, lon_org, & ! [IN]
1120  dims(1)+2, dims(2), dims(3) ) ! [IN]
1121 
1122  end if
1123 
1124  call intrpnest_interp_3d( velz(:,:,:), &
1125  velz_org(:,:,:), &
1126  hfact(:,:,:), &
1127  vfact(:,:,:,:,:), &
1128  kgrd(:,:,:,:,:), &
1129  igrd(:,:,:), &
1130  jgrd(:,:,:), &
1131  ia, ja, ks, ke-1 )
1132 
1133  call intrpnest_interp_3d( llvelx(:,:,:), &
1134  velx_org(:,:,:), &
1135  hfact(:,:,:), &
1136  vfact(:,:,:,:,:), &
1137  kgrd(:,:,:,:,:), &
1138  igrd(:,:,:), &
1139  jgrd(:,:,:), &
1140  ia, ja, ks, ke )
1141 
1142  call intrpnest_interp_3d( llvely(:,:,:), &
1143  vely_org(:,:,:), &
1144  hfact(:,:,:), &
1145  vfact(:,:,:,:,:), &
1146  kgrd(:,:,:,:,:), &
1147  igrd(:,:,:), &
1148  jgrd(:,:,:), &
1149  ia, ja, ks, ke )
1150 
1151  if ( rotate ) then
1152  ! convert from latlon coordinate to local mapping (x)
1153  do j = 1, ja
1154  do i = 1, ia
1155  do k = ks, ke
1156  work(k,i,j) = llvelx(k,i,j) * rotc(i,j,cosin) + llvely(k,i,j) * rotc(i,j,sine )
1157  end do
1158  end do
1159  end do
1160 
1161  ! from scalar point to staggered point
1162  do j = 1, ja
1163  do i = 1, ia-1
1164  do k = ks, ke
1165  velx(k,i,j) = ( work(k,i+1,j) + work(k,i,j) ) * 0.5_rp
1166  end do
1167  end do
1168  end do
1169  do j = 1, ja
1170  do k = ks, ke
1171  velx(k,ia,j) = work(k,ia,j)
1172  end do
1173  end do
1174  velx(ks-1,:,:) = 0.0_rp
1175  velx(ks-2,:,:) = 0.0_rp
1176  call comm_vars8( velx(:,:,:), 1 )
1177  call comm_wait ( velx(:,:,:), 1, .false. )
1178 
1179  ! convert from latlon coordinate to local mapping (y)
1180  do j = 1, ja
1181  do i = 1, ia
1182  do k = ks, ke
1183  work(k,i,j) = - llvelx(k,i,j) * rotc(i,j,sine ) + llvely(k,i,j) * rotc(i,j,cosin)
1184  end do
1185  end do
1186  end do
1187 
1188  do j = 1, ja-1
1189  do i = 1, ia
1190  do k = ks, ke
1191  vely(k,i,j) = ( work(k,i,j+1) + work(k,i,j) ) * 0.5_rp
1192  end do
1193  end do
1194  end do
1195  do i = 1, ia
1196  do k = ks, ke
1197  vely(k,i,ja) = work(k,i,ja)
1198  end do
1199  end do
1200  vely(ks-1,:,:) = 0.0_rp
1201  vely(ks-2,:,:) = 0.0_rp
1202  call comm_vars8( vely(:,:,:), 1 )
1203  call comm_wait ( vely(:,:,:), 1, .false. )
1204 
1205  else
1206  velx = llvelx
1207  vely = llvely
1208  end if
1209 
1210  if( trim(mptype_run)=='double' .and. mptype_parent <= 6 )then
1211  if( io_l ) write(io_fid_log,*) '--- Diagnose Number Concentration from Mixing Ratio'
1212  call hydrometeor_diagnose_number_concentration( qtrc_org(:,:,:,:) ) ! [inout]
1213  endif
1214 
1215  do j = 1, dims(3)
1216  do i = 1, dims(2)
1217  do k = 1, dims(1)+2
1218  do iq = 1, qa
1219  qtrc_org(k,i,j,iq) = max( qtrc_org(k,i,j,iq), 0.0_rp )
1220  end do
1221  end do
1222  end do
1223  end do
1224 
1225  call intrpnest_interp_3d( pott(:,:,:), &
1226  pott_org(:,:,:), &
1227  hfact(:,:,:), &
1228  vfact(:,:,:,:,:), &
1229  kgrd(:,:,:,:,:), &
1230  igrd(:,:,:), &
1231  jgrd(:,:,:), &
1232  ia, ja, ks, ke )
1233 
1234  do iq = 1, qa
1235  call intrpnest_interp_3d( qtrc(:,:,:,iq,nn), &
1236  qtrc_org(:,:,:,iq), &
1237  hfact(:,:,:), &
1238  vfact(:,:,:,:,:), &
1239  kgrd(:,:,:,:,:), &
1240  igrd(:,:,:), &
1241  jgrd(:,:,:), &
1242  ia, ja, ks, ke )
1243  end do
1244 
1245  if( use_file_density ) then
1246  ! use logarithmic density to interpolate more accurately
1247 
1248  dens_org = log( dens_org )
1249  call intrpnest_interp_3d( dens(:,:,:,nn), &
1250  dens_org(:,:,:), &
1251  hfact(:,:,:), &
1252  vfact(:,:,:,:,:), &
1253  kgrd(:,:,:,:,:), &
1254  igrd(:,:,:), &
1255  jgrd(:,:,:), &
1256  ia, ja, ks, ke, &
1257  logwegt=.true. )
1258  else
1259 
1260  pres_org = log( pres_org )
1261  call intrpnest_interp_3d( pres(:,:,:), &
1262  pres_org(:,:,:), &
1263  hfact(:,:,:), &
1264  vfact(:,:,:,:,:), &
1265  kgrd(:,:,:,:,:), &
1266  igrd(:,:,:), &
1267  jgrd(:,:,:), &
1268  ia, ja, ks, ke, &
1269  logwegt=.true. )
1270 
1271  qc = 0.0_rp
1272 #ifndef DRY
1273  if ( i_qc > 0 ) then
1274  do iq = qls, qle
1275  qc(:,:,:) = qc(:,:,:) + qtrc(:,:,:,iq,nn)
1276  enddo
1277  end if
1278 #endif
1279  ! make density & pressure profile in moist condition
1280  call hydrostatic_buildrho_real( dens(:,:,:,nn), & ! [OUT]
1281  temp(:,:,:), & ! [OUT]
1282  pres(:,:,:), & ! [INOUT]
1283  pott(:,:,:), & ! [IN]
1284  qtrc(:,:,:,i_qv,nn), & ! [IN]
1285  qc(:,:,:) ) ! [IN]
1286 
1287  call comm_vars8( dens(:,:,:,nn), 1 )
1288  call comm_wait ( dens(:,:,:,nn), 1 )
1289 
1290  end if
1291 
1292 
1293  do j = 1, ja
1294  do i = 1, ia
1295  do k = ks, ke-1
1296  momz(k,i,j,nn) = velz(k,i,j) * ( dens(k+1,i,j,nn) + dens(k,i,j,nn) ) * 0.5_rp
1297  end do
1298  end do
1299  end do
1300  do j = 1, ja
1301  do i = 1, ia
1302  momz(ke,i,j,nn) = 0.0_rp
1303  end do
1304  end do
1305  do j = 1, ja
1306  do i = 1, ia
1307  do k = ks, ke
1308  rhot(k,i,j,nn) = pott(k,i,j) * dens(k,i,j,nn)
1309  end do
1310  end do
1311  end do
1312  do j = 1, ja
1313  do i = 1, ia-1
1314  do k = ks, ke
1315  momx(k,i,j,nn) = velx(k,i,j) * ( dens(k,i+1,j,nn) + dens(k,i,j,nn) ) * 0.5_rp
1316  end do
1317  end do
1318  end do
1319  do j = 1, ja
1320  do k = ks, ke
1321  momx(k,ia,j,nn) = velx(k,ia,j) * dens(k,ia,j,nn)
1322  end do
1323  end do
1324  call comm_vars8( momx(:,:,:,nn), 1 )
1325 
1326  do j = 1, ja-1
1327  do i = 1, ia
1328  do k = ks, ke
1329  momy(k,i,j,nn) = vely(k,i,j) * ( dens(k,i,j+1,nn) + dens(k,i,j,nn) ) * 0.5_rp
1330  end do
1331  end do
1332  end do
1333  do i = 1, ia
1334  do k = ks, ke
1335  momy(k,i,ja,nn) = vely(k,i,ja) * dens(k,i,ja,nn)
1336  end do
1337  end do
1338  call comm_vars8( momy(:,:,:,nn), 2 )
1339 
1340  call comm_wait ( momx(:,:,:,nn), 1, .false. )
1341  call comm_wait ( momy(:,:,:,nn), 2, .false. )
1342 
1343  end do
1344 
1345  first = .false.
1346 
1347  return
1348  end subroutine parentatominput
1349 
1350  !-----------------------------------------------------------------------------
1352  subroutine parentatomboundary( &
1353  dens, &
1354  momz, &
1355  momx, &
1356  momy, &
1357  rhot, &
1358  qtrc, &
1359  numsteps, &
1360  update_dt, &
1361  basename, &
1362  title )
1363  use scale_comm, only: &
1364  comm_vars, &
1365  comm_wait
1366  use scale_fileio, only: &
1367  fileio_create, &
1368  fileio_def_var, &
1369  fileio_enddef, &
1370  fileio_write_var
1371  use scale_time, only: &
1372  time_nowdate
1373  use scale_atmos_phy_mp, only: &
1374  qa_mp, &
1375  qs_mp, &
1376  qe_mp
1377  implicit none
1378 
1379  real(RP), intent(in) :: dens(:,:,:,:)
1380  real(RP), intent(in) :: momz(:,:,:,:)
1381  real(RP), intent(in) :: momx(:,:,:,:)
1382  real(RP), intent(in) :: momy(:,:,:,:)
1383  real(RP), intent(in) :: rhot(:,:,:,:)
1384  real(RP), intent(in) :: qtrc(:,:,:,:,:)
1385  real(RP), intent(in) :: update_dt
1386  character(len=*), intent(in) :: basename
1387  character(len=*), intent(in) :: title
1388  integer, intent(in) :: numsteps ! total time steps
1389 
1390  character(len=H_SHORT) :: atmos_boundary_out_dtype = 'DEFAULT'
1391  real(RP), allocatable :: buffer(:,:,:,:)
1392  integer :: nowdate(6)
1393 
1394  integer :: fid, vid(5+QA_MP)
1395  integer :: k, i, j, n, iq
1396  integer :: ts, te
1397  !---------------------------------------------------------------------------
1398 
1399  ts = 1
1400  te = numsteps
1401 
1402  nowdate = time_nowdate
1403  nowdate(1) = nowdate(1)
1404 
1405  allocate( buffer(ka,ia,ja,te-ts+1) )
1406 
1407  if( io_l ) write(io_fid_log,*)
1408  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[RealinputAtmos]/Categ[Boundary]'
1409 
1410  call fileio_create( fid, basename, title, atmos_boundary_out_dtype, nowdate )
1411 
1412  call fileio_def_var( fid, vid(1), 'DENS', 'Reference Density', 'kg/m3', 'ZXYT', &
1413  atmos_boundary_out_dtype, update_dt, numsteps )
1414  call fileio_def_var( fid, vid(2), 'VELZ', 'Reference VELZ', 'm/s', 'ZXYT', &
1415  atmos_boundary_out_dtype, update_dt, numsteps )
1416  call fileio_def_var( fid, vid(3), 'VELX', 'Reference VELX', 'm/s', 'ZXYT', &
1417  atmos_boundary_out_dtype, update_dt, numsteps )
1418  call fileio_def_var( fid, vid(4), 'VELY', 'Reference VELY', 'm/s', 'ZXYT', &
1419  atmos_boundary_out_dtype, update_dt, numsteps )
1420  call fileio_def_var( fid, vid(5), 'POTT', 'Reference PT', 'K', 'ZXYT', &
1421  atmos_boundary_out_dtype, update_dt, numsteps )
1422  do iq = qs_mp, qe_mp
1423  call fileio_def_var( fid, vid(6+iq-qs_mp), tracer_name(iq), 'Reference '//tracer_name(iq), 'kg/kg', 'ZXYT', &
1424  atmos_boundary_out_dtype, update_dt, numsteps )
1425  end do
1426 
1427  call fileio_enddef( fid )
1428 
1429  call fileio_write_var( fid, vid(1), dens(:,:,:,ts:te), 'DENS', 'ZXYT', update_dt )
1430  do n = ts, te
1431  do j = 1, ja
1432  do i = 1, ia
1433  do k = ks, ke
1434  buffer(k,i,j,n-ts+1) = 2.0_rp * momz(k,i,j,n) / ( dens(k+1,i,j,n) + dens(k,i,j,n) )
1435  end do
1436  end do
1437  end do
1438  end do
1439  call fileio_write_var( fid, vid(2), buffer, 'VELZ', 'ZXYT', update_dt )
1440  do n = ts, te
1441  do j = 1, ja
1442  do i = 1, ia-1
1443  do k = ks, ke
1444  buffer(k,i,j,n-ts+1) = 2.0_rp * momx(k,i,j,n) / ( dens(k,i+1,j,n) + dens(k,i,j,n) )
1445  end do
1446  end do
1447  end do
1448  end do
1449  do n = ts, te
1450  buffer(:,ia,:,n-ts+1) = buffer(:,ia-1,:,n-ts+1)
1451  end do
1452  call fileio_write_var( fid, vid(3), buffer, 'VELX', 'ZXYT', update_dt )
1453  do n = ts, te
1454  do j = 1, ja-1
1455  do i = 1, ia
1456  do k = ks, ke
1457  buffer(k,i,j,n-ts+1) = 2.0_rp * momy(k,i,j,n) / ( dens(k,i,j+1,n) + dens(k,i,j,n) )
1458  end do
1459  end do
1460  end do
1461  end do
1462  do n = ts, te
1463  buffer(:,:,ja,n-ts+1) = buffer(:,:,ja-1,n-ts+1)
1464  end do
1465  call fileio_write_var( fid, vid(4), buffer, 'VELY', 'ZXYT', update_dt )
1466  do n = ts, te
1467  do j = 1, ja
1468  do i = 1, ia
1469  do k = ks, ke
1470  buffer(k,i,j,n-ts+1) = rhot(k,i,j,n) / dens(k,i,j,n)
1471  end do
1472  end do
1473  end do
1474  end do
1475  call fileio_write_var( fid, vid(5), buffer, 'POTT', 'ZXYT', update_dt )
1476 
1477  do iq = qs_mp, qe_mp
1478  do n = ts, te
1479  do j = 1, ja
1480  do i = 1, ia
1481  do k = ks, ke
1482  buffer(k,i,j,n-ts+1) = qtrc(k,i,j,iq,n)
1483  end do
1484  end do
1485  end do
1486  end do
1487  call fileio_write_var( fid, vid(6+iq-qs_mp), buffer, tracer_name(iq), 'ZXYT', update_dt )
1488  end do
1489 
1490  deallocate( buffer )
1491 
1492  return
1493  end subroutine parentatomboundary
1494 
1495  !-----------------------------------------------------------------------------
1497  subroutine parentsurfacesetup( &
1498  ldims, odims, &
1499  lmdlid, omdlid, &
1500  timelen, &
1501  basename_land, &
1502  basename_ocean, &
1503  filetype_land, &
1504  filetype_ocean, &
1505  use_file_landwater, &
1506  intrp_land_temp, &
1507  intrp_land_water, &
1508  intrp_land_sfc_temp, &
1509  intrp_ocean_temp, &
1510  intrp_ocean_sfc_temp )
1511  use scale_external_io, only: &
1512  iscale, &
1513  iwrfarw, &
1514  inicam, &
1515  igrads
1516  use mod_realinput_scale, only: &
1519  use mod_realinput_wrfarw, only: &
1522  use mod_realinput_nicam, only: &
1525  use mod_realinput_grads, only: &
1528  implicit none
1529 
1530  integer, intent(out) :: ldims(3) ! dims for land
1531  integer, intent(out) :: odims(2) ! dims for ocean
1532  integer, intent(out) :: lmdlid ! model id for land
1533  integer, intent(out) :: omdlid ! model id for ocean
1534  integer, intent(out) :: timelen ! number of time steps in ocean file
1535  character(len=*), intent(in) :: basename_land
1536  character(len=*), intent(in) :: basename_ocean
1537  character(len=*), intent(in) :: filetype_land
1538  character(len=*), intent(in) :: filetype_ocean
1539  logical, intent(in) :: use_file_landwater ! use land water data from files
1540  character(len=*), intent(in) :: intrp_land_temp
1541  character(len=*), intent(in) :: intrp_land_water
1542  character(len=*), intent(in) :: intrp_land_sfc_temp
1543  character(len=*), intent(in) :: intrp_ocean_temp
1544  character(len=*), intent(in) :: intrp_ocean_sfc_temp
1545  !---------------------------------------------------------------------------
1546 
1547  if( io_l ) write(io_fid_log,*)
1548  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[RealinputSurface]/Categ[Setup]'
1549 
1550  ! Land
1551 
1552  if( lkmax < 4 )then
1553  write(*,*) 'xxx LKMAX less than 4: ', lkmax
1554  write(*,*) 'xxx in Real Case, LKMAX should be set more than 4'
1555  call prc_mpistop
1556  endif
1557 
1558  if( serial_land ) then
1559  if( prc_ismaster ) then
1560  do_read_land = .true.
1561  else
1562  do_read_land = .false.
1563  endif
1564  else
1565  do_read_land = .true.
1566  endif
1567 
1568  select case(trim(filetype_land))
1569  case('SCALE-RM')
1570 
1571  lmdlid = iscale
1572  serial_land = .false.
1573  do_read_land = .true.
1574  call parentlandsetupscale( ldims ) ! (out)
1575  use_waterratio = .false.
1576 
1577  case('WRF-ARW')
1578 
1579  lmdlid = iwrfarw
1580  if ( do_read_land ) call parentlandsetupwrfarw( ldims, & ! (out)
1581  basename_land ) ! (in)
1582  use_waterratio = .false.
1583 
1584  case('NICAM-NETCDF')
1585 
1586  lmdlid = inicam
1587  if ( do_read_land ) call parentlandsetupnicam( ldims, & ! (out)
1588  basename_land ) ! (in)
1589  use_waterratio = .false.
1590 
1591  case('GrADS')
1592 
1593  lmdlid = igrads
1594  if ( do_read_land ) call parentlandsetupgrads( ldims, & ! (out)
1595  use_waterratio, & ! (out)
1596  use_file_landwater, & ! (in)
1597  basename_land ) ! (in)
1598 
1599  case default
1600 
1601  write(*,*) 'xxx Unsupported FILE TYPE:', trim(filetype_land)
1602  call prc_mpistop
1603 
1604  endselect
1605 
1606  if( serial_land ) then
1607  call comm_bcast( ldims(:), 3 )
1608  call comm_bcast( use_waterratio )
1609  endif
1610 
1611 
1612  select case( intrp_land_temp )
1613  case( 'off' )
1614  i_intrp_land_temp = i_intrp_off
1615  case( 'mask' )
1616  i_intrp_land_temp = i_intrp_mask
1617  case( 'fill' )
1618  i_intrp_land_temp = i_intrp_fill
1619  case default
1620  write(*,*) 'xxx INTRP_LAND_TEMP is invalid. ', intrp_land_temp
1621  call prc_mpistop
1622  end select
1623  select case( intrp_land_sfc_temp )
1624  case( 'off' )
1625  i_intrp_land_sfc_temp = i_intrp_off
1626  case( 'mask' )
1627  i_intrp_land_sfc_temp = i_intrp_mask
1628  case( 'fill' )
1629  i_intrp_land_sfc_temp = i_intrp_fill
1630  case default
1631  write(*,*) 'xxx INTRP_LAND_SFC_TEMP is invalid. ', intrp_land_sfc_temp
1632  call prc_mpistop
1633  end select
1634  select case( intrp_land_water )
1635  case( 'off' )
1636  i_intrp_land_water = i_intrp_off
1637  case( 'mask' )
1638  i_intrp_land_water = i_intrp_mask
1639  case( 'fill' )
1640  i_intrp_land_water = i_intrp_fill
1641  case default
1642  write(*,*) 'xxx INTRP_LAND_WATER is invalid. ', intrp_land_water
1643  call prc_mpistop
1644  end select
1645 
1646  select case( lmdlid )
1647  case( iscale, iwrfarw, inicam )
1648  i_intrp_land_temp = i_intrp_mask
1649  i_intrp_land_sfc_temp = i_intrp_mask
1650  i_intrp_land_water = i_intrp_mask
1651  end select
1652 
1653 
1654  ! Ocean
1655 
1656  if( serial_ocean ) then
1657  if( prc_ismaster ) then
1658  do_read_ocean = .true.
1659  else
1660  do_read_ocean = .false.
1661  endif
1662  else
1663  do_read_ocean = .true.
1664  endif
1665 
1666  select case(trim(filetype_ocean))
1667  case('SCALE-RM')
1668 
1669  timelen = -1
1670  omdlid = iscale
1671  serial_ocean = .false.
1672  do_read_ocean = .true.
1673  call parentoceansetupscale( odims )
1674  update_coord = .false.
1675 
1676  case('WRF-ARW')
1677 
1678  omdlid = iwrfarw
1679  if ( do_read_ocean ) call parentoceansetupwrfarw( odims, timelen, & ! (out)
1680  basename_ocean ) ! (in)
1681  update_coord = .true.
1682 
1683  case('NICAM-NETCDF')
1684 
1685  omdlid = inicam
1686  if ( do_read_ocean ) call parentoceansetupnicam( odims, timelen, & ! (out)
1687  basename_ocean ) ! (in)
1688  update_coord = .false.
1689 
1690  case('GrADS')
1691 
1692  omdlid = igrads
1693  if ( do_read_ocean ) call parentoceansetupgrads( odims, timelen, & ! (out)
1694  basename_ocean ) ! (out)
1695  update_coord = .false.
1696 
1697  case default
1698 
1699  write(*,*) 'xxx Unsupported FILE TYPE:', trim(filetype_ocean)
1700  call prc_mpistop
1701 
1702  endselect
1703 
1704  if( serial_ocean ) then
1705  call comm_bcast( odims(:), 2 )
1706  call comm_bcast( timelen )
1707  endif
1708 
1709 
1710  select case( intrp_ocean_temp )
1711  case( 'off' )
1712  i_intrp_ocean_temp = i_intrp_off
1713  case( 'mask' )
1714  i_intrp_ocean_temp = i_intrp_mask
1715  case( 'fill' )
1716  i_intrp_ocean_temp = i_intrp_fill
1717  case default
1718  write(*,*) 'xxx INTRP_OCEAN_TEMP is invalid. ', intrp_ocean_temp
1719  call prc_mpistop
1720  end select
1721  select case( intrp_ocean_sfc_temp )
1722  case( 'off' )
1723  i_intrp_ocean_sfc_temp = i_intrp_off
1724  case( 'mask' )
1725  i_intrp_ocean_sfc_temp = i_intrp_mask
1726  case( 'fill' )
1727  i_intrp_ocean_sfc_temp = i_intrp_fill
1728  case default
1729  write(*,*) 'xxx INTRP_OCEAN_SFC_TEMP is invalid. ', intrp_ocean_sfc_temp
1730  call prc_mpistop
1731  end select
1732 
1733  select case( omdlid )
1734  case( iscale, iwrfarw, inicam )
1735  i_intrp_ocean_temp = i_intrp_mask
1736  i_intrp_ocean_sfc_temp = i_intrp_mask
1737  end select
1738 
1739 
1740  allocate( tw_org( odims(1), odims(2) ) )
1741  allocate( sst_org( odims(1), odims(2) ) )
1742  allocate( albw_org( odims(1), odims(2), 2 ) )
1743  allocate( olon_org( odims(1), odims(2) ) )
1744  allocate( olat_org( odims(1), odims(2) ) )
1745  allocate( omask_org( odims(1), odims(2) ) )
1746 
1747  first = .true.
1748 
1749  return
1750  end subroutine parentsurfacesetup
1751 
1752  !-----------------------------------------------------------------------------
1754  subroutine parentsurfaceinput( &
1755  tg, &
1756  strg, &
1757  lst, &
1758  albg, &
1759  tc_urb, &
1760  qc_urb, &
1761  uc_urb, &
1762  ust, &
1763  albu, &
1764  tw, &
1765  sst, &
1766  albw, &
1767  z0w, &
1768  DENS, &
1769  MOMZ, &
1770  MOMX, &
1771  MOMY, &
1772  RHOT, &
1773  QTRC, &
1774  basename_land, &
1775  basename_ocean, &
1776  mdlid_land, &
1777  mdlid_ocean, &
1778  ldims, &
1779  odims, &
1780  use_file_landwater, &
1781  init_landwater_ratio, &
1782  init_ocean_alb_lw, &
1783  init_ocean_alb_sw, &
1784  init_ocean_z0w, &
1785  intrp_iter_max, &
1786  soilwater_ds2vc_flag, &
1787  elevation_collection, &
1788  timelen, &
1789  skiplen, &
1790  lit )
1791  use scale_comm, only: &
1792  comm_bcast, &
1793  comm_vars8, &
1794  comm_wait
1795  use scale_const, only: &
1796  eps => const_eps, &
1797  undef => const_undef, &
1798  i_sw => const_i_sw, &
1799  i_lw => const_i_lw
1800  use scale_interpolation_nest, only: &
1803  use scale_land_grid, only: &
1804  lcz => grid_lcz
1805  use scale_atmos_thermodyn, only: &
1806  thermodyn_temp_pres => atmos_thermodyn_temp_pres
1807  use scale_atmos_hydrometeor, only: &
1808  i_qv
1809  use scale_landuse, only: &
1810  lsmask_nest => landuse_frac_land
1811  use mod_realinput_scale, only: &
1815  use mod_realinput_wrfarw, only: &
1819  use mod_realinput_nicam, only: &
1823  use mod_realinput_grads, only: &
1827  implicit none
1828 
1829  real(RP), intent(inout) :: tg(LKMAX,IA,JA)
1830  real(RP), intent(inout) :: strg(LKMAX,IA,JA)
1831  real(RP), intent(inout) :: lst(IA,JA)
1832  real(RP), intent(inout) :: albg(IA,JA,2)
1833  real(RP), intent(inout) :: tc_urb(IA,JA)
1834  real(RP), intent(inout) :: qc_urb(IA,JA)
1835  real(RP), intent(inout) :: uc_urb(IA,JA)
1836  real(RP), intent(inout) :: ust(IA,JA)
1837  real(RP), intent(inout) :: albu(IA,JA,2)
1838  real(RP), intent(out) :: tw(:,:,:)
1839  real(RP), intent(out) :: sst(:,:,:)
1840  real(RP), intent(out) :: albw(:,:,:,:)
1841  real(RP), intent(out) :: z0w(:,:,:)
1842  real(RP), intent(in) :: DENS(KA,IA,JA)
1843  real(RP), intent(in) :: MOMZ(KA,IA,JA)
1844  real(RP), intent(in) :: MOMX(KA,IA,JA)
1845  real(RP), intent(in) :: MOMY(KA,IA,JA)
1846  real(RP), intent(in) :: RHOT(KA,IA,JA)
1847  real(RP), intent(in) :: QTRC(KA,IA,JA,QA)
1848  character(len=*), intent(in) :: basename_land
1849  character(len=*), intent(in) :: basename_ocean
1850  integer, intent(in) :: mdlid_land
1851  integer, intent(in) :: mdlid_ocean
1852  integer, intent(in) :: ldims(3)
1853  integer, intent(in) :: odims(2)
1854  logical, intent(in) :: use_file_landwater ! use land water data from files
1855  real(RP), intent(in) :: init_landwater_ratio ! Ratio of land water to storage is constant,
1856  ! if use_file_landwater is ".false."
1857  real(RP), intent(in) :: init_ocean_alb_lw
1858  real(RP), intent(in) :: init_ocean_alb_sw
1859  real(RP), intent(in) :: init_ocean_z0w
1860  integer, intent(in) :: intrp_iter_max
1861  logical, intent(in) :: soilwater_ds2vc_flag
1862  logical, intent(in) :: elevation_collection
1863  integer, intent(in) :: timelen ! time steps in one file
1864  integer, intent(in) :: skiplen ! skip steps
1865  integer, intent(in) :: lit
1866 
1867  ! land
1868  real(RP) :: tg_org (ldims(1),ldims(2),ldims(3))
1869  real(RP) :: strg_org (ldims(1),ldims(2),ldims(3))
1870  real(RP) :: smds_org (ldims(1),ldims(2),ldims(3))
1871 ! real(RP) :: skint_org( ldims(2),ldims(3))
1872  real(RP) :: lst_org ( ldims(2),ldims(3))
1873  real(RP) :: ust_org ( ldims(2),ldims(3))
1874  real(RP) :: albg_org ( ldims(2),ldims(3),2)
1875  real(RP) :: topo_org ( ldims(2),ldims(3))
1876  real(RP) :: lmask_org( ldims(2),ldims(3))
1877  real(RP) :: lz_org (ldims(1) )
1878  real(RP) :: llon_org ( ldims(2),ldims(3))
1879  real(RP) :: llat_org ( ldims(2),ldims(3))
1880 
1881  ! ocean
1882  real(RP) :: tw_org ( odims(1),odims(2))
1883  real(RP) :: sst_org ( odims(1),odims(2))
1884  real(RP) :: albw_org ( odims(1),odims(2),2)
1885  real(RP) :: z0w_org ( odims(1),odims(2))
1886  real(RP) :: olon_org ( odims(1),odims(2))
1887  real(RP) :: olat_org ( odims(1),odims(2))
1888  real(RP) :: omask_org( odims(1),odims(2))
1889  real(RP) :: omask ( odims(1),odims(2))
1890  real(RP) :: lst_ocean( odims(1),odims(2))
1891 
1892  real(RP) :: hfact_o(odims(1),odims(2),itp_nh)
1893  integer :: igrd_o (odims(1),odims(2),itp_nh)
1894  integer :: jgrd_o (odims(1),odims(2),itp_nh)
1895 
1896  real(RP) :: temp
1897  real(RP) :: pres
1898 
1899  integer :: i, j
1900  integer :: n, nn
1901  !---------------------------------------------------------------------------
1902 
1903  if( io_l ) write(io_fid_log,*)
1904  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[RealinputSurface]/Categ[Input]'
1905 
1906  first = .true.
1907 
1908  if ( first ) then ! read land data only once
1909 
1910  if ( do_read_land ) then
1911 
1912  select case( mdlid_land )
1913  case( iscale ) ! TYPE: SCALE-RM
1914 
1915  call parentlandinputscale( &
1916  tg_org, strg_org, & ! (out)
1917  lst_org, ust_org, albg_org, & ! (out)
1918  topo_org, lmask_org, & ! (out)
1919  llon_org, llat_org, lz_org, & ! (out)
1920  basename_land, ldims, & ! (in)
1921  use_file_landwater, lit ) ! (in)
1922 
1923  case( iwrfarw ) ! TYPE: WRF-ARW
1924 
1925  call parentlandinputwrfarw( &
1926  tg_org, strg_org, & ! (out)
1927  lst_org, ust_org, albg_org, & ! (out)
1928  topo_org, lmask_org, & ! (out)
1929  llon_org, llat_org, lz_org, & ! (out)
1930  basename_land, ldims, & ! (in)
1931  use_file_landwater, lit ) ! (in)
1932 
1933  case( inicam ) ! TYPE: NICAM-NETCDF
1934 
1935  call parentlandinputnicam( &
1936  tg_org, strg_org, & ! (out)
1937  lst_org, & ! (out)
1938  llon_org, llat_org, lz_org, & ! (out)
1939  topo_org, lmask_org, & ! (out)
1940  basename_land, ldims, & ! (in)
1941  use_file_landwater, lit ) ! (in)
1942  ust_org = undef
1943  albg_org = undef
1944 
1945  case( igrads ) ! TYPE: GrADS format
1946 
1947  call parentlandinputgrads( &
1948  tg_org, strg_org, smds_org, & ! (out)
1949  lst_org, & ! (out)
1950  llon_org, llat_org, lz_org, & ! (out)
1951  topo_org, lmask_org, & ! (out)
1952  basename_land, ldims, & ! (in)
1953  use_file_landwater, lit ) ! (in)
1954  ust_org = undef
1955  albg_org = undef
1956 
1957  end select
1958 
1959  end if
1960 
1961  if ( serial_land ) then
1962  call comm_bcast( tg_org, ldims(1), ldims(2), ldims(3) )
1963  if ( use_waterratio ) then
1964  call comm_bcast( smds_org, ldims(1), ldims(2), ldims(3) )
1965  else
1966  call comm_bcast( strg_org, ldims(1), ldims(2), ldims(3) )
1967  end if
1968  call comm_bcast( lst_org, ldims(2), ldims(3) )
1969  call comm_bcast( ust_org, ldims(2), ldims(3) )
1970  call comm_bcast( albg_org(:,:,i_lw), ldims(2), ldims(3) )
1971  call comm_bcast( albg_org(:,:,i_sw), ldims(2), ldims(3) )
1972  call comm_bcast( topo_org, ldims(2), ldims(3) )
1973  call comm_bcast( lmask_org, ldims(2), ldims(3) )
1974  call comm_bcast( llon_org, ldims(2), ldims(3) )
1975  call comm_bcast( llat_org, ldims(2), ldims(3) )
1976  call comm_bcast( lz_org, ldims(1) )
1977  end if
1978 
1979 
1980  ! urban data
1981 
1982  do j = 1, ja
1983  do i = 1, ia
1984  call thermodyn_temp_pres( temp, & ! [OUT]
1985  pres, & ! [OUT] not used
1986  dens(ks,i,j), & ! [IN]
1987  rhot(ks,i,j), & ! [IN]
1988  qtrc(ks,i,j,:), & ! [IN]
1989  tracer_cv(:), & ! [IN]
1990  tracer_r(:), & ! [IN]
1991  tracer_mass(:) ) ! [IN]
1992 
1993  tc_urb(i,j) = temp
1994 #ifdef DRY
1995  qc_urb(i,j) = 0.0_rp
1996 #else
1997  qc_urb(i,j) = qtrc(ks,i,j,i_qv)
1998 #endif
1999  enddo
2000  enddo
2001 
2002  do j = 1, ja-1
2003  do i = 1, ia-1
2004  uc_urb(i,j) = max(sqrt( ( momx(ks,i,j) / (dens(ks,i+1, j)+dens(ks,i,j)) * 2.0_rp )**2.0_rp &
2005  + ( momy(ks,i,j) / (dens(ks, i,j+1)+dens(ks,i,j)) * 2.0_rp )**2.0_rp ), &
2006  0.01_rp)
2007  enddo
2008  enddo
2009  do j = 1, ja-1
2010  uc_urb(ia,j) = max(sqrt( ( momx(ks,ia,j) / dens(ks,ia,j ) )**2.0_rp &
2011  + ( momy(ks,ia,j) / (dens(ks,ia,j+1)+dens(ks,ia,j)) * 2.0_rp )**2.0_rp ), &
2012  0.01_rp)
2013  enddo
2014  do i = 1, ia-1
2015  uc_urb(i,ja) = max(sqrt( ( momx(ks,i,ja) / (dens(ks,i+1,ja)+dens(ks,i,ja)) * 2.0_rp )**2.0_rp &
2016  + ( momy(ks,i,ja) / dens(ks,i ,ja) )**2.0_rp ), 0.01_rp)
2017  enddo
2018  uc_urb(ia,ja) = max(sqrt( ( momx(ks,ia,ja) / dens(ks,ia,ja) )**2.0_rp &
2019  + ( momy(ks,ia,ja) / dens(ks,ia,ja) )**2.0_rp ), 0.01_rp)
2020 
2021  call comm_vars8( uc_urb, 1 )
2022  call comm_wait ( uc_urb, 1, .false. )
2023 
2024 
2025  end if ! first
2026 
2027 
2028  if ( do_read_ocean ) then
2029 
2030  select case( mdlid_ocean )
2031  case( iscale ) ! TYPE: SCALE-RM
2032 
2033  call parentoceanopenscale( olon_org, olat_org, & ! (out)
2034  omask_org, & ! (out)
2035  basename_ocean, & ! (in)
2036  odims ) ! (in)
2037 
2038  case( iwrfarw ) ! TYPE: WRF-ARW
2039 
2041 
2042  case( inicam ) ! TYPE: NICAM-NETCDF
2043 
2044  call parentoceanopennicam( olon_org, olat_org, & ! (out)
2045  omask_org, & ! (out)
2046  basename_ocean, & ! (in)
2047  odims ) ! (in)
2048 
2049  case( igrads ) ! TYPE: GrADS format
2050 
2052 
2053  end select
2054 
2055  end if
2056 
2057 
2058  do n = skiplen+1, timelen
2059  nn = n - skiplen
2060 
2061  if ( do_read_ocean ) then
2062 
2063  select case( mdlid_ocean )
2064  case( iscale ) ! TYPE: SCALE-RM
2065 
2066  call parentoceaninputscale( &
2067  tw_org, sst_org, & ! (out)
2068  albw_org, z0w_org, & ! (out)
2069  omask_org, & ! (out)
2070  basename_ocean, odims, & ! (in)
2071  n ) ! (in)
2072 
2073  case( iwrfarw ) ! TYPE: WRF-ARW
2074 
2075  call parentoceaninputwrfarw( &
2076  tw_org, sst_org, & ! (out)
2077  albw_org, z0w_org, & ! (out)
2078  omask_org, & ! (out)
2079  olon_org, olat_org, & ! (out)
2080  basename_ocean, odims, & ! (in)
2081  n ) ! (in)
2082 
2083  case( inicam ) ! TYPE: NICAM-NETCDF
2084 
2085  call parentoceaninputnicam( &
2086  tw_org, sst_org, & ! (out)
2087  basename_ocean, odims, & ! (in)
2088  omask_org, & ! (in)
2089  n ) ! (in)
2090  albw_org = undef
2091  z0w_org = undef
2092 
2093  case( igrads ) ! TYPE: GrADS format
2094 
2095  call parentoceaninputgrads( &
2096  tw_org, sst_org, & ! (out)
2097  omask_org, & ! (out)
2098  olon_org, olat_org, & ! (out)
2099  basename_ocean, odims, & ! (in)
2100  n ) ! (in)
2101  albw_org = undef
2102  z0w_org = undef
2103 
2104  end select
2105 
2106  end if
2107 
2108  if ( serial_ocean ) then
2109  call comm_bcast( tw_org, odims(1), odims(2) )
2110  call comm_bcast( sst_org, odims(1), odims(2) )
2111  call comm_bcast( albw_org(:,:,i_lw), odims(1), odims(2) )
2112  call comm_bcast( albw_org(:,:,i_sw), odims(1), odims(2) )
2113  call comm_bcast( z0w_org, odims(1), odims(2) )
2114  call comm_bcast( omask_org, odims(1), odims(2) )
2115  if ( first .or. update_coord ) then
2116  call comm_bcast( olon_org, odims(1), odims(2) )
2117  call comm_bcast( olat_org, odims(1), odims(2) )
2118  end if
2119  end if
2120 
2121 
2122  if ( first .or. update_coord ) then
2123 
2124  ! interpolation facter between outer ocean grid
2125  call intrpnest_interp_fact_latlon( hfact_o(:,:,:), & ! [OUT]
2126  igrd_o(:,:,:), jgrd_o(:,:,:), & ! [OUT]
2127  olat_org(:,:), olon_org(:,:), & ! [IN]
2128  odims(1), odims(2), & ! [IN]
2129  llat_org(:,:), llon_org(:,:), & ! [IN]
2130  ldims(2), ldims(3) ) ! [IN]
2131 
2132  end if
2133 
2134  ! Ocean temp: interpolate over the land
2135  if ( i_intrp_ocean_temp .ne. i_intrp_off ) then
2136  select case( i_intrp_ocean_temp )
2137  case( i_intrp_mask )
2138  omask = omask_org
2139  case( i_intrp_fill )
2140  call make_mask( omask, tw_org, odims(1), odims(2), landdata=.false.)
2141  end select
2142  call interp_oceanland_data(tw_org, omask, odims(1), odims(2), .false., intrp_iter_max)
2143  end if
2144 
2145  ! SST: interpolate over the land
2146  if ( i_intrp_ocean_sfc_temp .ne. i_intrp_off ) then
2147  select case( i_intrp_ocean_sfc_temp )
2148  case( i_intrp_mask )
2149  omask = omask_org
2150  case( i_intrp_fill )
2151  call make_mask( omask, sst_org, odims(1), odims(2), landdata=.false.)
2152  end select
2153  call interp_oceanland_data(sst_org, omask, odims(1), odims(2), .false., intrp_iter_max)
2154  end if
2155 
2156  if ( first ) then ! interporate land data only once
2157 
2158  call land_interporation( &
2159  tg, strg, & ! (out)
2160  lst, albg, & ! (out)
2161  ust, albu, & ! (out)
2162  tg_org, strg_org, smds_org, & ! (inout)
2163  lst_org, albg_org, & ! (inout)
2164  ust_org, & ! (inout)
2165  sst_org, & ! (in)
2166  lmask_org, & ! (in)
2167  lsmask_nest, & ! (in)
2168  topo_org, & ! (in)
2169  lz_org, llon_org, llat_org, & ! (in)
2170  lcz, lon, lat, & ! (in)
2171  ldims, odims, & ! (in)
2172  maskval_tg, maskval_strg, & ! (in)
2173  init_landwater_ratio, & ! (in)
2174  use_file_landwater, & ! (in)
2175  use_waterratio, & ! (in)
2176  soilwater_ds2vc_flag, & ! (in)
2177  elevation_collection, & ! (in)
2178  intrp_iter_max ) ! (in)
2179 
2180  end if ! first
2181 
2182  if ( first .or. update_coord ) then
2183  ! land surface temperature at ocean grid
2184  call intrpnest_interp_2d( lst_ocean(:,:), lst_org(:,:), hfact_o(:,:,:), &
2185  igrd_o(:,:,:), jgrd_o(:,:,:), odims(1), odims(2) )
2186  end if
2187 
2188  call replace_misval_map( sst_org, lst_ocean, odims(1), odims(2), "SST")
2189  call replace_misval_map( tw_org, lst_ocean, odims(1), odims(2), "OCEAN_TEMP")
2190 
2191  do j = 1, odims(2)
2192  do i = 1, odims(1)
2193  if ( albw_org(i,j,i_lw) == undef ) albw_org(i,j,i_lw) = init_ocean_alb_lw
2194  if ( albw_org(i,j,i_sw) == undef ) albw_org(i,j,i_sw) = init_ocean_alb_sw
2195  if ( z0w_org(i,j) == undef ) z0w_org(i,j) = init_ocean_z0w
2196  end do
2197  end do
2198 
2199 
2200  if ( first .or. update_coord ) then
2201  ! interporation for ocean variables
2202  call intrpnest_interp_fact_latlon( hfact(:,:,:), & ! [OUT]
2203  igrd(:,:,:), jgrd(:,:,:), & ! [OUT]
2204  lat(:,:), lon(:,:), & ! [IN]
2205  ia, ja, & ! [IN]
2206  olat_org(:,:), olon_org(:,:), & ! [IN]
2207  odims(1), odims(2) ) ! [IN]
2208  end if
2209 
2210  call intrpnest_interp_2d( tw(:,:,nn), tw_org(:,:), hfact(:,:,:), &
2211  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2212  call intrpnest_interp_2d( sst(:,:,nn), sst_org(:,:), hfact(:,:,:), &
2213  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2214  call intrpnest_interp_2d( albw(:,:,i_lw,nn), albw_org(:,:,i_lw), hfact(:,:,:), &
2215  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2216  call intrpnest_interp_2d( albw(:,:,i_sw,nn), albw_org(:,:,i_sw), hfact(:,:,:), &
2217  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2218  call intrpnest_interp_2d( z0w(:,:,nn), z0w_org(:,:), hfact(:,:,:), &
2219  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2220 
2221  if ( first ) then
2222 
2223  ! replace values over the ocean ####
2224  do j = 1, ja
2225  do i = 1, ia
2226  if( abs(lsmask_nest(i,j)-0.0_rp) < eps ) then ! ocean grid
2227  lst(i,j) = sst(i,j,nn)
2228  ust(i,j) = sst(i,j,nn)
2229  endif
2230  enddo
2231  enddo
2232 
2233  end if
2234 
2235 
2236  first = .false.
2237 
2238  end do ! time loop
2239 
2240  return
2241  end subroutine parentsurfaceinput
2242 
2244  subroutine parentoceanboundary( &
2245  tw, &
2246  sst, &
2247  albw, &
2248  z0, &
2249  numsteps, &
2250  update_dt, &
2251  basename, &
2252  title )
2253  use scale_const, only: &
2254  i_sw => const_i_sw, &
2255  i_lw => const_i_lw
2256  use scale_fileio, only: &
2257  fileio_create, &
2258  fileio_def_var, &
2259  fileio_enddef, &
2260  fileio_write_var
2261  use scale_time, only: &
2262  time_nowdate
2263  implicit none
2264 
2265  real(RP), intent(in) :: tw(:,:,:)
2266  real(RP), intent(in) :: sst(:,:,:)
2267  real(RP), intent(in) :: albw(:,:,:,:)
2268  real(RP), intent(in) :: z0(:,:,:)
2269  real(RP), intent(in) :: update_dt
2270  character(len=*), intent(in) :: basename
2271  character(len=*), intent(in) :: title
2272  integer, intent(in) :: numsteps ! total time steps
2273 
2274  character(len=H_SHORT) :: ocean_boundary_out_dtype = 'DEFAULT'
2275  integer :: nowdate(6)
2276  integer :: fid, vid(5)
2277  integer :: ts, te
2278  !---------------------------------------------------------------------------
2279 
2280  ts = 1
2281  te = numsteps
2282 
2283  if( io_l ) write(io_fid_log,*)
2284  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[RealinputSurface]/Categ[Boundary]'
2285 
2286  nowdate = time_nowdate
2287  nowdate(1) = nowdate(1)
2288 
2289  call fileio_create( fid, basename, title, ocean_boundary_out_dtype, nowdate )
2290 
2291  call fileio_def_var( fid, vid(1), &
2292  'OCEAN_TEMP', 'Reference Ocean Temperature', 'K', 'XYT', &
2293  ocean_boundary_out_dtype, update_dt, numsteps )
2294  call fileio_def_var( fid, vid(2), &
2295  'OCEAN_SFC_TEMP', 'Reference Ocean Surface Temperature', 'K', 'XYT', &
2296  ocean_boundary_out_dtype, update_dt, numsteps )
2297  call fileio_def_var( fid, vid(3), &
2298  'OCEAN_ALB_LW', 'Reference Ocean Surface Albedo Long-wave', '1', 'XYT', &
2299  ocean_boundary_out_dtype, update_dt, numsteps )
2300  call fileio_def_var( fid, vid(4), &
2301  'OCEAN_ALB_SW', 'Reference Ocean Surface Albedo Short-wave', '1', 'XYT', &
2302  ocean_boundary_out_dtype, update_dt, numsteps )
2303  call fileio_def_var( fid, vid(5), &
2304  'OCEAN_SFC_Z0', 'Reference Ocean Surface Z0', 'm', 'XYT', &
2305  ocean_boundary_out_dtype, update_dt, numsteps )
2306 
2307  call fileio_enddef( fid )
2308 
2309  call fileio_write_var( fid, vid(1), tw(:,:,ts:te), 'OCEAN_TEMP', 'XYT', update_dt )
2310  call fileio_write_var( fid, vid(2), sst(:,:,ts:te), 'OCEAN_SFC_TEMP', 'XYT', update_dt )
2311  call fileio_write_var( fid, vid(3), albw(:,:,i_lw,ts:te), 'OCEAN_ALB_LW', 'XYT', update_dt )
2312  call fileio_write_var( fid, vid(4), albw(:,:,i_sw,ts:te), 'OCEAN_ALB_SW', 'XYT', update_dt )
2313  call fileio_write_var( fid, vid(5), z0(:,:,ts:te), 'OCEAN_SFC_Z0', 'XYT', update_dt )
2314 
2315  return
2316  end subroutine parentoceanboundary
2317 
2318 
2319  !-------------------------------
2320  subroutine land_interporation( &
2321  tg, &
2322  strg, &
2323  lst, &
2324  albg, &
2325  ust, &
2326  albu, &
2327  tg_org, &
2328  strg_org, &
2329  smds_org, &
2330  lst_org, &
2331  albg_org, &
2332  ust_org, &
2333  sst_org, &
2334  lmask_org, &
2335  lsmask_nest, &
2336  topo_org, &
2337  lz_org, &
2338  llon_org, &
2339  llat_org, &
2340  LCZ, &
2341  LON, &
2342  LAT, &
2343  ldims, &
2344  odims, &
2345  maskval_tg, &
2346  maskval_strg, &
2347  init_landwater_ratio, &
2348  use_file_landwater, &
2349  use_waterratio, &
2350  soilwater_ds2vc_flag, &
2351  elevation_collection, &
2352  intrp_iter_max )
2353  use scale_process, only: &
2354  prc_mpistop
2355  use scale_const, only: &
2356  undef => const_undef, &
2357  i_sw => const_i_sw, &
2358  i_lw => const_i_lw, &
2359  laps => const_laps
2360  use scale_interpolation_nest, only: &
2365  use scale_topography, only: &
2366  topo_zsfc
2367  use mod_land_vars, only: &
2369  implicit none
2370  real(RP), intent(out) :: tg(LKMAX,IA,JA)
2371  real(RP), intent(out) :: strg(LKMAX,IA,JA)
2372  real(RP), intent(out) :: lst(IA,JA)
2373  real(RP), intent(out) :: albg(IA,JA,2)
2374  real(RP), intent(out) :: ust(IA,JA)
2375  real(RP), intent(out) :: albu(IA,JA,2)
2376  real(RP), intent(inout) :: tg_org(:,:,:)
2377  real(RP), intent(inout) :: strg_org(:,:,:)
2378  real(RP), intent(inout) :: smds_org(:,:,:)
2379  real(RP), intent(inout) :: lst_org(:,:)
2380  real(RP), intent(inout) :: albg_org(:,:,:)
2381  real(RP), intent(inout) :: ust_org(:,:)
2382  real(RP), intent(inout) :: sst_org(:,:)
2383  real(RP), intent(in) :: lmask_org(:,:)
2384  real(RP), intent(in) :: lsmask_nest(:,:)
2385  real(RP), intent(in) :: topo_org(:,:)
2386  real(RP), intent(in) :: lz_org(:)
2387  real(RP), intent(in) :: llon_org(:,:)
2388  real(RP), intent(in) :: llat_org(:,:)
2389  real(RP), intent(in) :: LCZ(LKMAX)
2390  real(RP), intent(in) :: LON(IA,JA)
2391  real(RP), intent(in) :: LAT(IA,JA)
2392  integer, intent(in) :: ldims(3)
2393  integer, intent(in) :: odims(2)
2394  real(RP), intent(in) :: maskval_tg
2395  real(RP), intent(in) :: maskval_strg
2396  real(RP), intent(in) :: init_landwater_ratio
2397  logical, intent(in) :: use_file_landwater
2398  logical, intent(in) :: use_waterratio
2399  logical, intent(in) :: soilwater_ds2vc_flag
2400  logical, intent(in) :: elevation_collection
2401  integer, intent(in) :: intrp_iter_max
2402 
2403  real(RP) :: lmask(ldims(2), ldims(3))
2404  real(RP) :: smds(LKMAX,IA,JA)
2405 
2406  ! data for interporation
2407  real(RP) :: hfact_l(ldims(2), ldims(3), itp_nh)
2408  integer :: igrd_l (ldims(2), ldims(3), itp_nh)
2409  integer :: jgrd_l (ldims(2), ldims(3), itp_nh)
2410  real(RP) :: vfactl(LKMAX,IA,JA,itp_nh,itp_nv)
2411  integer :: kgrdl (LKMAX,IA,JA,itp_nh,itp_nv)
2412 
2413  real(RP) :: sst_land(ldims(2), ldims(3))
2414  real(RP) :: work(ldims(2), ldims(3))
2415 
2416  real(RP) :: lz3d_org(ldims(1),ldims(2),ldims(3))
2417  real(RP) :: lcz_3D(LKMAX,IA,JA)
2418 
2419  ! elevation collection
2420  real(RP) :: topo(IA,JA)
2421  real(RP) :: tdiff
2422 
2423  integer :: k, i, j
2424 
2425 
2426  ! Surface skin temp: interpolate over the ocean
2427  if ( i_intrp_land_sfc_temp .ne. i_intrp_off ) then
2428  select case( i_intrp_land_sfc_temp )
2429  case( i_intrp_mask )
2430  lmask = lmask_org
2431  case( i_intrp_fill )
2432  call make_mask( lmask, lst_org, ldims(2), ldims(3), landdata=.true.)
2433  case default
2434  write(*,*) 'xxx INTRP_LAND_SFC_TEMP is invalid.'
2435  call prc_mpistop
2436  end select
2437  call interp_oceanland_data(lst_org, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2438  end if
2439 
2440  ! Urban surface temp: interpolate over the ocean
2441  ! if ( i_INTRP_URB_SFC_TEMP .ne. i_intrp_off ) then
2442  ! select case( i_INTRP_URB_SFC_TEMP )
2443  ! case( i_intrp_mask )
2444  ! lmask = lmask_org
2445  ! case( i_intrp_fill )
2446  ! call make_mask( lmask, ust_org, ldims(2), ldims(3), landdata=.true.)
2447  ! case default
2448  ! write(*,*) 'xxx INTRP_URB_SFC_TEMP is invalid.'
2449  ! call PRC_MPIstop
2450  ! end select
2451  ! call interp_OceanLand_data(ust_org, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2452  !end if
2453 
2454  ! interpolation facter between outer land grid and ocean grid
2455  call intrpnest_interp_fact_latlon( hfact_l(:,:,:), & ! [OUT]
2456  igrd_l(:,:,:), jgrd_l(:,:,:), & ! [OUT]
2457  llat_org(:,:), llon_org(:,:), & ! [IN]
2458  ldims(2), ldims(3), & ! [IN]
2459  olat_org(:,:), olon_org(:,:), & ! [IN]
2460  odims(1), odims(2) ) ! [IN]
2461 
2462 
2463  ! sst on land grid
2464  call intrpnest_interp_2d( sst_land(:,:), sst_org(:,:), hfact_l(:,:,:), &
2465  igrd_l(:,:,:), jgrd_l(:,:,:), ldims(2), ldims(3) )
2466  call replace_misval_map( lst_org, sst_land, ldims(2), ldims(3), "SKINT")
2467 
2468  ! replace missing value
2469  do j = 1, ldims(3)
2470  do i = 1, ldims(2)
2471  if ( ust_org(i,j) == undef ) ust_org(i,j) = lst_org(i,j)
2472 ! if ( skinw_org(i,j) == UNDEF ) skinw_org(i,j) = 0.0_RP
2473 ! if ( snowq_org(i,j) == UNDEF ) snowq_org(i,j) = 0.0_RP
2474 ! if ( snowt_org(i,j) == UNDEF ) snowt_org(i,j) = TEM00
2475  if ( albg_org(i,j,i_lw) == undef ) albg_org(i,j,i_lw) = 0.03_rp ! emissivity of general ground surface : 0.95-0.98
2476  if ( albg_org(i,j,i_sw) == undef ) albg_org(i,j,i_sw) = 0.22_rp
2477  end do
2478  end do
2479 
2480  ! Land temp: interpolate over the ocean
2481  if ( i_intrp_land_temp .ne. i_intrp_off ) then
2482  do k = 1, ldims(1)
2483  work(:,:) = tg_org(k,:,:)
2484  select case( i_intrp_land_temp )
2485  case( i_intrp_mask )
2486  lmask = lmask_org
2487  case( i_intrp_fill )
2488  call make_mask( lmask, work, ldims(2), ldims(3), landdata=.true.)
2489  end select
2490  call interp_oceanland_data( work, lmask, ldims(2), ldims(3), .true., intrp_iter_max )
2491  !replace land temp using skin temp
2492  call replace_misval_map( work, lst_org, ldims(2), ldims(3), "STEMP")
2493  tg_org(k,:,:) = work(:,:)
2494  end do
2495  end if
2496 
2497 
2498  ! fill grid data
2499  do j = 1, ldims(3)
2500  do i = 1, ldims(2)
2501  lz3d_org(:,i,j) = lz_org(:)
2502  end do
2503  end do
2504 
2505  do j = 1, ja
2506  do i = 1, ia
2507  lcz_3d(:,i,j) = lcz(:)
2508  enddo
2509  enddo
2510 
2511  call intrpnest_interp_fact_llz( hfact(:,:,:), & ! [OUT]
2512  vfactl(:,:,:,:,:), & ! [OUT]
2513  kgrdl(:,:,:,:,:), & ! [OUT]
2514  igrd(:,:,:), jgrd(:,:,:), & ! [OUT]
2515  ncopy(:,:,:), & ! [OUT]
2516  lcz_3d(:,:,:), & ! [IN]
2517  lat(:,:), lon(:,:), & ! [IN]
2518  1, lkmax, ia, ja, & ! [IN]
2519  lz3d_org(:,:,:), & ! [IN]
2520  llat_org(:,:), llon_org(:,:), & ! [IN]
2521  ldims(1), ldims(2), ldims(3), & ! [IN]
2522  landgrid=.true. ) ! [IN]
2523 
2524  call intrpnest_interp_2d( lst(:,:), lst_org(:,:), hfact(:,:,:), &
2525  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2526  call intrpnest_interp_2d( ust(:,:), ust_org(:,:), hfact(:,:,:), &
2527  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2528 ! call INTRPNEST_interp_2d( skinw(:,:), skinw_org(:,:), hfact(:,:,:), &
2529 ! igrd(:,:,:), jgrd(:,:,:), IA, JA )
2530 ! call INTRPNEST_interp_2d( snowq(:,:), snowq_org(:,:), hfact(:,:,:), &
2531 ! igrd(:,:,:), jgrd(:,:,:), IA, JA )
2532 ! call INTRPNEST_interp_2d( snowt(:,:), snowt_org(:,:), hfact(:,:,:), &
2533 ! igrd(:,:,:), jgrd(:,:,:), IA, JA )
2534  call intrpnest_interp_2d( albg(:,:,i_lw), albg_org(:,:,i_lw), hfact(:,:,:), &
2535  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2536  call intrpnest_interp_2d( albg(:,:,i_sw), albg_org(:,:,i_sw), hfact(:,:,:), &
2537  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2538 
2539  call intrpnest_interp_3d( tg(:,:,:), &
2540  tg_org(:,:,:), &
2541  hfact(:,:,:), &
2542  vfactl(:,:,:,:,:), &
2543  kgrdl(:,:,:,:,:), &
2544  igrd(:,:,:), &
2545  jgrd(:,:,:), &
2546  ia, ja, 1, lkmax-1 )
2547 
2548  do j = 1, ja
2549  do i = 1, ia
2550  tg(lkmax,i,j) = tg(lkmax-1,i,j)
2551  enddo ! i
2552  enddo ! j
2553 
2554  ! replace values over the ocean
2555  do k = 1, lkmax
2556  call replace_misval_const( tg(k,:,:), maskval_tg, lsmask_nest )
2557  enddo
2558 
2559 
2560  ! elevation collection
2561  if ( elevation_collection ) then
2562  call intrpnest_interp_2d( topo(:,:), topo_org(:,:), hfact(:,:,:), &
2563  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2564 
2565  do j = 1, ja
2566  do i = 1, ia
2567  if ( topo(i,j) > 0.0_rp ) then ! ignore UNDEF value
2568  tdiff = ( topo_zsfc(i,j) - topo(i,j) ) * laps
2569  lst(i,j) = lst(i,j) - tdiff
2570  ust(i,j) = ust(i,j) - tdiff
2571  do k = 1, lkmax
2572  tg(k,i,j) = tg(k,i,j) - tdiff
2573  end do
2574  end if
2575  end do
2576  end do
2577  end if
2578 
2579 
2580 
2581  ! Land water: interpolate over the ocean
2582  if( use_file_landwater )then
2583 
2584  if ( use_waterratio ) then
2585 
2586  if ( i_intrp_land_water .ne. i_intrp_off ) then
2587  do k = 1, ldims(1)
2588  work(:,:) = smds_org(k,:,:)
2589  select case( i_intrp_land_water )
2590  case( i_intrp_mask )
2591  lmask = lmask_org
2592  case( i_intrp_fill )
2593  call make_mask( lmask, work, ldims(2), ldims(3), landdata=.true.)
2594  end select
2595  call interp_oceanland_data(work, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2596  lmask(:,:) = init_landwater_ratio
2597  !replace missing value to init_landwater_ratio
2598  call replace_misval_map( work, lmask, ldims(2), ldims(3), "SMOISDS")
2599  smds_org(k,:,:) = work(:,:)
2600  enddo
2601  end if
2602 
2603  call intrpnest_interp_3d( smds(:,:,:), &
2604  smds_org(:,:,:), &
2605  hfact(:,:,:), &
2606  vfactl(:,:,:,:,:), &
2607  kgrdl(:,:,:,:,:), &
2608  igrd(:,:,:), &
2609  jgrd(:,:,:), &
2610  ia, ja, 1, lkmax-1 )
2611  do k = 1, lkmax-1
2612  strg(k,:,:) = convert_ws2vwc( smds(k,:,:), critical=soilwater_ds2vc_flag )
2613  end do
2614 
2615  else
2616 
2617  if ( i_intrp_land_water .ne. i_intrp_off ) then
2618  do k = 1, ldims(1)
2619  work(:,:) = strg_org(k,:,:)
2620  select case( i_intrp_land_water )
2621  case( i_intrp_mask )
2622  lmask = lmask_org
2623  case( i_intrp_fill )
2624  call make_mask( lmask, work, ldims(2), ldims(3), landdata=.true.)
2625  end select
2626  call interp_oceanland_data(work, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2627  lmask(:,:) = maskval_strg
2628  !replace missing value to init_landwater_ratio
2629  call replace_misval_map( work, lmask, ldims(2), ldims(3), "SMOIS")
2630  strg_org(k,:,:) = work(:,:)
2631  enddo
2632  end if
2633 
2634  call intrpnest_interp_3d( strg(:,:,:), &
2635  strg_org(:,:,:), &
2636  hfact(:,:,:), &
2637  vfactl(:,:,:,:,:), &
2638  kgrdl(:,:,:,:,:), &
2639  igrd(:,:,:), &
2640  jgrd(:,:,:), &
2641  ia, ja, 1, lkmax-1 )
2642  end if
2643 
2644  ! replace values over the ocean
2645  do k = 1, lkmax-1
2646  call replace_misval_const( strg(k,:,:), maskval_strg, lsmask_nest )
2647  enddo
2648 
2649  do j = 1, ja
2650  do i = 1, ia
2651  strg(lkmax,i,j) = strg(lkmax-1,i,j)
2652  enddo
2653  enddo
2654 
2655  else ! not read from boundary file
2656 
2657  smds(:,:,:) = init_landwater_ratio
2658  ! conversion from water saturation [fraction] to volumetric water content [m3/m3]
2659  do k = 1, lkmax
2660  strg(k,:,:) = convert_ws2vwc( smds(k,:,:), critical=.true. )
2661  end do
2662 
2663  endif ! use_file_waterratio
2664 
2665 
2666  ! copy albedo of land to urban
2667  do j = 1, ja
2668  do i = 1, ia
2669  albu(i,j,:) = albg(i,j,:)
2670  enddo
2671  enddo
2672 
2673 
2674  return
2675  end subroutine land_interporation
2676 
2677  !-------------------------------
2678  subroutine make_mask( &
2679  gmask, & ! (out)
2680  data, & ! (in)
2681  nx, & ! (in)
2682  ny, & ! (in)
2683  landdata ) ! (in)
2684  use scale_const, only: &
2685  eps => const_eps, &
2686  undef => const_undef
2687  implicit none
2688  real(RP), intent(out) :: gmask(:,:)
2689  real(RP), intent(in) :: data(:,:)
2690  integer, intent(in) :: nx
2691  integer, intent(in) :: ny
2692  logical, intent(in) :: landdata ! .true. => land data , .false. => ocean data
2693 
2694  real(RP) :: dd
2695  integer :: i,j
2696 
2697  if( landdata )then
2698  gmask(:,:) = 1.0_rp ! gmask=1 will be skip in "interp_OceanLand_data"
2699  dd = 0.0_rp
2700  else
2701  gmask(:,:) = 0.0_rp ! gmask=0 will be skip in "interp_OceanLand_data"
2702  dd = 1.0_rp
2703  endif
2704 
2705  do j = 1, ny
2706  do i = 1, nx
2707  if( abs(data(i,j) - undef) < sqrt(eps) )then
2708  gmask(i,j) = dd
2709  endif
2710  enddo
2711  enddo
2712 
2713  return
2714  end subroutine make_mask
2715  !-----------------------------------------------------------------------------
2716  subroutine interp_oceanland_data( &
2717  data, &
2718  lsmask, &
2719  nx, &
2720  ny, &
2721  landdata, &
2722  iter_max )
2723  use scale_const, only: &
2724  eps => const_eps
2725  implicit none
2726 
2727  integer, intent(in) :: nx
2728  integer, intent(in) :: ny
2729  real(RP), intent(inout) :: data (nx,ny)
2730  real(RP), intent(in) :: lsmask(nx,ny)
2731  logical, intent(in) :: landdata ! .true. => land data , .false. => ocean data
2732  integer, intent(in) :: iter_max
2733 
2734  integer :: mask (nx,ny)
2735  integer :: mask_prev(nx,ny)
2736  real(RP) :: data_prev(nx,ny)
2737  real(RP) :: tmp, cnt, sw
2738  integer :: mask_target
2739 
2740  integer :: num_land, num_ocean, num_replaced
2741  integer :: istr, iend, jstr, jend
2742  integer :: i, j, ii, jj, ite
2743  !---------------------------------------------------------------------------
2744 
2745  if( io_l ) write(io_fid_log,*)
2746  if( io_l ) write(io_fid_log,*) '*** [interp_OceanLand_data]/Categ[realinit]'
2747 
2748  if ( landdata ) then
2749  if( io_l ) write(io_fid_log,*) '*** target mask : LAND'
2750  mask_target = 1 ! interpolation for land data
2751  else
2752  if( io_l ) write(io_fid_log,*) '*** target mask : OCEAN'
2753  mask_target = 0 ! interpolation for ocean data
2754  endif
2755 
2756  ! search target cell for interpolation
2757  num_land = 0
2758  num_ocean = 0
2759  do j = 1, ny
2760  do i = 1, nx
2761  mask(i,j) = int( 0.5_rp - sign(0.5_rp,abs(lsmask(i,j)-1.0_rp)-eps) ) ! 1 for land, 0 for ocean
2762  num_land = num_land + ( mask(i,j) )
2763  num_ocean = num_ocean + ( 1-mask(i,j) )
2764  enddo
2765  enddo
2766 
2767  if( io_l ) write(io_fid_log,'(1x,A,I3.3,A,3I8,A,I8)') '*** ite = ', 0, &
2768  ', (land,ocean,replaced) = ', num_land, num_ocean, 0, ' / ', nx*ny
2769 
2770  ! start interpolation
2771  do ite = 1, iter_max
2772  ! save previous state
2773  mask_prev(:,:) = mask(:,:)
2774  data_prev(:,:) = data(:,:)
2775  num_replaced = 0
2776 
2777  do j = 1, ny
2778  do i = 1, nx
2779 
2780  if( mask(i,j) == mask_target ) cycle ! already filled
2781 
2782  ! collect neighbor grid
2783  istr = max(i-1,1 )
2784  iend = min(i+1,nx)
2785  jstr = max(j-1,1 )
2786  jend = min(j+1,ny)
2787 
2788  tmp = 0.0_rp
2789  cnt = 0.0_rp
2790  do jj = jstr, jend
2791  do ii = istr, iend
2792  sw = 0.5_rp - sign(0.5_rp,real(abs(mask_prev(ii,jj)-mask_target),kind=rp)-eps)
2793 
2794  tmp = tmp + sw * data_prev(ii,jj)
2795  cnt = cnt + sw
2796  enddo
2797  enddo
2798 
2799  if ( cnt >= 3.0_rp ) then ! replace by average of neighbor grid value
2800  data(i,j) = tmp / cnt
2801  mask(i,j) = mask_target
2802 
2803  num_replaced = num_replaced + 1
2804  endif
2805 
2806  enddo
2807  enddo
2808 
2809  if ( landdata ) then
2810  num_land = num_land + num_replaced
2811  num_ocean = num_ocean - num_replaced
2812  else
2813  num_land = num_land - num_replaced
2814  num_ocean = num_ocean + num_replaced
2815  endif
2816  if( io_l ) write(io_fid_log,'(1x,A,I3.3,A,3I8,A,I8)') '*** ite = ', ite, &
2817  ', (land,ocean,replaced) = ', num_land, num_ocean, num_replaced, ' / ', nx*ny
2818 
2819  if( num_replaced == 0 ) exit
2820 
2821  enddo ! itelation
2822 
2823 
2824  return
2825  end subroutine interp_oceanland_data
2826 
2827  !-----------------------------------------------------------------------------
2828  subroutine replace_misval_const( data, maskval, frac_land )
2829  use scale_const, only: &
2830  eps => const_eps
2831  implicit none
2832  real(RP), intent(inout) :: data(:,:)
2833  real(RP), intent(in) :: maskval
2834  real(RP), intent(in) :: frac_land(:,:)
2835  integer :: i, j
2836 
2837  do j = 1, ja
2838  do i = 1, ia
2839  if( abs(frac_land(i,j)-0.0_rp) < eps )then ! ocean grid
2840  data(i,j) = maskval
2841  endif
2842  enddo
2843  enddo
2844 
2845  end subroutine replace_misval_const
2846 
2847  !-----------------------------------------------------------------------------
2848  subroutine replace_misval_map( data, maskval, nx, ny, elem)
2849  use scale_const, only: &
2850  eps => const_eps, &
2851  undef => const_undef
2852  implicit none
2853 
2854  real(RP), intent(inout) :: data(:,:)
2855  real(RP), intent(in) :: maskval(:,:)
2856  integer, intent(in) :: nx, ny
2857  character(len=*), intent(in) :: elem
2858 
2859  integer :: i, j
2860 
2861  do j = 1, ny
2862  do i = 1, nx
2863  if( abs(data(i,j) - undef) < sqrt(eps) )then
2864  if( abs(maskval(i,j) - undef) < sqrt(eps) )then
2865  write(*,*) "xxx data for mask of "//trim(elem)//"(",i,",",j,") includes missing value."
2866  write(*,*) "xxx Please check input data of SKINTEMP or SST. "
2867  call prc_mpistop
2868  else
2869  data(i,j) = maskval(i,j)
2870  endif
2871  endif
2872  enddo
2873  enddo
2874 
2875  end subroutine replace_misval_map
2876 
2877 end module mod_realinput
module ATMOS admin
subroutine, public parentatomsetupgrads(dims, basename)
Atmos Setup.
subroutine replace_misval_map(data, maskval, nx, ny, elem)
real(rp), dimension(:,:,:), allocatable, target, public momz
subroutine make_mask(gmask, data, nx, ny, landdata)
subroutine, public intrpnest_domain_compatibility(lon_org, lat_org, lev_org, lon_loc, lat_loc, lev_loc, skip_x, skip_y, skip_z)
subroutine, public parentoceaninputgrads(tw_org, sst_org, omask_org, olon_org, olat_org, basename_num, odims, nt)
logical, public prc_ismaster
master process in local communicator?
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
subroutine, public parentatomopenwrfarw
subroutine, public parentoceanopengrads
real(rp), dimension(:,:,:), allocatable, public ocean_sfc_albedo
ocean surface albedo (0-1)
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:95
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public urban_qc
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
subroutine, public parentlandinputwrfarw(tg_org, sh2o_org, lst_org, ust_org, albg_org, topo_org, lmask_org, llon_org, llat_org, lz_org, basename, ldims, use_file_landwater, it)
module REAL input WRF-ARW
subroutine, public parentoceanopennicam(olon_org, olat_org, omask_org, basename_num, odims)
module GRID (nesting system)
subroutine, public parentatomopennicam(lon_org, lat_org, cz_org, basename_num, dims)
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(qa_max), public tracer_r
subroutine, public parentatominputnicam(velz_org, velx_org, vely_org, pres_org, temp_org, qtrc_org, basename_num, dims, it)
integer, parameter, public inicam
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
subroutine, public parentlandsetupwrfarw(ldims, basename_land)
Land Setup.
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0e
ocean surface roughness length for vapor [m]
subroutine, public parentoceanopenscale(olon_org, olat_org, omask_org, basename_ocean, odims)
real(rp), dimension(:,:), allocatable, public ocean_temp
temperature at uppermost ocean layer [K]
module ATMOSPHERE / Physics Cloud Microphysics
subroutine land_interporation(tg, strg, lst, albg, ust, albu, tg_org, strg_org, smds_org, lst_org, albg_org, ust_org, sst_org, lmask_org, lsmask_nest, topo_org, lz_org, llon_org, llat_org, LCZ, LON, LAT, ldims, odims, maskval_tg, maskval_strg, init_landwater_ratio, use_file_landwater, use_waterratio, soilwater_ds2vc_flag, elevation_collection, intrp_iter_max)
module ATMOSPHERIC Variables
subroutine parentatomsetup(dims, timelen, mdlid, basename_org, filetype, use_file_density_in, serial_in)
Atmos Setup.
real(rp), dimension(:,:,:), allocatable, target, public momx
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
integer, public qa
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
integer, parameter, public igrads
procedure(intrpnest_intfc_interp_2d), pointer, public intrpnest_interp_2d
subroutine, public parentlandsetupscale(ldims)
Land Setup.
real(rp), dimension(:,:), allocatable, public urban_tb
integer, public nest_interp_level
horizontal interpolation level
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:60
subroutine, public parentoceaninputwrfarw(tw_org, sst_org, albw_org, z0w_org, omask_org, olon_org, olat_org, basename, odims, it)
module GRID (cartesian) for land
integer, parameter, public iwrfarw
integer, parameter, public iscale
subroutine, public intrpnest_interp_fact_latlon(hfact, igrd, jgrd, mylat, mylon, myIA, myJA, inlat, inlon, inIA, inJA)
module URBAN Variables
real(rp), dimension(:,:), allocatable, public urban_raing
subroutine, public atmos_hydrometeor_diagnose_number_concentration(QTRC)
character(len=h_short), dimension(qa_max), public tracer_name
subroutine, public parentlandsetupgrads(ldims, use_waterratio, use_file_landwater, basename)
Land Setup.
real(rp), dimension(qa_max), public tracer_cv
module FILE I/O (netcdf)
real(rp), dimension(:,:), allocatable, public urban_uc
module FILE I/O (netcdf)
subroutine, public intrpnest_interp_fact_llz(hfact, vfact, kgrd, igrd, jgrd, ncopy, myhgt, mylat, mylon, myKS, myKE, myIA, myJA, inhgt, inlat, inlon, inKA, inIA, inJA, landgrid)
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
real(rp), dimension(:,:), allocatable, public land_sfc_temp
land surface skin temperature [K]
subroutine, public parentatomopenscale(lon_org, lat_org, cz_org, basename_org, dims)
module ATMOSPHERIC Surface Variables
subroutine, public parentoceansetupnicam(odims, timelen, basename_org)
Ocean Setup.
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:62
subroutine, public parentoceansetupwrfarw(odims, timelen, basename_org)
Ocean Setup.
module TRACER
module Index
Definition: scale_index.F90:14
real(rp), dimension(:,:), allocatable, public urban_tr
module REAL input NICAM
module REAL input
integer, public ia
of whole cells: x, local, with HALO
real(rp), dimension(:,:,:), allocatable, public urban_tgl
module GRIDTRANS
subroutine, public parentoceaninputscale(tw_org, sst_org, albw_org, z0w_org, omask_org, basename_ocean, odims, it)
module GRID (real space)
module LANDUSE
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp) function, dimension(ia, ja), public convert_ws2vwc(WS, critical)
conversion from water saturation [fraction] to volumetric water content [m3/m3]
subroutine, public parentlandinputnicam(tg_org, strg_org, lst_org, llon_org, llat_org, lz_org, topo_org, lmask_org, basename_num, ldims, use_file_landwater, it)
integer, public ka
of whole cells: z, local, with HALO
subroutine, public parentatomsetupscale(dims)
Atmos Setup.
subroutine, public parentatomsetupnicam(dims, timelen, basename_org)
Atmos Setup.
module REAL input SCALE
subroutine, public parentlandsetupnicam(ldims, basename_org)
Land Setup.
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
subroutine, public fileio_create(fid, basename, title, datatype, date, subsec, append, nozcoord)
Create/open a netCDF file.
real(rp), dimension(:,:,:), allocatable, public land_temp
temperature of each soil layer [K]
subroutine, public parentlandinputscale(tg_org, strg_org, lst_org, ust_org, albg_org, topo_org, lmask_org, llon_org, llat_org, lz_org, basename_land, ldims, use_file_landwater, it)
module COMMUNICATION
Definition: scale_comm.F90:23
module ATMOSPHERE / Hydrostatic barance
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
module TIME
Definition: scale_time.F90:15
real(rp), dimension(:,:), allocatable, public urban_roff
module PROCESS
subroutine, public realinput_atmos(flg_intrp)
subroutine, public parentlandinputgrads(tg_org, strg_org, smds_org, lst_org, llon_org, llat_org, lz_org, topo_org, lmask_org, basename_num, ldims, use_file_landwater, nt)
real(rp), dimension(:,:,:), allocatable, public gtrans_rotc
rotation coefficient
real(rp), dimension(:,:,:), allocatable, public land_sfc_albedo
land surface albedo (0-1)
procedure(intrpnest_intfc_interp_3d), pointer, public intrpnest_interp_3d
module LAND Variables
subroutine, public parentoceansetupgrads(odims, timelen, basename)
Ocean Setup.
subroutine, public parentatominputwrfarw(velz_org, llvelx_org, llvely_org, pres_org, temp_org, qtrc_org, lon_org, lat_org, cz_org, basename, mptype_parent, dims, it)
subroutine, public parentatomsetupwrfarw(dims, timelen, basename_org)
Atmos Setup.
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(:,:,:), allocatable, public urban_sfc_albedo
character(len=h_short), public atmos_phy_mp_type
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:,:), allocatable, public urban_tc
real(rp), dimension(:,:,:), allocatable, target, public momy
subroutine, public fileio_enddef(fid)
Exit netCDF file define mode.
real(rp), dimension(:,:), allocatable, public urban_rainr
module INTERPOLATION (nesting system)
real(rp), dimension(:,:,:), allocatable, public land_water
moisture of each soil layer [m3/m3]
subroutine, public parentatominputscale(velz_org, velx_org, vely_org, pres_org, dens_org, pott_org, qtrc_org, cz_org, flg_bin, flg_intrp, basename_org, mptype_parent, dims, it)
subroutine, public parentatominputgrads(velz_org, velx_org, vely_org, pres_org, dens_org, temp_org, qtrc_org, lon_org, lat_org, cz_org, basename_num, dims, nt)
real(rp), public const_eps
small number
Definition: scale_const.F90:36
module ATMOSPHERE / Thermodynamics
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0m
ocean surface roughness length for momentum [m]
real(rp), dimension(:,:), allocatable, public real_lon
longitude [rad,0-2pi]
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:96
subroutine, public parentoceansetupscale(odims)
Ocean Setup.
module PRECISION
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
subroutine, public fileio_def_var(fid, vid, varname, desc, unit, axistype, datatype, timeintv, nsteps)
Define a variable to file.
real(rp), dimension(:), allocatable, public grid_lcz
center coordinate [m]: z, local=global
module REAL input GrADS
module TOPOGRAPHY
real(rp), dimension(:,:), allocatable, public urban_tg
real(rp), dimension(:,:), allocatable, public ocean_sfc_temp
ocean surface skin temperature [K]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
module land grid index
real(rp), dimension(:,:,:), allocatable, public urban_trl
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
subroutine, public parentoceanopenwrfarw
real(rp), dimension(:,:), allocatable, public real_lat
latitude [rad,-pi,pi]
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public realinput_surface
integer, parameter, public rp
real(rp), dimension(:,:,:), allocatable, public urban_tbl
real(rp), dimension(:,:), allocatable, public real_lonx
longitude at staggered point (uy) [rad,0-2pi]
real(rp), dimension(:,:), allocatable, public landuse_frac_land
land fraction
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0h
ocean surface roughness length for heat [m]
module OCEAN Variables
integer, public io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
real(rp), dimension(:,:), allocatable, public urban_rainb
subroutine, public parentatomopengrads
subroutine replace_misval_const(data, maskval, frac_land)
real(rp), dimension(qa_max), public tracer_mass
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
subroutine, public parentoceaninputnicam(tw_org, sst_org, basename_num, odims, omask_org, it)
module urban grid index
real(rp), dimension(:,:), allocatable, public real_laty
latitude at staggered point (xv) [rad,-pi,pi]
integer, public ja
of whole cells: y, local, with HALO