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]
    INIT_OCEAN_ALB_LW real(RP)
    INIT_OCEAN_ALB_SW real(RP)
    INIT_OCEAN_Z0W real(RP)
    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 162 of file mod_realinput.f90.

References mod_atmos_admin::atmos_phy_mp_type, mod_atmos_vars::dens, scale_grid_index::ia, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, 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, and mod_atmos_vars::rhot.

Referenced by mod_mkinit::interporation_fact().

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
module ATMOS admin
real(rp), dimension(:,:,:), allocatable, target, public momz
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:,:), allocatable, target, public rhot
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:), allocatable, target, public momx
real(rp), dimension(:,:,:), allocatable, target, public dens
character(len=h_short), public atmos_phy_mp_type
real(rp), dimension(:,:,:), allocatable, target, public momy
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
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_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, 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  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
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:95
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
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
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]
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:96
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
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
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
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 772 of file mod_realinput.f90.

References scale_atmos_hydrometeor::atmos_hydrometeor_diagnose_number_concentration(), mod_atmos_admin::atmos_phy_mp_type, scale_const::const_eps, scale_const::const_i_lw, scale_const::const_i_sw, scale_const::const_undef, scale_fileio::fileio_create(), scale_fileio::fileio_def_var(), scale_fileio::fileio_enddef(), scale_land_grid::grid_lcz, scale_gridtrans::gtrans_rotc, scale_atmos_hydrometeor::i_qc, scale_atmos_hydrometeor::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, 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_atmos_phy_mp::qa_mp, scale_atmos_phy_mp::qe_mp, scale_atmos_hydrometeor::qle, scale_atmos_hydrometeor::qls, scale_atmos_phy_mp::qs_mp, replace_misval_map(), scale_time::time_nowdate, scale_tracer::tracer_cv, scale_tracer::tracer_mass, scale_tracer::tracer_name, and scale_tracer::tracer_r.

Referenced by realinput_atmos().

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
subroutine, public parentatomsetupgrads(dims, basename)
Atmos Setup.
subroutine, public prc_mpistop
Abort MPI.
module REAL input WRF-ARW
integer, parameter, public inicam
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
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
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 2353 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_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(), make_mask(), scale_process::prc_mpistop(), replace_misval_const(), replace_misval_map(), and scale_topography::topo_zsfc.

Referenced by parentatomsetup().

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
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:95
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
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:96
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
module TOPOGRAPHY
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 2684 of file mod_realinput.f90.

References scale_const::const_eps, scale_const::const_undef, scale_stdio::io_fid_log, and scale_stdio::io_l.

Referenced by land_interporation(), and parentatomsetup().

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
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 2829 of file mod_realinput.f90.

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

Referenced by land_interporation().

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 
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_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(len=*), intent(in)  elem 
)

Definition at line 2849 of file mod_realinput.f90.

References scale_const::const_eps, scale_const::const_undef, and scale_process::prc_mpistop().

Referenced by land_interporation(), and parentatomsetup().

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 
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: