SCALE-RM
Functions/Subroutines
mod_realinput Module Reference

module REAL input More...

Functions/Subroutines

subroutine, public realinput_atmos (flg_intrp)
 
subroutine, public realinput_surface
 
subroutine parentatomsetup (dims, timelen, mdlid, basename_org, filetype, use_file_density_in, serial_in)
 Atmos Setup. More...
 
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)
 
subroutine make_mask (gmask, data, nx, ny, landdata)
 
subroutine replace_misval_const (data, maskval, frac_land)
 
subroutine replace_misval_map (data, maskval, nx, ny, elem)
 

Detailed Description

module REAL input

Description
read data from file for real atmospheric simulations
Author
Team SCALE
History
  • 2014-04-28 (R.Yoshida) [new]
NAMELIST
  • PARAM_MKINIT_REAL_ATMOS
    nametypedefault valuecomment
    NUMBER_OF_FILES integer 1
    NUMBER_OF_TSTEPS integer 1 num of time steps in one file
    NUMBER_OF_SKIP_TSTEPS integer 0 num of skipped first several data
    FILETYPE_ORG character(len=H_LONG) ''
    BASENAME_ORG character(len=*)
    BASENAME_ADD_NUM logical .false.
    BASENAME_BOUNDARY character(len=H_LONG) 'boundary_atmos'
    BOUNDARY_TITLE character(len=H_LONG) 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
    BOUNDARY_UPDATE_DT real(RP) 0.0_RP inteval time of boudary data update [s]
    PARENT_MP_TYPE integer 6 microphysics type of the parent model (number of classes)
    SERIAL_PROC_READ logical .true. read by one MPI process and broadcast
    USE_FILE_DENSITY logical .false. use density data from files

  • PARAM_MKINIT_REAL_LAND
    nametypedefault valuecomment
    NUMBER_OF_FILES integer 1
    NUMBER_OF_TSTEPS integer 1 num of time steps in one file
    NUMBER_OF_SKIP_TSTEPS integer 0 num of skipped first several data
    FILETYPE_ORG character(len=H_LONG) ''
    BASENAME_ORG character(len=*)
    BASENAME_ADD_NUM logical .false.
    USE_FILE_LANDWATER logical
    INIT_LANDWATER_RATIO real(RP)
    INTRP_LAND_TEMP character(len=*)
    INTRP_LAND_WATER character(len=*)
    INTRP_LAND_SFC_TEMP character(len=*)
    INTRP_ITER_MAX integer
    SOILWATER_DS2VC character(len=H_SHORT) 'limit'
    ELEVATION_COLLECTION logical
    SERIAL_PROC_READ logical .true. read by one MPI process and broadcast

  • PARAM_MKINIT_REAL_OCEAN
    nametypedefault valuecomment
    NUMBER_OF_FILES integer 1
    NUMBER_OF_TSTEPS integer 1 num of time steps in one file
    NUMBER_OF_SKIP_TSTEPS integer 0 num of skipped first several data
    FILETYPE_ORG character(len=H_LONG) ''
    BASENAME_ORG character(len=*)
    BASENAME_ADD_NUM logical .false.
    BASENAME_BOUNDARY character(len=H_LONG) 'boundary_atmos'
    BOUNDARY_TITLE character(len=H_LONG) 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
    BOUNDARY_UPDATE_DT real(RP) 0.0_RP inteval time of boudary data update [s]
    INTRP_OCEAN_TEMP character(len=*)
    INTRP_OCEAN_SFC_TEMP character(len=*)
    INTRP_ITER_MAX integer
    SERIAL_PROC_READ logical .true. read by one MPI process and broadcast

History Output
No history output

Function/Subroutine Documentation

◆ realinput_atmos()

subroutine, public mod_realinput::realinput_atmos ( logical, intent(in)  flg_intrp)

Definition at line 164 of file mod_realinput.f90.

References mod_atmos_vars::dens, scale_grid_index::ia, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::ja, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, mod_atmos_vars::momx, mod_atmos_vars::momy, mod_atmos_vars::momz, parentatomsetup(), scale_process::prc_mpistop(), scale_tracer::qa, mod_atmos_vars::qtrc, mod_atmos_vars::rhot, and scale_tracer::tracer_type.

Referenced by mod_mkinit::interporation_fact().

164  use mod_atmos_vars, only: &
165  dens, &
166  momz, &
167  momx, &
168  momy, &
169  rhot, &
170  qtrc
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_lnml ) write(io_fid_log,nml=param_mkinit_real_atmos)
227 
228  if( tracer_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
real(rp), dimension(:,:,:), allocatable, target, public momz
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:,:), allocatable, target, public rhot
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:), allocatable, target, public momx
integer, public ke
end point of inner domain: z, local
integer, public qa
real(rp), dimension(:,:,:), allocatable, target, public dens
character(len=h_short), public tracer_type
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, target, public momy
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ realinput_surface()

subroutine, public mod_realinput::realinput_surface ( )

Definition at line 354 of file mod_realinput.f90.

References mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_albedo, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_temp, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0e, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0h, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0m, scale_const::const_i_lw, scale_const::const_i_sw, mod_atmos_vars::dens, scale_grid_index::ia, scale_external_io::igrads, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::ja, mod_land_vars::land_sfc_albedo, mod_land_vars::land_sfc_temp, mod_land_vars::land_temp, mod_land_vars::land_water, scale_landuse::landuse_fact_land, scale_landuse::landuse_fact_ocean, scale_landuse::landuse_fact_urban, scale_land_grid_index::lkmax, mod_atmos_vars::momx, mod_atmos_vars::momy, mod_atmos_vars::momz, mod_ocean_vars::ocean_sfc_albedo, mod_ocean_vars::ocean_sfc_temp, mod_ocean_vars::ocean_sfc_z0e, mod_ocean_vars::ocean_sfc_z0h, mod_ocean_vars::ocean_sfc_z0m, mod_ocean_vars::ocean_temp, scale_process::prc_mpistop(), mod_atmos_vars::qtrc, mod_atmos_vars::rhot, scale_urban_grid_index::uke, scale_urban_grid_index::uks, mod_urban_vars::urban_qc, mod_urban_vars::urban_rainb, mod_urban_vars::urban_raing, mod_urban_vars::urban_rainr, mod_urban_vars::urban_roff, mod_urban_vars::urban_sfc_albedo, mod_urban_vars::urban_sfc_temp, mod_urban_vars::urban_tb, mod_urban_vars::urban_tbl, mod_urban_vars::urban_tc, mod_urban_vars::urban_tg, mod_urban_vars::urban_tgl, mod_urban_vars::urban_tr, mod_urban_vars::urban_trl, and mod_urban_vars::urban_uc.

Referenced by mod_mkinit::interporation_fact().

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  character(len=H_SHORT) :: intrp_land_temp = 'off'
407  character(len=H_SHORT) :: intrp_land_water = 'off'
408  character(len=H_SHORT) :: intrp_land_sfc_temp = 'off'
409  character(len=H_SHORT) :: intrp_ocean_temp = 'off'
410  character(len=H_SHORT) :: intrp_ocean_sfc_temp = 'off'
411  integer :: intrp_iter_max = 20
412  character(len=H_SHORT) :: soilwater_ds2vc = 'limit'
413  logical :: soilwater_ds2vc_flag ! true: 'critical', false: 'limit'
414  logical :: elevation_collection = .true.
415 
416  namelist / param_mkinit_real_land / &
417  number_of_files, &
418  number_of_tsteps, &
419  number_of_skip_tsteps, &
420  filetype_org, &
421  basename_org, &
422  basename_add_num, &
423  use_file_landwater, &
424  init_landwater_ratio, &
425  intrp_land_temp, &
426  intrp_land_water, &
427  intrp_land_sfc_temp, &
428  intrp_iter_max, &
429  soilwater_ds2vc, &
430  elevation_collection, &
431  serial_proc_read
432 
433  namelist / param_mkinit_real_ocean / &
434  number_of_files, &
435  number_of_tsteps, &
436  number_of_skip_tsteps, &
437  filetype_org, &
438  basename_org, &
439  basename_add_num, &
440  basename_boundary, &
441  boundary_title, &
442  boundary_update_dt, &
443  intrp_ocean_temp, &
444  intrp_ocean_sfc_temp, &
445  intrp_iter_max, &
446  serial_proc_read
447 
448  character(len=H_LONG) :: filetype_land
449  character(len=H_LONG) :: filetype_ocean
450  character(len=H_LONG) :: basename_land
451  character(len=H_LONG) :: basename_ocean
452  character(len=H_LONG) :: basename_withnum = ''
453  character(len=5) :: num = ''
454  logical :: serial_proc_read_land
455  logical :: serial_proc_read_ocean
456 
457  ! land
458  real(RP) :: land_temp_org(lkmax,ia,ja)
459  real(RP) :: land_water_org(lkmax,ia,ja)
460  real(RP) :: land_sfc_temp_org(ia,ja)
461  real(RP) :: land_sfc_albedo_org(ia,ja,2)
462 
463  ! urban
464  real(RP) :: urban_tc_org(ia,ja)
465  real(RP) :: urban_qc_org(ia,ja)
466  real(RP) :: urban_uc_org(ia,ja)
467  real(RP) :: urban_sfc_temp_org(ia,ja)
468  real(RP) :: urban_sfc_albedo_org(ia,ja,2)
469 
470  ! ocean
471  real(RP), allocatable :: ocean_temp_org(:,:,:)
472  real(RP), allocatable :: ocean_sfc_temp_org(:,:,:)
473  real(RP), allocatable :: ocean_sfc_albedo_org(:,:,:,:)
474  real(RP), allocatable :: ocean_sfc_z0_org(:,:,:)
475 
476  integer :: mdlid_land, mdlid_ocean
477  integer :: ldims(3), odims(2)
478 
479  integer :: totaltimesteps = 1
480  integer :: timelen
481  integer :: skip_steps
482  integer :: lit
483  integer :: lfn
484  integer :: ierr
485 
486  integer :: k, i, j, n, ns, ne
487  !---------------------------------------------------------------------------
488 
489  if( io_l ) write(io_fid_log,*)
490  if( io_l ) write(io_fid_log,*) '+++ Module[RealCaseSurface]/Categ[MKINIT]'
491 
492 
493  ! LAND/URBAN
494 
495  !--- read namelist
496  rewind(io_fid_conf)
497  read(io_fid_conf,nml=param_mkinit_real_land,iostat=ierr)
498  if( ierr < 0 ) then !--- missing
499  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
500  elseif( ierr > 0 ) then !--- fatal error
501  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_LAND. Check!'
502  call prc_mpistop
503  endif
504  if( io_lnml ) write(io_fid_log,nml=param_mkinit_real_land)
505 
506  filetype_land = filetype_org
507 
508  lfn = number_of_skip_tsteps / number_of_tsteps
509  if ( filetype_land .ne. "GrADS" .and. ( number_of_files > 1 .or. basename_add_num ) ) then
510  write(num,'(I5.5)') lfn
511  basename_land = trim(basename_org)//"_"//num
512  else
513  basename_land = trim(basename_org)
514  end if
515 
516  serial_land = serial_proc_read
517 
518  lit = mod(number_of_skip_tsteps,number_of_tsteps)+1
519 
520  !--- read external file
521  if( io_l ) write(io_fid_log,*) ' '
522  if( io_l ) write(io_fid_log,*) '+++ Target File Name (Land): ',trim(basename_land)
523  if( io_l ) write(io_fid_log,*) ' Time Steps in One File: ', number_of_tsteps
524  if( io_l ) write(io_fid_log,*) ' Time Step to read: ', lit
525 
526 
527 
528  ! OCEAN
529 
530  !--- read namelist
531  rewind(io_fid_conf)
532  read(io_fid_conf,nml=param_mkinit_real_ocean,iostat=ierr)
533  if( ierr < 0 ) then !--- missing
534  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
535  elseif( ierr > 0 ) then !--- fatal error
536  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_OCEAN. Check!'
537  call prc_mpistop
538  endif
539  if( io_lnml ) write(io_fid_log,nml=param_mkinit_real_ocean)
540 
541  filetype_ocean = filetype_org
542 
543  if ( filetype_ocean .ne. "GrADS" .and. ( number_of_files > 1 .or. basename_add_num ) ) then
544  basename_ocean = trim(basename_org)//"_00000"
545  else
546  basename_ocean = trim(basename_org)
547  end if
548 
549  select case ( soilwater_ds2vc )
550  case ( 'critical' )
551  soilwater_ds2vc_flag = .true.
552  case ('limit' )
553  soilwater_ds2vc_flag = .false.
554  case default
555  write(*,*) ' xxx Unsupported SOILWATER_DS2CV TYPE:', trim(soilwater_ds2vc)
556  call prc_mpistop
557  end select
558 
559  serial_ocean = serial_proc_read
560 
561 
562  call parentsurfacesetup( ldims, odims, & ![OUT]
563  mdlid_land, & ![OUT]
564  mdlid_ocean, & ![OUT]
565  timelen, & ![OUT]
566  basename_land, & ![IN]
567  basename_ocean, & ![IN]
568  filetype_land, & ![IN]
569  filetype_ocean, & ![IN]
570  use_file_landwater, & ![IN]
571  intrp_land_temp, & ![IN]
572  intrp_land_water, & ![IN]
573  intrp_land_sfc_temp, & ![IN]
574  intrp_ocean_temp, & ![IN]
575  intrp_ocean_sfc_temp ) ![IN]
576 
577  if ( timelen > 0 ) then
578  number_of_tsteps = timelen ! read from file
579  endif
580 
581  totaltimesteps = number_of_files * number_of_tsteps
582 
583  allocate( ocean_temp_org(ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
584  allocate( ocean_sfc_temp_org(ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
585  allocate( ocean_sfc_albedo_org(ia,ja,2,1+number_of_skip_tsteps:totaltimesteps ) )
586  allocate( ocean_sfc_z0_org(ia,ja,1+number_of_skip_tsteps:totaltimesteps ) )
587 
588  if ( mdlid_land == igrads .and. ( number_of_files > 1 .or. basename_add_num ) ) then
589  write(num,'(I5.5)') lfn
590  basename_land = "_"//num
591  end if
592 
593  if ( mdlid_ocean == igrads ) then
594  basename_org = ""
595  end if
596 
597  !--- read external file
598  do n = 1, number_of_files
599 
600  if ( number_of_files > 1 .or. basename_add_num ) then
601  write(num,'(I5.5)') n-1
602  basename_ocean = trim(basename_org)//"_"//num
603  else
604  basename_ocean = trim(basename_org)
605  end if
606 
607  if( io_l ) write(io_fid_log,*) ' '
608  if( io_l ) write(io_fid_log,*) '+++ Target File Name (Ocean): ', trim(basename_ocean)
609  if( io_l ) write(io_fid_log,*) ' Time Steps in One File: ', number_of_tsteps
610 
611  ns = number_of_tsteps * (n - 1) + 1
612  ne = ns + (number_of_tsteps - 1)
613 
614  if ( ne <= number_of_skip_tsteps ) then
615  if( io_l ) write(io_fid_log,*) ' SKIP'
616  cycle
617  end if
618 
619  skip_steps = max(number_of_skip_tsteps - ns + 1, 0)
620  ns = max(ns, number_of_skip_tsteps+1)
621 
622  ! read all prepared data
623  call parentsurfaceinput( land_temp_org, &
624  land_water_org, &
625  land_sfc_temp_org, &
626  land_sfc_albedo_org, &
627  urban_tc_org, &
628  urban_qc_org, &
629  urban_uc_org, &
630  urban_sfc_temp_org, &
631  urban_sfc_albedo_org, &
632  ocean_temp_org(:,:,ns:ne), &
633  ocean_sfc_temp_org(:,:,ns:ne), &
634  ocean_sfc_albedo_org(:,:,:,ns:ne), &
635  ocean_sfc_z0_org(:,:,ns:ne), &
636  dens, &
637  momz, &
638  momx, &
639  momy, &
640  rhot, &
641  qtrc, &
642  basename_land, &
643  basename_ocean, &
644  mdlid_land, mdlid_ocean, &
645  ldims, odims, &
646  use_file_landwater, &
647  init_landwater_ratio, &
648  intrp_iter_max, &
649  soilwater_ds2vc_flag, &
650  elevation_collection, &
651  number_of_tsteps, &
652  skip_steps, lit )
653 
654  enddo
655 
656 
657  !--- input initial data
658  do j = 1, ja
659  do i = 1, ia
660  land_sfc_temp(i,j) = land_sfc_temp_org(i,j)
661  land_sfc_albedo(i,j,i_lw) = land_sfc_albedo_org(i,j,i_lw)
662  land_sfc_albedo(i,j,i_sw) = land_sfc_albedo_org(i,j,i_sw)
663  do k = 1, lkmax
664  land_temp(k,i,j) = land_temp_org(k,i,j)
665  land_water(k,i,j) = land_water_org(k,i,j)
666  end do
667 
668  urban_sfc_temp(i,j) = urban_sfc_temp_org(i,j)
669  urban_sfc_albedo(i,j,i_lw) = urban_sfc_albedo_org(i,j,i_lw)
670  urban_sfc_albedo(i,j,i_sw) = urban_sfc_albedo_org(i,j,i_sw)
671  do k = uks, uke
672  urban_trl(k,i,j) = urban_sfc_temp_org(i,j)
673  urban_tbl(k,i,j) = urban_sfc_temp_org(i,j)
674  urban_tgl(k,i,j) = urban_sfc_temp_org(i,j)
675  end do
676  urban_tc(i,j) = urban_tc_org(i,j)
677  urban_qc(i,j) = urban_qc_org(i,j)
678  urban_uc(i,j) = urban_uc_org(i,j)
679  urban_tr(i,j) = urban_sfc_temp_org(i,j)
680  urban_tb(i,j) = urban_sfc_temp_org(i,j)
681  urban_tg(i,j) = urban_sfc_temp_org(i,j)
682  urban_rainr(i,j) = 0.0_rp
683  urban_rainb(i,j) = 0.0_rp
684  urban_raing(i,j) = 0.0_rp
685  urban_roff(i,j) = 0.0_rp
686  enddo
687  enddo
688 
689 
690  ns = number_of_skip_tsteps + 1 ! skip first several data
691  do j = 1, ja
692  do i = 1, ia
693  ocean_temp(i,j) = ocean_temp_org(i,j,ns)
694  ocean_sfc_temp(i,j) = ocean_sfc_temp_org(i,j,ns)
695  ocean_sfc_albedo(i,j,i_lw) = ocean_sfc_albedo_org(i,j,i_lw,ns)
696  ocean_sfc_albedo(i,j,i_sw) = ocean_sfc_albedo_org(i,j,i_sw,ns)
697  ocean_sfc_z0m(i,j) = ocean_sfc_z0_org(i,j,ns)
698  ocean_sfc_z0h(i,j) = ocean_sfc_z0_org(i,j,ns)
699  ocean_sfc_z0e(i,j) = ocean_sfc_z0_org(i,j,ns)
700  enddo
701  enddo
702 
703  do j = 1, ja
704  do i = 1, ia
705  atmos_phy_sf_sfc_temp(i,j) = fact_ocean(i,j) * ocean_sfc_temp(i,j) &
706  + fact_land(i,j) * land_sfc_temp(i,j) &
707  + fact_urban(i,j) * urban_sfc_temp(i,j)
708  atmos_phy_sf_sfc_albedo(i,j,i_lw) = fact_ocean(i,j) * ocean_sfc_albedo(i,j,i_lw) &
709  + fact_land(i,j) * land_sfc_albedo(i,j,i_lw) &
710  + fact_urban(i,j) * urban_sfc_albedo(i,j,i_lw)
711  atmos_phy_sf_sfc_albedo(i,j,i_sw) = fact_ocean(i,j) * ocean_sfc_albedo(i,j,i_sw) &
712  + fact_land(i,j) * land_sfc_albedo(i,j,i_sw) &
713  + fact_urban(i,j) * urban_sfc_albedo(i,j,i_sw)
717  end do
718  end do
719 
720 
721  !--- output boundary data
722  totaltimesteps = totaltimesteps - number_of_skip_tsteps ! skip first several data
723  if ( totaltimesteps > 1 ) then
724  if ( boundary_update_dt <= 0.0_rp ) then
725  write(*,*) 'xxx BOUNDARY_UPDATE_DT is necessary in real case preprocess'
726  call prc_mpistop
727  endif
728 
729  call parentoceanboundary( ocean_temp_org(:,:,ns:ne), &
730  ocean_sfc_temp_org(:,:,ns:ne), &
731  ocean_sfc_albedo_org(:,:,:,ns:ne), &
732  ocean_sfc_z0_org(:,:,ns:ne), &
733  totaltimesteps, &
734  boundary_update_dt, &
735  basename_boundary, &
736  boundary_title )
737 
738  end if
739 
740  deallocate( ocean_temp_org )
741  deallocate( ocean_sfc_temp_org )
742  deallocate( ocean_sfc_albedo_org )
743  deallocate( ocean_sfc_z0_org )
744 
745  return
real(rp), dimension(:,:,:), allocatable, target, public momz
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
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:98
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public urban_qc
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
real(rp), dimension(:,:,:), allocatable, target, public rhot
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0e
ocean surface roughness length for vapor [m]
real(rp), dimension(:,:), allocatable, public ocean_temp
temperature at uppermost ocean layer [K]
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:), allocatable, target, public momx
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
real(rp), dimension(:,:), allocatable, public urban_tb
real(rp), dimension(:,:,:), allocatable, target, public dens
module URBAN Variables
real(rp), dimension(:,:), allocatable, public urban_raing
real(rp), dimension(:,:), allocatable, public urban_uc
integer, parameter, public i_lw
real(rp), dimension(:,:), allocatable, public land_sfc_temp
land surface skin temperature [K]
module ATMOSPHERIC Surface Variables
integer, parameter, public i_sw
real(rp), dimension(:,:), allocatable, public urban_tr
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:,:,:), allocatable, public urban_tgl
module LANDUSE
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
real(rp), dimension(:,:,:), allocatable, public land_temp
temperature of each soil layer [K]
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
real(rp), dimension(:,:), allocatable, public urban_roff
real(rp), dimension(:,:,:), allocatable, public land_sfc_albedo
land surface albedo [0-1]
module LAND Variables
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(:,:,:), allocatable, public urban_sfc_albedo
real(rp), dimension(:,:), allocatable, public urban_tc
real(rp), dimension(:,:,:), allocatable, target, public momy
real(rp), dimension(:,:), allocatable, public urban_rainr
real(rp), dimension(:,:,:), allocatable, public land_water
moisture of each soil layer [m3/m3]
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0m
ocean surface roughness length for momentum [m]
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
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
real(rp), dimension(:,:,:), allocatable, public urban_trl
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), dimension(:,:,:), allocatable, public urban_tbl
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0h
ocean surface roughness length for heat [m]
module OCEAN Variables
real(rp), dimension(:,:), allocatable, public urban_rainb
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentatomsetup()

subroutine mod_realinput::parentatomsetup ( integer, dimension(6), intent(out)  dims,
integer, intent(out)  timelen,
integer, intent(out)  mdlid,
character(len=*), intent(in)  basename_org,
character(len=*), intent(in)  filetype,
logical, intent(in)  use_file_density_in,
logical, intent(in)  serial_in 
)

Atmos Setup.

Definition at line 759 of file mod_realinput.f90.

References scale_tracer::aq_name, scale_tracer::aq_unit, scale_const::const_eps, scale_const::const_i_lw, scale_const::const_i_sw, scale_const::const_undef, scale_land_grid::grid_lcz, scale_gridtrans::gtrans_rotc, scale_tracer::i_qc, scale_tracer::i_qv, scale_grid_index::ia, scale_external_io::igrads, scale_external_io::inicam, scale_interpolation_nest::intrpnest_domain_compatibility(), scale_interpolation_nest::intrpnest_interp_2d, scale_interpolation_nest::intrpnest_interp_3d, scale_interpolation_nest::intrpnest_interp_fact_latlon(), scale_interpolation_nest::intrpnest_interp_fact_llz(), scale_stdio::io_fid_log, scale_stdio::io_l, scale_external_io::iscale, scale_external_io::iwrfarw, scale_grid_index::ja, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, land_interporation(), scale_landuse::landuse_frac_land, scale_land_grid_index::lkmax, dc_log::log(), make_mask(), scale_grid_nest::nest_interp_level, mod_realinput_grads::parentatominputgrads(), mod_realinput_nicam::parentatominputnicam(), mod_realinput_scale::parentatominputscale(), mod_realinput_wrfarw::parentatominputwrfarw(), mod_realinput_grads::parentatomopengrads(), mod_realinput_nicam::parentatomopennicam(), mod_realinput_scale::parentatomopenscale(), mod_realinput_wrfarw::parentatomopenwrfarw(), mod_realinput_grads::parentatomsetupgrads(), mod_realinput_nicam::parentatomsetupnicam(), mod_realinput_scale::parentatomsetupscale(), mod_realinput_wrfarw::parentatomsetupwrfarw(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_nicam::parentlandinputnicam(), mod_realinput_scale::parentlandinputscale(), mod_realinput_wrfarw::parentlandinputwrfarw(), mod_realinput_grads::parentlandsetupgrads(), mod_realinput_nicam::parentlandsetupnicam(), mod_realinput_scale::parentlandsetupscale(), mod_realinput_wrfarw::parentlandsetupwrfarw(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_nicam::parentoceaninputnicam(), mod_realinput_scale::parentoceaninputscale(), mod_realinput_wrfarw::parentoceaninputwrfarw(), mod_realinput_grads::parentoceanopengrads(), mod_realinput_nicam::parentoceanopennicam(), mod_realinput_scale::parentoceanopenscale(), mod_realinput_wrfarw::parentoceanopenwrfarw(), mod_realinput_grads::parentoceansetupgrads(), mod_realinput_nicam::parentoceansetupnicam(), mod_realinput_scale::parentoceansetupscale(), mod_realinput_wrfarw::parentoceansetupwrfarw(), scale_process::prc_ismaster, scale_process::prc_mpistop(), scale_tracer::qa, scale_tracer::qwe, scale_tracer::qws, replace_misval_map(), scale_time::time_nowdate, and scale_tracer::tracer_type.

Referenced by realinput_atmos().

759  use scale_external_io, only: &
760  iscale, &
761  iwrfarw, &
762  inicam, &
763  igrads
764  use mod_realinput_scale, only: &
766  use mod_realinput_wrfarw, only: &
768  use mod_realinput_nicam, only: &
770  use mod_realinput_grads, only: &
772  implicit none
773 
774  integer, intent(out) :: dims(6)
775  integer, intent(out) :: timelen
776  integer, intent(out) :: mdlid
777  character(len=*), intent(in) :: basename_org
778  character(len=*), intent(in) :: filetype
779  logical, intent(in) :: serial_in ! read by a serial process
780  logical, intent(in) :: use_file_density_in ! use density data from files
781 
782  !---------------------------------------------------------------------------
783 
784  if( io_l ) write(io_fid_log,*)
785  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[RealinputAtmos]/Categ[Setup]'
786 
787  serial = serial_in
788  if( serial ) then
789  if( prc_ismaster ) then
790  do_read_atom = .true.
791  else
792  do_read_atom = .false.
793  endif
794  else
795  do_read_atom = .true.
796  endif
797 
798  select case(trim(filetype))
799  case('SCALE-RM')
800 
801  mdlid = iscale
802  do_read_atom = .true.
803  serial = .false.
804  call parentatomsetupscale( dims ) ! (out)
805  update_coord = .false.
806  use_file_density = use_file_density_in
807  use_temp = .false.
808  rotate = .false.
809  timelen = -1
810 
811  case('WRF-ARW')
812 
813  mdlid = iwrfarw
814  if ( do_read_atom ) call parentatomsetupwrfarw( dims, timelen, & ! (out)
815  basename_org ) ! (in)
816  update_coord = .true.
817  use_file_density = .false.
818  use_temp = .true.
819  rotate = .true.
820 
821  case('NICAM-NETCDF')
822 
823  mdlid = inicam
824  if ( do_read_atom ) call parentatomsetupnicam( dims, timelen, & ! (out)
825  basename_org ) ! (in)
826  update_coord = .false.
827  use_file_density = .false.
828  use_temp = .true.
829  rotate = .true.
830 
831  case('GrADS')
832 
833  mdlid = igrads
834  if ( do_read_atom ) call parentatomsetupgrads( dims, & ! (out)
835  basename_org ) ! (in)
836  update_coord = .true.
837  use_file_density = .false.
838  use_temp = .true.
839  rotate = .true.
840  timelen = -1
841 
842  case default
843 
844  write(*,*) ' xxx Unsupported FILE TYPE:', trim(filetype)
845  call prc_mpistop
846 
847  endselect
848 
849  if( serial ) then
850  call comm_bcast( dims(:), 6 )
851  call comm_bcast( timelen )
852  endif
853 
854  if( io_l ) write(io_fid_log,*) '+++ Horizontal Interpolation Level:', &
856  itp_nh = int( nest_interp_level )
857  itp_nv = 2
858 
859  allocate( hfact( ia, ja, itp_nh ) )
860  allocate( vfact( ka, ia, ja, itp_nh, itp_nv ) )
861  allocate( igrd( ia, ja, itp_nh ) )
862  allocate( jgrd( ia, ja, itp_nh ) )
863  allocate( kgrd( ka, ia, ja, itp_nh, itp_nv ) )
864  allocate( ncopy( ia, ja, itp_nh ) )
865 
866  allocate( lon_org( dims(2), dims(3) ) )
867  allocate( lat_org( dims(2), dims(3) ) )
868  allocate( cz_org( dims(1)+2, dims(2), dims(3) ) )
869 
870  allocate( velz_org( dims(1)+2, dims(2), dims(3) ) )
871  allocate( velx_org( dims(1)+2, dims(2), dims(3) ) )
872  allocate( vely_org( dims(1)+2, dims(2), dims(3) ) )
873  allocate( pott_org( dims(1)+2, dims(2), dims(3) ) )
874  allocate( temp_org( dims(1)+2, dims(2), dims(3) ) )
875  allocate( pres_org( dims(1)+2, dims(2), dims(3) ) )
876  allocate( qtrc_org( dims(1)+2, dims(2), dims(3), qa ) )
877  allocate( dens_org( dims(1)+2, dims(2), dims(3) ) )
878 
879  return
subroutine, public parentatomsetupgrads(dims, basename)
Atmos Setup.
subroutine, public prc_mpistop
Abort MPI.
module REAL input WRF-ARW
integer, parameter, public inicam
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public qa
integer, parameter, public igrads
integer, public nest_interp_level
horizontal interpolation level
integer, parameter, public iwrfarw
integer, parameter, public iscale
module FILE I/O (netcdf)
module REAL input NICAM
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
subroutine, public parentatomsetupscale(dims)
Atmos Setup.
subroutine, public parentatomsetupnicam(dims, timelen, basename_org)
Atmos Setup.
module REAL input SCALE
subroutine, public parentatomsetupwrfarw(dims, timelen, basename_org)
Atmos Setup.
module REAL input GrADS
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ land_interporation()

subroutine mod_realinput::land_interporation ( real(rp), dimension(lkmax,ia,ja), intent(out)  tg,
real(rp), dimension(lkmax,ia,ja), intent(out)  strg,
real(rp), dimension(ia,ja), intent(out)  lst,
real(rp), dimension(ia,ja,2), intent(out)  albg,
real(rp), dimension(ia,ja), intent(out)  ust,
real(rp), dimension(ia,ja,2), intent(out)  albu,
real(rp), dimension(:,:,:), intent(inout)  tg_org,
real(rp), dimension(:,:,:), intent(inout)  strg_org,
real(rp), dimension(:,:,:), intent(inout)  smds_org,
real(rp), dimension(:,:), intent(inout)  lst_org,
real(rp), dimension(:,:,:), intent(inout)  albg_org,
real(rp), dimension(:,:), intent(inout)  ust_org,
real(rp), dimension(:,:), intent(inout)  sst_org,
real(rp), dimension(:,:), intent(in)  lmask_org,
real(rp), dimension(:,:), intent(in)  lsmask_nest,
real(rp), dimension(:,:), intent(in)  topo_org,
real(rp), dimension(:), intent(in)  lz_org,
real(rp), dimension(:,:), intent(in)  llon_org,
real(rp), dimension(:,:), intent(in)  llat_org,
real(rp), dimension(lkmax), intent(in)  LCZ,
real(rp), dimension(ia,ja), intent(in)  LON,
real(rp), dimension(ia,ja), intent(in)  LAT,
integer, dimension(3), intent(in)  ldims,
integer, dimension(2), intent(in)  odims,
real(rp), intent(in)  maskval_tg,
real(rp), intent(in)  maskval_strg,
real(rp), intent(in)  init_landwater_ratio,
logical, intent(in)  use_file_landwater,
logical, intent(in)  use_waterratio,
logical, intent(in)  soilwater_ds2vc_flag,
logical, intent(in)  elevation_collection,
integer, intent(in)  intrp_iter_max 
)

Definition at line 2301 of file mod_realinput.f90.

References scale_const::const_i_lw, scale_const::const_i_sw, scale_const::const_laps, scale_const::const_undef, mod_land_vars::convert_ws2vwc(), scale_grid_index::ia, scale_interpolation_nest::intrpnest_interp_2d, scale_interpolation_nest::intrpnest_interp_3d, scale_interpolation_nest::intrpnest_interp_fact_latlon(), scale_interpolation_nest::intrpnest_interp_fact_llz(), scale_grid_index::ja, scale_land_grid_index::lkmax, make_mask(), scale_process::prc_mpistop(), replace_misval_const(), replace_misval_map(), and scale_topography::topo_zsfc.

Referenced by parentatomsetup().

2301  use scale_process, only: &
2302  prc_mpistop
2303  use scale_const, only: &
2304  undef => const_undef, &
2305  i_sw => const_i_sw, &
2306  i_lw => const_i_lw, &
2307  laps => const_laps
2308  use scale_interpolation_nest, only: &
2313  use scale_topography, only: &
2314  topo_zsfc
2315  use mod_land_vars, only: &
2317  implicit none
2318  real(RP), intent(out) :: tg(lkmax,ia,ja)
2319  real(RP), intent(out) :: strg(lkmax,ia,ja)
2320  real(RP), intent(out) :: lst(ia,ja)
2321  real(RP), intent(out) :: albg(ia,ja,2)
2322  real(RP), intent(out) :: ust(ia,ja)
2323  real(RP), intent(out) :: albu(ia,ja,2)
2324  real(RP), intent(inout) :: tg_org(:,:,:)
2325  real(RP), intent(inout) :: strg_org(:,:,:)
2326  real(RP), intent(inout) :: smds_org(:,:,:)
2327  real(RP), intent(inout) :: lst_org(:,:)
2328  real(RP), intent(inout) :: albg_org(:,:,:)
2329  real(RP), intent(inout) :: ust_org(:,:)
2330  real(RP), intent(inout) :: sst_org(:,:)
2331  real(RP), intent(in) :: lmask_org(:,:)
2332  real(RP), intent(in) :: lsmask_nest(:,:)
2333  real(RP), intent(in) :: topo_org(:,:)
2334  real(RP), intent(in) :: lz_org(:)
2335  real(RP), intent(in) :: llon_org(:,:)
2336  real(RP), intent(in) :: llat_org(:,:)
2337  real(RP), intent(in) :: lcz(lkmax)
2338  real(RP), intent(in) :: lon(ia,ja)
2339  real(RP), intent(in) :: lat(ia,ja)
2340  integer, intent(in) :: ldims(3)
2341  integer, intent(in) :: odims(2)
2342  real(RP), intent(in) :: maskval_tg
2343  real(RP), intent(in) :: maskval_strg
2344  real(RP), intent(in) :: init_landwater_ratio
2345  logical, intent(in) :: use_file_landwater
2346  logical, intent(in) :: use_waterratio
2347  logical, intent(in) :: soilwater_ds2vc_flag
2348  logical, intent(in) :: elevation_collection
2349  integer, intent(in) :: intrp_iter_max
2350 
2351  real(RP) :: lmask(ldims(2), ldims(3))
2352  real(RP) :: smds(lkmax,ia,ja)
2353 
2354  ! data for interporation
2355  real(RP) :: hfact_l(ldims(2), ldims(3), itp_nh)
2356  integer :: igrd_l (ldims(2), ldims(3), itp_nh)
2357  integer :: jgrd_l (ldims(2), ldims(3), itp_nh)
2358  real(RP) :: vfactl(lkmax,ia,ja,itp_nh,itp_nv)
2359  integer :: kgrdl (lkmax,ia,ja,itp_nh,itp_nv)
2360 
2361  real(RP) :: sst_land(ldims(2), ldims(3))
2362  real(RP) :: work(ldims(2), ldims(3))
2363 
2364  real(RP) :: lz3d_org(ldims(1),ldims(2),ldims(3))
2365  real(RP) :: lcz_3d(lkmax,ia,ja)
2366 
2367  ! elevation collection
2368  real(RP) :: topo(ia,ja)
2369  real(RP) :: tdiff
2370 
2371  integer :: k, i, j
2372 
2373 
2374  ! Surface skin temp: interpolate over the ocean
2375  if ( i_intrp_land_sfc_temp .ne. i_intrp_off ) then
2376  select case ( i_intrp_land_sfc_temp )
2377  case ( i_intrp_mask )
2378  lmask = lmask_org
2379  case ( i_intrp_fill )
2380  call make_mask( lmask, lst_org, ldims(2), ldims(3), landdata=.true.)
2381  case default
2382  write(*,*) 'xxx INTRP_LAND_SFC_TEMP is invalid.'
2383  call prc_mpistop
2384  end select
2385  call interp_oceanland_data(lst_org, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2386  end if
2387 
2388  ! Urban surface temp: interpolate over the ocean
2389  ! if ( i_INTRP_URB_SFC_TEMP .ne. i_intrp_off ) then
2390  ! select case ( i_INTRP_URB_SFC_TEMP )
2391  ! case ( i_intrp_mask )
2392  ! lmask = lmask_org
2393  ! case ( i_intrp_fill )
2394  ! call make_mask( lmask, ust_org, ldims(2), ldims(3), landdata=.true.)
2395  ! case default
2396  ! write(*,*) 'xxx INTRP_URB_SFC_TEMP is invalid.'
2397  ! call PRC_MPIstop
2398  ! end select
2399  ! call interp_OceanLand_data(ust_org, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2400  !end if
2401 
2402  ! interpolation facter between outer land grid and ocean grid
2403  call intrpnest_interp_fact_latlon( hfact_l(:,:,:), & ! [OUT]
2404  igrd_l(:,:,:), jgrd_l(:,:,:), & ! [OUT]
2405  llat_org(:,:), llon_org(:,:), & ! [IN]
2406  ldims(2), ldims(3), & ! [IN]
2407  olat_org(:,:), olon_org(:,:), & ! [IN]
2408  odims(1), odims(2) ) ! [IN]
2409 
2410 
2411  ! sst on land grid
2412  call intrpnest_interp_2d( sst_land(:,:), sst_org(:,:), hfact_l(:,:,:), &
2413  igrd_l(:,:,:), jgrd_l(:,:,:), ldims(2), ldims(3) )
2414  call replace_misval_map( lst_org, sst_land, ldims(2), ldims(3), "SKINT")
2415 
2416  ! replace missing value
2417  do j = 1, ldims(3)
2418  do i = 1, ldims(2)
2419  if ( ust_org(i,j) == undef ) ust_org(i,j) = lst_org(i,j)
2420 ! if ( skinw_org(i,j) == UNDEF ) skinw_org(i,j) = 0.0_RP
2421 ! if ( snowq_org(i,j) == UNDEF ) snowq_org(i,j) = 0.0_RP
2422 ! if ( snowt_org(i,j) == UNDEF ) snowt_org(i,j) = TEM00
2423  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
2424  if ( albg_org(i,j,i_sw) == undef ) albg_org(i,j,i_sw) = 0.22_rp
2425  end do
2426  end do
2427 
2428  ! Land temp: interpolate over the ocean
2429  if ( i_intrp_land_temp .ne. i_intrp_off ) then
2430  do k = 1, ldims(1)
2431  work(:,:) = tg_org(k,:,:)
2432  select case ( i_intrp_land_temp )
2433  case ( i_intrp_mask )
2434  lmask = lmask_org
2435  case ( i_intrp_fill )
2436  call make_mask( lmask, work, ldims(2), ldims(3), landdata=.true.)
2437  end select
2438  call interp_oceanland_data( work, lmask, ldims(2), ldims(3), .true., intrp_iter_max )
2439  !replace land temp using skin temp
2440  call replace_misval_map( work, lst_org, ldims(2), ldims(3), "STEMP")
2441  tg_org(k,:,:) = work(:,:)
2442  end do
2443  end if
2444 
2445 
2446  ! fill grid data
2447  do j = 1, ldims(3)
2448  do i = 1, ldims(2)
2449  lz3d_org(:,i,j) = lz_org(:)
2450  end do
2451  end do
2452 
2453  do j = 1, ja
2454  do i = 1, ia
2455  lcz_3d(:,i,j) = lcz(:)
2456  enddo
2457  enddo
2458 
2459  call intrpnest_interp_fact_llz( hfact(:,:,:), & ! [OUT]
2460  vfactl(:,:,:,:,:), & ! [OUT]
2461  kgrdl(:,:,:,:,:), & ! [OUT]
2462  igrd(:,:,:), jgrd(:,:,:), & ! [OUT]
2463  ncopy(:,:,:), & ! [OUT]
2464  lcz_3d(:,:,:), & ! [IN]
2465  lat(:,:), lon(:,:), & ! [IN]
2466  1, lkmax, ia, ja, & ! [IN]
2467  lz3d_org(:,:,:), & ! [IN]
2468  llat_org(:,:), llon_org(:,:), & ! [IN]
2469  ldims(1), ldims(2), ldims(3), & ! [IN]
2470  landgrid=.true. ) ! [IN]
2471 
2472  call intrpnest_interp_2d( lst(:,:), lst_org(:,:), hfact(:,:,:), &
2473  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2474  call intrpnest_interp_2d( ust(:,:), ust_org(:,:), hfact(:,:,:), &
2475  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2476 ! call INTRPNEST_interp_2d( skinw(:,:), skinw_org(:,:), hfact(:,:,:), &
2477 ! igrd(:,:,:), jgrd(:,:,:), IA, JA )
2478 ! call INTRPNEST_interp_2d( snowq(:,:), snowq_org(:,:), hfact(:,:,:), &
2479 ! igrd(:,:,:), jgrd(:,:,:), IA, JA )
2480 ! call INTRPNEST_interp_2d( snowt(:,:), snowt_org(:,:), hfact(:,:,:), &
2481 ! igrd(:,:,:), jgrd(:,:,:), IA, JA )
2482  call intrpnest_interp_2d( albg(:,:,i_lw), albg_org(:,:,i_lw), hfact(:,:,:), &
2483  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2484  call intrpnest_interp_2d( albg(:,:,i_sw), albg_org(:,:,i_sw), hfact(:,:,:), &
2485  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2486 
2487  call intrpnest_interp_3d( tg(:,:,:), &
2488  tg_org(:,:,:), &
2489  hfact(:,:,:), &
2490  vfactl(:,:,:,:,:), &
2491  kgrdl(:,:,:,:,:), &
2492  igrd(:,:,:), &
2493  jgrd(:,:,:), &
2494  ia, ja, 1, lkmax-1 )
2495 
2496  do j = 1, ja
2497  do i = 1, ia
2498  tg(lkmax,i,j) = tg(lkmax-1,i,j)
2499  enddo ! i
2500  enddo ! j
2501 
2502  ! replace values over the ocean
2503  do k = 1, lkmax
2504  call replace_misval_const( tg(k,:,:), maskval_tg, lsmask_nest )
2505  enddo
2506 
2507 
2508  ! elevation collection
2509  if ( elevation_collection ) then
2510  call intrpnest_interp_2d( topo(:,:), topo_org(:,:), hfact(:,:,:), &
2511  igrd(:,:,:), jgrd(:,:,:), ia, ja )
2512 
2513  do j = 1, ja
2514  do i = 1, ia
2515  if ( topo(i,j) > 0.0_rp ) then ! ignore UNDEF value
2516  tdiff = ( topo_zsfc(i,j) - topo(i,j) ) * laps
2517  lst(i,j) = lst(i,j) - tdiff
2518  ust(i,j) = ust(i,j) - tdiff
2519  do k = 1, lkmax
2520  tg(k,i,j) = tg(k,i,j) - tdiff
2521  end do
2522  end if
2523  end do
2524  end do
2525  end if
2526 
2527 
2528 
2529  ! Land water: interpolate over the ocean
2530  if( use_file_landwater )then
2531 
2532  if ( use_waterratio ) then
2533 
2534  if ( i_intrp_land_water .ne. i_intrp_off ) then
2535  do k = 1, ldims(1)
2536  work(:,:) = smds_org(k,:,:)
2537  select case ( i_intrp_land_water )
2538  case ( i_intrp_mask )
2539  lmask = lmask_org
2540  case ( i_intrp_fill )
2541  call make_mask( lmask, work, ldims(2), ldims(3), landdata=.true.)
2542  end select
2543  call interp_oceanland_data(work, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2544  lmask(:,:) = init_landwater_ratio
2545  !replace missing value to init_landwater_ratio
2546  call replace_misval_map( work, lmask, ldims(2), ldims(3), "SMOISDS")
2547  smds_org(k,:,:) = work(:,:)
2548  enddo
2549  end if
2550 
2551  call intrpnest_interp_3d( smds(:,:,:), &
2552  smds_org(:,:,:), &
2553  hfact(:,:,:), &
2554  vfactl(:,:,:,:,:), &
2555  kgrdl(:,:,:,:,:), &
2556  igrd(:,:,:), &
2557  jgrd(:,:,:), &
2558  ia, ja, 1, lkmax-1 )
2559  do k = 1, lkmax
2560  strg(k,:,:) = convert_ws2vwc( smds(k,:,:), critical=soilwater_ds2vc_flag )
2561  end do
2562 
2563  else
2564 
2565  if ( i_intrp_land_water .ne. i_intrp_off ) then
2566  do k = 1, ldims(1)
2567  work(:,:) = strg_org(k,:,:)
2568  select case ( i_intrp_land_water )
2569  case ( i_intrp_mask )
2570  lmask = lmask_org
2571  case ( i_intrp_fill )
2572  call make_mask( lmask, work, ldims(2), ldims(3), landdata=.true.)
2573  end select
2574  call interp_oceanland_data(work, lmask, ldims(2), ldims(3), .true., intrp_iter_max)
2575  lmask(:,:) = maskval_strg
2576  !replace missing value to init_landwater_ratio
2577  call replace_misval_map( work, lmask, ldims(2), ldims(3), "SMOIS")
2578  strg_org(k,:,:) = work(:,:)
2579  enddo
2580  end if
2581 
2582  call intrpnest_interp_3d( strg(:,:,:), &
2583  strg_org(:,:,:), &
2584  hfact(:,:,:), &
2585  vfactl(:,:,:,:,:), &
2586  kgrdl(:,:,:,:,:), &
2587  igrd(:,:,:), &
2588  jgrd(:,:,:), &
2589  ia, ja, 1, lkmax-1 )
2590  ! interpolation
2591  do j = 1, ja
2592  do i = 1, ia
2593  strg(lkmax,i,j) = strg(lkmax-1,i,j)
2594  enddo
2595  enddo
2596 
2597  end if
2598 
2599  ! replace values over the ocean
2600  do k = 1, lkmax
2601  call replace_misval_const( strg(k,:,:), maskval_strg, lsmask_nest )
2602  enddo
2603 
2604  else ! not read from boundary file
2605 
2606  smds(:,:,:) = init_landwater_ratio
2607  ! conversion from water saturation [fraction] to volumetric water content [m3/m3]
2608  do k = 1, lkmax
2609  strg(k,:,:) = convert_ws2vwc( smds(k,:,:), critical=.true. )
2610  end do
2611 
2612  endif ! use_file_waterratio
2613 
2614 
2615  ! copy albedo of land to urban
2616  do j = 1, ja
2617  do i = 1, ia
2618  albu(i,j,:) = albg(i,j,:)
2619  enddo
2620  enddo
2621 
2622 
2623  return
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
subroutine, public prc_mpistop
Abort MPI.
procedure(intrpnest_intfc_interp_2d), pointer, public intrpnest_interp_2d
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:60
subroutine, public intrpnest_interp_fact_latlon(hfact, igrd, jgrd, mylat, mylon, myIA, myJA, inlat, inlon, inIA, inJA)
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)
integer, parameter, public i_lw
real(rp), public const_undef
Definition: scale_const.F90:43
integer, parameter, public i_sw
integer, public ia
of x whole cells (local, with HALO)
real(rp) function, dimension(ia, ja), public convert_ws2vwc(WS, critical)
conversion from water saturation [fraction] to volumetric water content [m3/m3]
module PROCESS
procedure(intrpnest_intfc_interp_3d), pointer, public intrpnest_interp_3d
module LAND Variables
module CONSTANT
Definition: scale_const.F90:14
module INTERPOLATION (nesting system)
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
module TOPOGRAPHY
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ make_mask()

subroutine mod_realinput::make_mask ( real(rp), dimension(:,:), intent(out)  gmask,
real(rp), dimension(:,:), intent(in)  data,
integer, intent(in)  nx,
integer, intent(in)  ny,
logical, intent(in)  landdata 
)

Definition at line 2633 of file mod_realinput.f90.

References scale_const::const_eps, and scale_const::const_undef.

Referenced by land_interporation(), and parentatomsetup().

2633  use scale_const, only: &
2634  eps => const_eps, &
2635  undef => const_undef
2636  implicit none
2637  real(RP), intent(out) :: gmask(:,:)
2638  real(RP), intent(in) :: data(:,:)
2639  integer, intent(in) :: nx
2640  integer, intent(in) :: ny
2641  logical, intent(in) :: landdata ! .true. => land data , .false. => ocean data
2642 
2643  real(RP) :: dd
2644  integer :: i,j
2645 
2646  if( landdata )then
2647  gmask(:,:) = 1.0_rp ! gmask=1 will be skip in "interp_OceanLand_data"
2648  dd = 0.0_rp
2649  else
2650  gmask(:,:) = 0.0_rp ! gmask=0 will be skip in "interp_OceanLand_data"
2651  dd = 1.0_rp
2652  endif
2653 
2654  do j = 1, ny
2655  do i = 1, nx
2656  if( abs(data(i,j) - undef) < sqrt(eps) )then
2657  gmask(i,j) = dd
2658  endif
2659  enddo
2660  enddo
2661 
2662  return
real(rp), public const_undef
Definition: scale_const.F90:43
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
Here is the caller graph for this function:

◆ replace_misval_const()

subroutine mod_realinput::replace_misval_const ( real(rp), dimension(:,:), intent(inout)  data,
real(rp), intent(in)  maskval,
real(rp), dimension(:,:), intent(in)  frac_land 
)

Definition at line 2774 of file mod_realinput.f90.

References scale_const::const_eps, scale_grid_index::ia, and scale_grid_index::ja.

Referenced by land_interporation().

2774  use scale_const, only: &
2775  eps => const_eps
2776  implicit none
2777  real(RP), intent(inout) :: data(:,:)
2778  real(RP), intent(in) :: maskval
2779  real(RP), intent(in) :: frac_land(:,:)
2780  integer :: i, j
2781 
2782  do j = 1, ja
2783  do i = 1, ia
2784  if( abs(frac_land(i,j)-0.0_rp) < eps )then ! ocean grid
2785  data(i,j) = maskval
2786  endif
2787  enddo
2788  enddo
2789 
integer, public ia
of x whole cells (local, with HALO)
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

◆ replace_misval_map()

subroutine mod_realinput::replace_misval_map ( real(rp), dimension(:,:), intent(inout)  data,
real(rp), dimension(:,:), intent(in)  maskval,
integer, intent(in)  nx,
integer, intent(in)  ny,
character(*), intent(in)  elem 
)

Definition at line 2794 of file mod_realinput.f90.

References scale_const::const_eps, scale_const::const_pi, scale_const::const_undef, scale_tracer::i_nc, scale_tracer::i_ng, scale_tracer::i_ni, scale_tracer::i_nr, scale_tracer::i_ns, scale_tracer::i_qc, scale_tracer::i_qg, scale_tracer::i_qi, scale_tracer::i_qr, scale_tracer::i_qs, and scale_process::prc_mpistop().

Referenced by land_interporation(), and parentatomsetup().

2794  use scale_const, only: &
2795  eps => const_eps, &
2796  undef => const_undef
2797  implicit none
2798  real(RP), intent(inout) :: data(:,:)
2799  real(RP), intent(in) :: maskval(:,:)
2800  integer, intent(in) :: nx, ny
2801  character(*), intent(in) :: elem
2802  integer :: i, j
2803 
2804  do j = 1, ny
2805  do i = 1, nx
2806  if( abs(data(i,j) - undef) < sqrt(eps) )then
2807  if( abs(maskval(i,j) - undef) < sqrt(eps) )then
2808  write(*,*) "Data for mask has missing value. ",trim(elem),i,j
2809  call prc_mpistop
2810  else
2811  data(i,j) = maskval(i,j)
2812  endif
2813  endif
2814  enddo
2815  enddo
2816 
subroutine, public prc_mpistop
Abort MPI.
real(rp), public const_undef
Definition: scale_const.F90:43
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
Here is the call graph for this function:
Here is the caller graph for this function: