SCALE-RM
Functions/Subroutines | Variables
mod_realinput Module Reference

module REAL input More...

Functions/Subroutines

subroutine, public realinput_atmos
 
subroutine, public realinput_surface
 
subroutine land_interporation (kmax, imax, jmax, oimax, ojmax, tg, strg, lst, albg, tg_org, strg_org, smds_org, lst_org, albg_org, sst_org, lmask_org, lsmask_nest, topo_org, lz_org, llon_org, llat_org, LCZ, CX, CY, LON, LAT, maskval_tg, maskval_strg, init_landwater_ratio, use_file_landwater, use_waterratio, soilwater_ds2vc_flag, elevation_correction, intrp_iter_max, ol_interp)
 
subroutine ocean_interporation (imax, jmax, sst_org, tw_org, albw_org, z0w_org, CX, CY, elevation_correction_ocean, init_ocean_alb_lw, init_ocean_alb_sw, init_ocean_z0w, first_surface, update_coord, sst, tw, albw, z0w)
 
subroutine urban_input (lst, albg, tc_urb, qc_urb, uc_urb, ust, albu)
 
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)
 

Variables

integer, parameter, public iscale = 1
 
integer, parameter, public iwrfarw = 2
 
integer, parameter, public igrads = 4
 

Detailed Description

module REAL input

Description
read data from file for real atmospheric simulations
Author
Team SCALE
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
    SERIAL_PROC_READ logical .true. read by one MPI process and broadcast
    FILETYPE_ORG character(len=H_LONG) ''
    BASENAME_ORG character(len=H_LONG) ''
    BASENAME_ADD_NUM logical .false.
    BASENAME_BOUNDARY character(len=H_LONG) ''
    BOUNDARY_POSTFIX_TIMELABEL logical .false.
    BOUNDARY_TITLE character(len=H_LONG) 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
    BOUNDARY_DTYPE character(len=H_SHORT) 'DEFAULT'
    BOUNDARY_UPDATE_DT real(DP) 0.0_DP inteval time of boudary data update [s]
    FILTER_ORDER integer 8 order of the hyper-diffusion (must be even)
    FILTER_NITER integer 0 times for hyper-diffusion iteration
    USE_FILE_DENSITY logical .false. use density data from files
    USE_NONHYDRO_DENS_BOUNDARY logical .false. > use non-hydrostatic density for boundary data
    USE_SFC_DIAGNOSES logical .false. > use surface diagnoses
    USE_DATA_UNDER_SFC logical .true. > use data under the surface
    SAME_MP_TYPE logical .false. microphysics type of the parent model is same as it in this model
    INTRP_TYPE character(len=H_SHORT) "LINEAR" ! "LINEAR" or "DIST-WEIGHT"
    SKIP_VERTICAL_RANGE_CHECK logical .false. > skip chkecking if the domain top does not exceed that of the parent in the vertical direction

  • 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=H_LONG) ''
    BASENAME_ADD_NUM logical .false.
    BASENAME_BOUNDARY character(len=H_LONG) ''
    BOUNDARY_POSTFIX_TIMELABEL logical .false.
    BOUNDARY_TITLE character(len=H_LONG) 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
    BOUNDARY_UPDATE_DT real(DP) 0.0_DP inteval time of boudary data update [s]
    USE_FILE_LANDWATER logical
    INIT_LANDWATER_RATIO real(RP)

  • 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=H_LONG) ''
    BASENAME_ADD_NUM logical .false.
    BASENAME_BOUNDARY character(len=H_LONG) ''
    BOUNDARY_POSTFIX_TIMELABEL logical .false.
    BOUNDARY_TITLE character(len=H_LONG) 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
    BOUNDARY_UPDATE_DT real(DP) 0.0_DP 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_TYPE character(len=H_SHORT) "LINEAR" ! "LINEAR" or "DIST-WEIGHT"
    INTRP_OCEAN_TEMP character(len=*)
    INTRP_OCEAN_SFC_TEMP character(len=*)
    INTRP_ITER_MAX integer
    FILTER_ORDER integer 8 order of the hyper-diffusion (must be even)
    FILTER_NITER integer 0 times for hyper-diffusion iteration
    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

Definition at line 193 of file mod_realinput.F90.

193  use scale_const, only: &
194  p00 => const_pre00
195  use scale_time, only: &
197  use mod_atmos_vars, only: &
198  dens, &
199  momz, &
200  momx, &
201  momy, &
202  rhot, &
203  qtrc
204  use mod_atmos_admin, only: &
206  use scale_atmos_thermodyn, only: &
207  atmos_thermodyn_specific_heat
208  implicit none
209 
210  logical :: USE_SFC_DIAGNOSES = .false.
211  logical :: USE_DATA_UNDER_SFC = .true.
212  logical :: USE_NONHYDRO_DENS_BOUNDARY = .false.
213  logical :: SKIP_VERTICAL_RANGE_CHECK = .false.
214 
215 
216  namelist / param_mkinit_real_atmos / &
217  number_of_files, &
218  number_of_tsteps, &
219  number_of_skip_tsteps, &
220  serial_proc_read, &
221  filetype_org, &
222  basename_org, &
223  basename_add_num, &
224  basename_boundary, &
225  boundary_postfix_timelabel, &
226  boundary_title, &
227  boundary_dtype, &
228  boundary_update_dt, &
229  filter_order, &
230  filter_niter, &
231  use_file_density, &
232  use_nonhydro_dens_boundary, &
233  use_sfc_diagnoses, &
234  use_data_under_sfc, &
235  same_mp_type, &
236  intrp_type, &
237  skip_vertical_range_check
238 
239  character(len=H_LONG) :: basename_mod
240  character(len=H_LONG) :: basename_out_mod
241  character(len=19) :: timelabel
242 
243  integer :: dims(6) ! dims 1-3: normal, 4-6: staggerd
244  integer :: timelen
245 
246  integer :: fid_atmos
247  integer :: vid_atmos(5+QA)
248 
249  real(RP) :: DENS_in(KA,IA,JA)
250  real(RP) :: MOMZ_in(KA,IA,JA) ! staggered point
251  real(RP) :: MOMX_in(KA,IA,JA) ! staggered point
252  real(RP) :: MOMY_in(KA,IA,JA) ! staggered point
253  real(RP) :: RHOT_in(KA,IA,JA)
254  real(RP) :: QTRC_in(KA,IA,JA,QA)
255 
256  real(RP) :: VELZ_in(KA,IA,JA) ! staggered point
257  real(RP) :: VELX_in(KA,IA,JA) ! staggered point
258  real(RP) :: VELY_in(KA,IA,JA) ! staggered point
259  real(RP) :: POTT_in(KA,IA,JA)
260  real(RP) :: PRES_in(KA,IA,JA)
261 
262  real(RP) :: Qdry (KA,IA,JA)
263  real(RP) :: Rtot (KA,IA,JA)
264  real(RP) :: CPtot(KA,IA,JA)
265  real(RP) :: CVtot(KA,IA,JA)
266 
267  integer :: ifile, istep, t, tall
268  integer :: k, i, j, iq
269  integer :: ierr
270  !---------------------------------------------------------------------------
271 
272  log_newline
273  log_info('REALINPUT_atmos',*) 'Setup'
274 
275  !--- read namelist
276  rewind(io_fid_conf)
277  read(io_fid_conf,nml=param_mkinit_real_atmos,iostat=ierr)
278  if( ierr < 0 ) then !--- missing
279  log_info("REALINPUT_atmos",*) 'Not found namelist. Default used.'
280  elseif( ierr > 0 ) then !--- fatal error
281  log_error("REALINPUT_atmos",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_ATMOS. Check!'
282  call prc_abort
283  endif
284  log_nml(param_mkinit_real_atmos)
285 
286  if ( boundary_update_dt <= 0.0_dp ) then
287  log_error("REALINPUT_atmos",*) 'BOUNDARY_UPDATE_DT is necessary in real case preprocess'
288  call prc_abort
289  endif
290 
291  if ( filetype_org == 'GrADS' ) then
292  basename_mod = trim(basename_org) ! namelist file name
293  else
294  if ( number_of_files > 1 .OR. basename_add_num ) then
295  basename_mod = trim(basename_org)//'_00000'
296  else
297  basename_mod = trim(basename_org)
298  endif
299  endif
300 
301  select case( intrp_type )
302  case ( "LINEAR" )
303  itp_nh_a = 4
304  itp_type_a = i_intrp_linear
305  case ( "DIST-WEIGHT" )
306  itp_nh_a = comm_cartesc_nest_interp_level
307  itp_type_a = i_intrp_dstwgt
308  case default
309  log_error("REALINPUT_atmos",*) 'Unsupported type of INTRP_TYPE : ', trim(intrp_type)
310  log_error_cont(*) ' It must be "LINEAR" or "DIST-WEIGHT"'
311  call prc_abort
312  end select
313 
314  call parentatmossetup( filetype_org, & ![IN]
315  basename_mod, & ![IN]
316  serial_proc_read, & ![IN]
317  use_file_density, & ![IN]
318  dims(:), & ![OUT]
319  timelen ) ![OUT]
320 
321  if ( timelen > 0 ) then
322  number_of_tsteps = timelen ! read from file
323  endif
324 
325  log_newline
326  log_info("REALINPUT_atmos",*) 'Number of temporal data in each file : ', number_of_tsteps
327 
328  do ifile = 1, number_of_files
329 
330  if ( filetype_org == 'GrADS' ) then
331  if ( number_of_files > 1 .OR. basename_add_num ) then
332  write(basename_mod,'(A,I5.5)') '_', ifile-1 ! only the number postfix
333  else
334  basename_mod = ''
335  endif
336  else
337  if ( number_of_files > 1 .OR. basename_add_num ) then
338  write(basename_mod,'(A,A,I5.5)') trim(basename_org), '_', ifile-1
339  else
340  basename_mod = trim(basename_org)
341  endif
342  endif
343 
344  log_newline
345  log_info("REALINPUT_atmos",*) 'read external data from : ', trim(basename_mod)
346 
347  call parentatmosopen( filetype_org, & ![IN]
348  basename_mod, & ![IN]
349  dims(:) ) ![IN]
350 
351  do istep = 1, number_of_tsteps
352 
353  tall = number_of_tsteps * (ifile-1) + istep ! consecutive time step (input)
354  t = tall - number_of_skip_tsteps ! time step (output)
355 
356  if ( t <= 0 ) then
357  log_progress('(1x,A,I4,A,I5,A,I6,A)') &
358  '[file,step,cons.] = [', ifile, ',', istep, ',', tall, '] ...skip.'
359  cycle
360  endif
361 
362  if ( t == 1 .OR. basename_boundary /= '' ) then
363 
364  log_progress('(1x,A,I4,A,I5,A,I6,A)') &
365  '[file,step,cons.] = [', ifile, ',', istep, ',', tall, ']'
366 
367  ! read prepared data
368  call parentatmosinput( filetype_org, & ! [IN]
369  basename_mod, & ! [IN]
370  dims(:), & ! [IN]
371  istep, & ! [IN]
372  use_sfc_diagnoses, & ! [IN]
373  use_data_under_sfc, & ! [IN]
374  same_mp_type, & ! [IN]
375  skip_vertical_range_check, & ! [IN]
376  dens_in(:,:,:), & ! [OUT]
377  momz_in(:,:,:), & ! [OUT]
378  momx_in(:,:,:), & ! [OUT]
379  momy_in(:,:,:), & ! [OUT]
380  rhot_in(:,:,:), & ! [OUT]
381  qtrc_in(:,:,:,:), & ! [OUT]
382  velz_in(:,:,:), & ! [OUT]
383  velx_in(:,:,:), & ! [OUT]
384  vely_in(:,:,:), & ! [OUT]
385  pott_in(:,:,:), & ! [OUT]
386  pres_in(:,:,:) ) ! [OUT]
387  else
388  log_progress('(1x,A,I4,A,I5,A,I6,A)') &
389  '[file,step,cons.] = [', ifile, ',', istep, ',', tall, '] ...skip.'
390  endif
391 
392  !--- store prognostic variables as initial
393  if ( t == 1 ) then
394  log_newline
395  log_info("REALINPUT_atmos",*) 'store initial state.'
396 
397  !$omp parallel do collapse(2)
398  do j = 1, ja
399  do i = 1, ia
400  do k = 1, ka
401  dens(k,i,j) = dens_in(k,i,j)
402  momz(k,i,j) = momz_in(k,i,j)
403  momx(k,i,j) = momx_in(k,i,j)
404  momy(k,i,j) = momy_in(k,i,j)
405  rhot(k,i,j) = rhot_in(k,i,j)
406  enddo
407  enddo
408  enddo
409 
410  !$omp parallel do collapse(3)
411  do iq = 1, qa
412  do j = 1, ja
413  do i = 1, ia
414  do k = 1, ka
415  qtrc(k,i,j,iq) = qtrc_in(k,i,j,iq)
416  enddo
417  enddo
418  enddo
419  enddo
420 
421  endif
422 
423  !--- output boundary data
424  if ( basename_boundary /= '' ) then
425 
426  if ( t == 1 ) then
427  if ( boundary_postfix_timelabel ) then
428  call time_gettimelabel( timelabel )
429  basename_out_mod = trim(basename_boundary)//'_'//trim(timelabel)
430  else
431  basename_out_mod = trim(basename_boundary)
432  endif
433 
434  call boundaryatmossetup( basename_out_mod, & ! [IN]
435  boundary_title, & ! [IN]
436  boundary_dtype, & ! [IN]
437  boundary_update_dt, & ! [IN]
438  fid_atmos, & ! [OUT]
439  vid_atmos(:) ) ! [OUT]
440  endif
441 
442  if ( use_nonhydro_dens_boundary ) then
443  call atmos_thermodyn_specific_heat( ka, ks, ke, ia, 1, ia, ja, 1, ja, qa, &
444  qtrc_in(:,:,:,:), & ! [IN]
445  tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), & ! [IN]
446  qdry(:,:,:), rtot(:,:,:), cvtot(:,:,:), cptot(:,:,:) ) ! [OUT]
447  !$omp parallel do collapse(2)
448  do j = 1, ja
449  do i = 1, ia
450  do k = ks, ke
451  dens_in(k,i,j) = ( pres_in(k,i,j) / p00 )**( cvtot(k,i,j) / cptot(k,i,j) ) * p00 / ( rtot(k,i,j) * pott_in(k,i,j) )
452  end do
453  end do
454  end do
455  end if
456 
457  call boundaryatmosoutput( dens_in(:,:,:), & ! [IN]
458  velz_in(:,:,:), & ! [IN]
459  velx_in(:,:,:), & ! [IN]
460  vely_in(:,:,:), & ! [IN]
461  pott_in(:,:,:), & ! [IN]
462  qtrc_in(:,:,:,:), & ! [IN]
463  fid_atmos, & ! [IN]
464  vid_atmos(:), & ! [IN]
465  boundary_update_dt, & ! [IN]
466  t ) ! [IN]
467  endif
468 
469  enddo ! istep loop
470  enddo ! ifile loop
471 
472  return

References mod_atmos_admin::atmos_phy_mp_type, scale_comm_cartesc_nest::comm_cartesc_nest_interp_level, scale_const::const_pre00, mod_atmos_vars::dens, scale_atmos_grid_cartesc_index::ia, scale_io::io_fid_conf, scale_atmos_grid_cartesc_index::ja, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, mod_atmos_vars::momx, mod_atmos_vars::momy, mod_atmos_vars::momz, scale_prc::prc_abort(), scale_tracer::qa, mod_atmos_vars::qtrc, mod_atmos_vars::rhot, scale_time::time_gettimelabel(), scale_tracer::tracer_cp, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by mod_mkinit::read_sounding().

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 477 of file mod_realinput.F90.

477  use scale_const, only: &
478  tem00 => const_tem00
479  use scale_time, only: &
481  use scale_landuse, only: &
482  fact_ocean => landuse_fact_ocean, &
483  fact_land => landuse_fact_land, &
484  fact_urban => landuse_fact_urban, &
486  use mod_atmos_phy_sf_vars, only: &
492  use mod_ocean_admin, only: &
493  ocean_do
494  use scale_ocean_phy_ice_simple, only: &
496  use mod_ocean_vars, only: &
497  ice_flag, &
498  ocean_temp, &
499  ocean_salt, &
500  ocean_uvel, &
501  ocean_vvel, &
502  ocean_ocn_z0m, &
503  ocean_ice_temp, &
504  ocean_ice_mass, &
505  ocean_sfc_temp, &
507  ocean_sfc_z0m, &
508  ocean_sfc_z0h, &
510  use mod_land_admin, only: &
511  land_do
512  use mod_land_vars, only: &
513  land_temp, &
514  land_water, &
515  land_ice, &
516  land_sfc_temp, &
518  use mod_urban_admin, only: &
519  urban_do
520  use mod_urban_vars, only: &
521  urban_tc, &
522  urban_qc, &
523  urban_uc, &
524  urban_tr, &
525  urban_tb, &
526  urban_tg, &
527  urban_trl, &
528  urban_tbl, &
529  urban_tgl, &
530  urban_rainr, &
531  urban_rainb, &
532  urban_raing, &
533  urban_roff, &
534  urban_sfc_temp, &
536  implicit none
537 
538  logical :: USE_FILE_LANDWATER = .true. ! use land water data from files
539  real(RP) :: INIT_LANDWATER_RATIO = 0.5_rp ! Ratio of land water to storage is constant, if USE_FILE_LANDWATER is ".false." (all PFT)
540 ! real(RP) :: INIT_LANDWATER_RATIO_EACH(LANDUSE_PFT_nmax) ! Ratio of land water to storage is constant, if USE_FILE_LANDWATER is ".false." (each PFT)
541  real(RP) :: INIT_OCEAN_ALB_LW = 0.04_rp ! initial LW albedo on the ocean
542  real(RP) :: INIT_OCEAN_ALB_SW = 0.10_rp ! initial SW albedo on the ocean
543  real(RP) :: INIT_OCEAN_Z0W = 1.0e-3_rp ! initial surface roughness on the ocean
544  character(len=H_SHORT) :: INTRP_LAND_TEMP = 'off'
545  character(len=H_SHORT) :: INTRP_LAND_WATER = 'off'
546  character(len=H_SHORT) :: INTRP_LAND_SFC_TEMP = 'off'
547  character(len=H_SHORT) :: INTRP_OCEAN_TEMP = 'off'
548  character(len=H_SHORT) :: INTRP_OCEAN_SFC_TEMP = 'off'
549  integer :: INTRP_ITER_MAX = 100
550  character(len=H_SHORT) :: SOILWATER_DS2VC = 'limit'
551  logical :: soilwater_DS2VC_flag ! true: 'critical', false: 'limit'
552  logical :: elevation_correction = .true.
553  logical :: elevation_correction_land
554  logical :: elevation_correction_ocean
555 
556  namelist / param_mkinit_real_land / &
557  number_of_files, &
558  number_of_tsteps, &
559  number_of_skip_tsteps, &
560  filetype_org, &
561  basename_org, &
562  basename_add_num, &
563  basename_boundary, &
564  boundary_postfix_timelabel, &
565  boundary_title, &
566  boundary_update_dt, &
567  use_file_landwater, &
568  init_landwater_ratio, &
569 ! INIT_LANDWATER_RATIO_EACH, &
570  intrp_type, &
571  intrp_land_temp, &
572  intrp_land_water, &
573  intrp_land_sfc_temp, &
574  intrp_iter_max, &
575  filter_order, &
576  filter_niter, &
577  soilwater_ds2vc, &
578  elevation_correction, &
579  serial_proc_read
580 
581  namelist / param_mkinit_real_ocean / &
582  number_of_files, &
583  number_of_tsteps, &
584  number_of_skip_tsteps, &
585  filetype_org, &
586  basename_org, &
587  basename_add_num, &
588  basename_boundary, &
589  boundary_postfix_timelabel, &
590  boundary_title, &
591  boundary_update_dt, &
592  init_ocean_alb_lw, &
593  init_ocean_alb_sw, &
594  init_ocean_z0w, &
595  intrp_type, &
596  intrp_ocean_temp, &
597  intrp_ocean_sfc_temp, &
598  intrp_iter_max, &
599  filter_order, &
600  filter_niter, &
601  serial_proc_read
602 
603  character(len=H_LONG) :: FILETYPE_LAND
604  character(len=H_LONG) :: FILETYPE_OCEAN
605  character(len=H_LONG) :: BASENAME_LAND
606  character(len=H_LONG) :: BASENAME_OCEAN
607  character(len=5) :: NUM = ''
608 
609  ! land
610  real(RP), allocatable :: LAND_TEMP_org (:,:,:,:)
611  real(RP), allocatable :: LAND_WATER_org (:,:,:,:)
612  real(RP), allocatable :: LAND_SFC_TEMP_org (:,:,:)
613  real(RP), allocatable :: LAND_SFC_albedo_org(:,:,:,:,:)
614 
615  ! urban
616  real(RP) :: URBAN_TC_ORG(IA,JA)
617  real(RP) :: URBAN_QC_ORG(IA,JA)
618  real(RP) :: URBAN_UC_ORG(IA,JA)
619  real(RP) :: URBAN_SFC_TEMP_ORG(IA,JA)
620  real(RP) :: URBAN_SFC_albedo_ORG(IA,JA,N_RAD_DIR,N_RAD_RGN)
621 
622  ! ocean
623  real(RP), allocatable :: OCEAN_TEMP_org (:,:,:,:)
624  real(RP), allocatable :: OCEAN_SFC_TEMP_org (:,:,:)
625  real(RP), allocatable :: OCEAN_SFC_albedo_org(:,:,:,:,:)
626  real(RP), allocatable :: OCEAN_SFC_Z0_org (:,:,:)
627 
628  integer :: NUMBER_OF_FILES_LAND = 1
629  integer :: NUMBER_OF_FILES_OCEAN = 1
630  integer :: NUMBER_OF_TSTEPS_LAND = 1 ! num of time steps in one file
631  integer :: NUMBER_OF_TSTEPS_OCEAN = 1 ! num of time steps in one file
632  integer :: NUMBER_OF_SKIP_TSTEPS_LAND = 0 ! num of skipped first several data
633  integer :: NUMBER_OF_SKIP_TSTEPS_OCEAN = 0 ! num of skipped first several data
634 
635  character(len=H_LONG) :: BASENAME_BOUNDARY_LAND = ''
636  character(len=H_LONG) :: BASENAME_BOUNDARY_OCEAN = ''
637  logical :: BOUNDARY_POSTFIX_TIMELABEL_LAND = .false.
638  logical :: BOUNDARY_POSTFIX_TIMELABEL_OCEAN = .false.
639  character(len=H_LONG) :: BOUNDARY_TITLE_LAND = 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
640  character(len=H_LONG) :: BOUNDARY_TITLE_OCEAN = 'SCALE-RM BOUNDARY CONDITION for REAL CASE'
641  real(DP) :: BOUNDARY_UPDATE_DT_LAND = 0.0_dp ! inteval time of boudary data update [s]
642  real(DP) :: BOUNDARY_UPDATE_DT_OCEAN = 0.0_dp ! inteval time of boudary data update [s]
643  logical :: BASENAME_ADD_NUM_LAND
644  logical :: BASENAME_ADD_NUM_OCEAN
645 
646  integer :: mdlid_land, mdlid_ocean
647  integer :: ldims(3), odims(2)
648 
649  integer :: totaltimesteps = 1
650  integer :: timelen
651  integer :: skip_steps, skip_steps_land
652  integer :: ierr
653 
654  character(len=H_LONG) :: basename_out_mod
655  character(len=19) :: timelabel
656 
657  logical :: land_flag
658  logical :: multi_land
659  logical :: multi_ocean
660 
661  integer :: ns, ne, nsl, nel
662  integer :: idir, irgn
663 
664  integer :: k, i, j, n
665  !---------------------------------------------------------------------------
666 
667  if ( land_do .or. urban_do ) then
668  land_flag = .true.
669  else
670  land_flag = .false.
671  end if
672 
673  if ( .not. land_flag .or. .not. ocean_do ) then
674  log_error("REALINPUT_surface",*) 'OCEAN_ and LAND_DYN_TYPE must be set'
675  end if
676 
677 
678  log_newline
679  log_info('REALINPUT_surface',*) 'Setup LAND'
680 
681  ! LAND/URBAN
682 ! INIT_LANDWATER_RATIO_EACH(:) = -1.0_RP
683 
684  !--- read namelist
685  rewind(io_fid_conf)
686  read(io_fid_conf,nml=param_mkinit_real_land,iostat=ierr)
687  if( ierr < 0 ) then !--- missing
688  log_info("REALINPUT_surface",*) 'Not found namelist. Default used.'
689  elseif( ierr > 0 ) then !--- fatal error
690  log_error("REALINPUT_surface",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_LAND. Check!'
691  call prc_abort
692  endif
693  log_nml(param_mkinit_real_land)
694 
695  number_of_files_land = number_of_files
696  number_of_tsteps_land = number_of_tsteps
697  number_of_skip_tsteps_land = number_of_skip_tsteps
698  filetype_land = filetype_org
699  basename_add_num_land = basename_add_num
700  basename_boundary_land = basename_boundary
701  boundary_postfix_timelabel_land = boundary_postfix_timelabel
702  boundary_title_land = boundary_title
703  boundary_update_dt_land = boundary_update_dt
704  elevation_correction_land = elevation_correction
705 
706  if ( filetype_land .ne. "GrADS" .and. ( number_of_files > 1 .OR. basename_add_num_land ) ) then
707  basename_land = trim(basename_org)//"_00000"
708  else
709  basename_land = trim(basename_org)
710  endif
711 
712 !!$ if( .NOT. USE_FILE_LANDWATER ) then
713 !!$ if( all( INIT_LANDWATER_RATIO_EACH(:) < 0.0_RP ) ) then
714 !!$ LOG_INFO("REALINPUT_surface",*) 'Applied INIT_LANDWATER_RATIO, instead of INIT_LANDWATER_RATIO_EACH.'
715 !!$ INIT_LANDWATER_RATIO_EACH(:) = INIT_LANDWATER_RATIO
716 !!$ else
717 !!$ if( any( INIT_LANDWATER_RATIO_EACH(:) < 0.0_RP ) ) then
718 !!$ LOG_ERROR("REALINPUT_surface",*) 'Insufficient elemtents of array (INIT_LANDWATER_RATIO_EACH):', INIT_LANDWATER_RATIO_EACH(:)
719 !!$ call PRC_abort
720 !!$ endif
721 !!$ endif
722 !!$ endif
723 
724  select case( soilwater_ds2vc )
725  case( 'critical' )
726  soilwater_ds2vc_flag = .true.
727  case('limit' )
728  soilwater_ds2vc_flag = .false.
729  case default
730  log_error("REALINPUT_surface",*) 'Unsupported SOILWATER_DS2CV TYPE:', trim(soilwater_ds2vc)
731  call prc_abort
732  end select
733 
734  serial_land = serial_proc_read
735 
736  select case( intrp_type )
737  case ( "LINEAR" )
738  itp_nh_l = 4
739  itp_type_l = i_intrp_linear
740  case ( "DIST-WEIGHT" )
741  itp_nh_l = comm_cartesc_nest_interp_level
742  itp_type_l = i_intrp_dstwgt
743  case default
744  log_error("REALINPUT_surface",*) 'Unsupported type of INTRP_TYPE : ', trim(intrp_type)
745  log_error_cont(*) ' It must be "LINEAR" or "DIST-WEIGHT"'
746  call prc_abort
747  end select
748 
749 
750 
751  log_newline
752  log_info('REALINPUT_surface',*) 'Setup OCEAN'
753 
754  !--- read namelist
755  rewind(io_fid_conf)
756  read(io_fid_conf,nml=param_mkinit_real_ocean,iostat=ierr)
757  if( ierr < 0 ) then !--- missing
758  log_info("REALINPUT_surface",*) 'Not found namelist. Default used.'
759  elseif( ierr > 0 ) then !--- fatal error
760  log_error("REALINPUT_surface",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_OCEAN. Check!'
761  call prc_abort
762  endif
763  log_nml(param_mkinit_real_ocean)
764 
765  number_of_files_ocean = number_of_files
766  number_of_tsteps_ocean = number_of_tsteps
767  number_of_skip_tsteps_ocean = number_of_skip_tsteps
768  filetype_ocean = filetype_org
769  basename_add_num_ocean = basename_add_num
770  basename_boundary_ocean = basename_boundary
771  boundary_postfix_timelabel_ocean = boundary_postfix_timelabel
772  boundary_title_ocean = boundary_title
773  boundary_update_dt_ocean = boundary_update_dt
774  elevation_correction_ocean = elevation_correction
775 
776  if ( filetype_ocean .ne. "GrADS" .and. ( number_of_files > 1 .OR. basename_add_num_ocean ) ) then
777  basename_ocean = trim(basename_org)//"_00000"
778  else
779  basename_ocean = trim(basename_org)
780  endif
781 
782  serial_ocean = serial_proc_read
783 
784  select case( intrp_type )
785  case ( "LINEAR" )
786  itp_nh_o = 4
787  itp_type_o = i_intrp_linear
788  case ( "DIST-WEIGHT" )
789  itp_nh_o = comm_cartesc_nest_interp_level
790  itp_type_o = i_intrp_dstwgt
791  case default
792  log_error("REALINPUT_surface",*) 'Unsupported type of INTRP_TYPE : ', trim(intrp_type)
793  log_error_cont(*) ' It must be "LINEAR" or "DIST-WEIGHT"'
794  call prc_abort
795  end select
796 
797  itp_nh_ol = comm_cartesc_nest_interp_level
798 
799  multi_land = ( number_of_files_land * number_of_tsteps_land - number_of_skip_tsteps_land ) > 1
800  multi_ocean = basename_boundary_ocean .ne. ''
801 
802  if ( ( multi_land .and. multi_ocean ) .AND. &
803  ( ( number_of_files_land .NE. number_of_files_ocean ) .OR. &
804  ( number_of_tsteps_land .NE. number_of_tsteps_ocean ) .OR. &
805  ( number_of_skip_tsteps_land .NE. number_of_skip_tsteps_ocean ) .OR. &
806  ( basename_boundary_land .NE. basename_boundary_ocean ) .OR. &
807  ( boundary_postfix_timelabel_land .NEQV. boundary_postfix_timelabel_ocean ) .OR. &
808  ( boundary_title_land .NE. boundary_title_ocean ) .OR. &
809  ( boundary_update_dt_land .NE. boundary_update_dt_ocean ) ) ) then
810  log_error("REALINPUT_surface",*) 'The following LAND/OCEAN parameters must be consistent due to technical problem:'
811  log_error_cont(*) ' NUMBER_OF_FILES, NUMBER_OF_TSTEPS, NUMBER_OF_SKIP_TSTEPS,'
812  log_error_cont(*) ' BASENAME_BOUNDARY, BOUNDARY_POSTFIX_TIMELABEL, BOUNDARY_TITLE, BOUNDARY_UPDATE_DT.'
813  call prc_abort
814  end if
815 
816  call parentsurfacesetup( ldims, odims, & ![OUT]
817  mdlid_land, & ![OUT]
818  mdlid_ocean, & ![OUT]
819  timelen, & ![OUT]
820  basename_land, & ![IN]
821  basename_ocean, & ![IN]
822  filetype_land, & ![IN]
823  filetype_ocean, & ![IN]
824  use_file_landwater, & ![IN]
825  intrp_land_temp, & ![IN]
826  intrp_land_water, & ![IN]
827  intrp_land_sfc_temp, & ![IN]
828  intrp_ocean_temp, & ![IN]
829  intrp_ocean_sfc_temp ) ![IN]
830 
831  if ( timelen > 0 ) then
832  number_of_tsteps_ocean = timelen ! read from file
833  endif
834 
835  totaltimesteps = number_of_files_ocean * number_of_tsteps_ocean
836 
837  if ( multi_land ) then
838  allocate( land_temp_org(lkmax,ia,ja, 1+number_of_skip_tsteps_land:totaltimesteps) )
839  allocate( land_water_org(lkmax,ia,ja, 1+number_of_skip_tsteps_land:totaltimesteps) )
840  allocate( land_sfc_temp_org( ia,ja, 1+number_of_skip_tsteps_land:totaltimesteps) )
841  allocate( land_sfc_albedo_org( ia,ja,n_rad_dir,n_rad_rgn,1+number_of_skip_tsteps_land:totaltimesteps) )
842  else
843  allocate( land_temp_org(lkmax,ia,ja, 1) )
844  allocate( land_water_org(lkmax,ia,ja, 1) )
845  allocate( land_sfc_temp_org( ia,ja, 1) )
846  allocate( land_sfc_albedo_org( ia,ja,n_rad_dir,n_rad_rgn,1) )
847  end if
848 
849  allocate( ocean_temp_org(okmax,ia,ja, 1+number_of_skip_tsteps_ocean:totaltimesteps) )
850  allocate( ocean_sfc_temp_org( ia,ja, 1+number_of_skip_tsteps_ocean:totaltimesteps) )
851  allocate( ocean_sfc_albedo_org( ia,ja,n_rad_dir,n_rad_rgn,1+number_of_skip_tsteps_ocean:totaltimesteps) )
852  allocate( ocean_sfc_z0_org( ia,ja, 1+number_of_skip_tsteps_ocean:totaltimesteps) )
853 
854  if ( mdlid_ocean == igrads ) then
855  basename_org = ""
856  endif
857 
858  !--- read external file
859  do n = 1, number_of_files_ocean
860 
861  if ( number_of_files_land > 1 .OR. basename_add_num_land ) then
862  write(num,'(I5.5)') n-1
863  basename_land = trim(basename_org)//"_"//num
864  else
865  basename_land = trim(basename_org)
866  endif
867  if ( number_of_files_ocean > 1 .OR. basename_add_num_ocean ) then
868  write(num,'(I5.5)') n-1
869  basename_ocean = trim(basename_org)//"_"//num
870  else
871  basename_ocean = trim(basename_org)
872  endif
873 
874  log_newline
875  log_info("REALINPUT_surface",*) 'Target File Name (Land) : ', trim(basename_land)
876  log_info("REALINPUT_surface",*) 'Target File Name (Ocean): ', trim(basename_ocean)
877  log_info("REALINPUT_surface",*) 'Time Steps in One File : ', number_of_tsteps
878 
879  ns = number_of_tsteps_ocean * (n - 1) + 1
880  ne = ns + (number_of_tsteps_ocean - 1)
881 
882  if ( ne <= number_of_skip_tsteps_ocean ) then
883  log_info("REALINPUT_surface",*) ' SKIP'
884  cycle
885  endif
886 
887  skip_steps = max(number_of_skip_tsteps_ocean - ns + 1, 0)
888  ns = max(ns, number_of_skip_tsteps_ocean+1)
889 
890  skip_steps_land = max(number_of_skip_tsteps_land - ns + 1, 0)
891 
892  if ( multi_land ) then
893  nsl = ns
894  nel = ne
895  else
896  nsl = 1
897  nel = 1
898  end if
899 
900  ! read all prepared data
901  call parentsurfaceinput( land_temp_org(:,:,:, nsl:nel), &
902  land_water_org(:,:,:, nsl:nel), &
903  land_sfc_temp_org(:,:, nsl:nel), &
904  land_sfc_albedo_org(:,:,:,:,nsl:nel), &
905  urban_tc_org(:,:), &
906  urban_qc_org(:,:), &
907  urban_uc_org(:,:), &
908  urban_sfc_temp_org(:,:), &
909  urban_sfc_albedo_org(:,:,:,:), &
910  ocean_temp_org(oks,:,:, ns:ne), &
911  ocean_sfc_temp_org( :,:, ns:ne), &
912  ocean_sfc_albedo_org( :,:,:,:,ns:ne), &
913  ocean_sfc_z0_org( :,:, ns:ne), &
914  basename_land, basename_ocean, &
915  mdlid_land, mdlid_ocean, &
916  ldims, odims, &
917  use_file_landwater, &
918  init_landwater_ratio, &
919 ! INIT_LANDWATER_RATIO_EACH(:), &
920  init_ocean_alb_lw, init_ocean_alb_sw, &
921  init_ocean_z0w, &
922  intrp_iter_max, &
923  soilwater_ds2vc_flag, &
924  elevation_correction_land, &
925  elevation_correction_ocean, &
926  multi_land, multi_ocean, &
927  number_of_tsteps_ocean, &
928  skip_steps_land, skip_steps, &
929  urban_do )
930 
931  ! required one-step data only
932  if( .not. ( multi_land .or. multi_ocean ) ) exit
933 
934  enddo
935 
936 
937  !--- input initial data
938  ns = number_of_skip_tsteps_ocean + 1 ! skip first several data
939  if ( multi_land ) then
940  nsl = ns
941  else
942  nsl = 1
943  end if
944 
945  !$omp parallel do
946  do j = 1, ja
947  do i = 1, ia
948  ocean_sfc_temp(i,j) = ocean_sfc_temp_org(i,j,ns)
949  ocean_sfc_z0m(i,j) = ocean_sfc_z0_org(i,j,ns)
950  ocean_sfc_z0h(i,j) = ocean_sfc_z0_org(i,j,ns)
951  ocean_sfc_z0e(i,j) = ocean_sfc_z0_org(i,j,ns)
952  do irgn = i_r_ir, i_r_vis
953  do idir = i_r_direct, i_r_diffuse
954  ocean_sfc_albedo(i,j,idir,irgn) = ocean_sfc_albedo_org(i,j,idir,irgn,ns)
955  enddo
956  enddo
957  do k = 1, okmax
958  ocean_temp(k,i,j) = ocean_temp_org(oks,i,j,ns)
959  ocean_salt(k,i,j) = 0.0_rp
960  ocean_uvel(k,i,j) = 0.0_rp
961  ocean_vvel(k,i,j) = 0.0_rp
962  enddo
963  ocean_ocn_z0m(i,j) = ocean_sfc_z0_org(i,j,ns)
964  if ( ice_flag ) then
965  ocean_ice_temp(i,j) = min( ocean_sfc_temp_org(i,j,ns), ocean_phy_ice_freezetemp )
966  ocean_ice_mass(i,j) = 0.0_rp
967  end if
968 
969  land_sfc_temp(i,j) = land_sfc_temp_org(i,j, nsl)
970  do irgn = i_r_ir, i_r_vis
971  do idir = i_r_direct, i_r_diffuse
972  land_sfc_albedo(i,j,idir,irgn) = land_sfc_albedo_org(i,j,idir,irgn,nsl)
973  enddo
974  enddo
975  do k = 1, lkmax
976  land_temp(k,i,j) = land_temp_org(k,i,j,nsl)
977  if ( land_temp(k,i,j) >= tem00 ) then
978  land_water(k,i,j) = land_water_org(k,i,j,nsl)
979  land_ice(k,i,j) = 0.0_rp
980  else
981  land_water(k,i,j) = 0.0_rp
982  land_ice(k,i,j) = land_water_org(k,i,j,nsl)
983  end if
984  enddo
985 
986  if ( urban_do ) then
987  urban_sfc_temp(i,j) = urban_sfc_temp_org(i,j)
988  do irgn = i_r_ir, i_r_vis
989  do idir = i_r_direct, i_r_diffuse
990  urban_sfc_albedo(i,j,idir,irgn) = urban_sfc_albedo_org(i,j,idir,irgn)
991  enddo
992  enddo
993  do k = 1, ukmax
994  urban_trl(k,i,j) = urban_sfc_temp_org(i,j)
995  urban_tbl(k,i,j) = urban_sfc_temp_org(i,j)
996  urban_tgl(k,i,j) = urban_sfc_temp_org(i,j)
997  enddo
998  urban_tc(i,j) = urban_tc_org(i,j)
999  urban_qc(i,j) = urban_qc_org(i,j)
1000  urban_uc(i,j) = urban_uc_org(i,j)
1001  urban_tr(i,j) = urban_sfc_temp_org(i,j)
1002  urban_tb(i,j) = urban_sfc_temp_org(i,j)
1003  urban_tg(i,j) = urban_sfc_temp_org(i,j)
1004  urban_rainr(i,j) = 0.0_rp
1005  urban_rainb(i,j) = 0.0_rp
1006  urban_raing(i,j) = 0.0_rp
1007  urban_roff(i,j) = 0.0_rp
1008  end if
1009 
1013 
1014  if ( urban_do ) then
1015  atmos_phy_sf_sfc_temp(i,j) = fact_ocean(i,j) * ocean_sfc_temp(i,j) &
1016  + fact_land(i,j) * land_sfc_temp(i,j) &
1017  + fact_urban(i,j) * urban_sfc_temp(i,j)
1018  do irgn = i_r_ir, i_r_vis
1019  do idir = i_r_direct, i_r_diffuse
1020  atmos_phy_sf_sfc_albedo(i,j,idir,irgn) = fact_ocean(i,j) * ocean_sfc_albedo(i,j,idir,irgn) &
1021  + fact_land(i,j) * land_sfc_albedo(i,j,idir,irgn) &
1022  + fact_urban(i,j) * urban_sfc_albedo(i,j,idir,irgn)
1023  enddo
1024  enddo
1025  else
1026  atmos_phy_sf_sfc_temp(i,j) = fact_ocean(i,j) * ocean_sfc_temp(i,j) &
1027  + fact_land(i,j) * land_sfc_temp(i,j)
1028  do irgn = i_r_ir, i_r_vis
1029  do idir = i_r_direct, i_r_diffuse
1030  atmos_phy_sf_sfc_albedo(i,j,idir,irgn) = fact_ocean(i,j) * ocean_sfc_albedo(i,j,idir,irgn) &
1031  + fact_land(i,j) * land_sfc_albedo(i,j,idir,irgn)
1032  enddo
1033  enddo
1034  endif
1035  enddo
1036  enddo
1037 
1038 
1039  !--- output boundary data
1040  if( basename_boundary_ocean /= '' ) then
1041  totaltimesteps = totaltimesteps - number_of_skip_tsteps_ocean ! skip first several data
1042  if ( totaltimesteps > 1 ) then
1043  if ( boundary_update_dt_ocean <= 0.0_dp ) then
1044  log_error("REALINPUT_surface",*) 'BOUNDARY_UPDATE_DT is necessary in real case preprocess'
1045  call prc_abort
1046  endif
1047 
1048  if ( boundary_postfix_timelabel_ocean ) then
1049  call time_gettimelabel( timelabel )
1050  basename_out_mod = trim(basename_boundary_ocean)//'_'//trim(timelabel)
1051  else
1052  basename_out_mod = trim(basename_boundary_ocean)
1053  endif
1054 
1055  if ( multi_land ) then
1056  nsl = ns
1057  nel = ne
1058  else
1059  nsl = 1
1060  nel = 1
1061  end if
1062 
1063  call parentsurfaceboundary( land_temp_org(:,:,:,nsl:nel), &
1064  land_water_org(:,:,:,nsl:nel), &
1065  land_sfc_temp_org( :,:,nsl:nel), &
1066  ocean_temp_org(:,:,:,ns:ne), &
1067  ocean_sfc_temp_org( :,:,ns:ne), &
1068  ocean_sfc_z0_org( :,:,ns:ne), &
1069  totaltimesteps, &
1070  boundary_update_dt_ocean, &
1071  basename_out_mod, &
1072  boundary_title_ocean, &
1073  multi_land )
1074 
1075  endif
1076  endif
1077 
1078  deallocate( land_temp_org )
1079  deallocate( land_water_org )
1080  deallocate( land_sfc_temp_org )
1081  deallocate( land_sfc_albedo_org )
1082  deallocate( ocean_temp_org )
1083  deallocate( ocean_sfc_temp_org )
1084  deallocate( ocean_sfc_albedo_org )
1085  deallocate( ocean_sfc_z0_org )
1086 
1087  return

References scale_atmos_grid_cartesc::atmos_grid_cartesc_cx, scale_atmos_grid_cartesc::atmos_grid_cartesc_cy, scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon, scale_atmos_hydrometeor::atmos_hydrometeor_dry, mod_atmos_admin::atmos_phy_ch_type, mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc(), 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_comm_cartesc_nest::comm_cartesc_nest_interp_level, scale_const::const_eps, scale_const::const_i_lw, scale_const::const_i_sw, scale_const::const_laps, scale_const::const_pi, scale_const::const_tem00, scale_const::const_undef, scale_precision::dp, scale_file_cartesc::file_cartesc_create(), scale_file_cartesc::file_cartesc_def_var(), scale_file_cartesc::file_cartesc_enddef(), scale_atmos_hydrometeor::i_qv, scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_atmos_grid_cartesc_index::ia, mod_ocean_vars::ice_flag, scale_atmos_grid_cartesc_index::ieb, igrads, scale_interp::interp_domain_compatibility(), scale_interp::interp_interp2d(), scale_interp::interp_interp3d(), scale_io::io_fid_conf, scale_atmos_grid_cartesc_index::isb, iscale, iwrfarw, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::jsb, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, mod_land_admin::land_do, scale_land_grid_cartesc::land_grid_cartesc_cz, mod_land_vars::land_ice, land_interporation(), 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_landuse::landuse_frac_land, scale_landuse::landuse_pft_nmax, scale_land_grid_cartesc_index::lkmax, make_mask(), scale_atmos_hydrometeor::n_hyd, scale_cpl_sfc_index::n_rad_dir, scale_cpl_sfc_index::n_rad_rgn, mod_ocean_admin::ocean_do, mod_ocean_vars::ocean_ice_mass, mod_ocean_vars::ocean_ice_temp, ocean_interporation(), mod_ocean_vars::ocean_ocn_z0m, scale_ocean_phy_ice_simple::ocean_phy_ice_freezetemp, mod_ocean_vars::ocean_salt, 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, mod_ocean_vars::ocean_uvel, mod_ocean_vars::ocean_vvel, scale_ocean_grid_cartesc_index::okmax, scale_ocean_grid_cartesc_index::oks, mod_realinput_grads::parentatmosinputgrads(), mod_realinput_scale::parentatmosinputscale(), mod_realinput_wrfarw::parentatmosinputwrfarw(), mod_realinput_grads::parentatmosopengrads(), mod_realinput_scale::parentatmosopenscale(), mod_realinput_wrfarw::parentatmosopenwrfarw(), mod_realinput_grads::parentatmossetupgrads(), mod_realinput_scale::parentatmossetupscale(), mod_realinput_wrfarw::parentatmossetupwrfarw(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_scale::parentlandinputscale(), mod_realinput_wrfarw::parentlandinputwrfarw(), mod_realinput_grads::parentlandsetupgrads(), mod_realinput_scale::parentlandsetupscale(), mod_realinput_wrfarw::parentlandsetupwrfarw(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_scale::parentoceaninputscale(), mod_realinput_wrfarw::parentoceaninputwrfarw(), mod_realinput_grads::parentoceanopengrads(), mod_realinput_scale::parentoceanopenscale(), mod_realinput_wrfarw::parentoceanopenwrfarw(), mod_realinput_grads::parentoceansetupgrads(), mod_realinput_scale::parentoceansetupscale(), mod_realinput_wrfarw::parentoceansetupwrfarw(), scale_prc::prc_abort(), scale_prc::prc_ismaster, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, mod_atmos_phy_ch_vars::qe_ch, mod_atmos_phy_mp_vars::qe_mp, scale_atmos_hydrometeor::qle, scale_atmos_hydrometeor::qls, mod_atmos_phy_ch_vars::qs_ch, mod_atmos_phy_mp_vars::qs_mp, replace_misval_map(), scale_time::time_gettimelabel(), scale_time::time_nowdate, scale_topography::topography_zsfc, scale_tracer::tracer_cp, scale_tracer::tracer_mass, scale_tracer::tracer_name, scale_tracer::tracer_r, scale_urban_grid_cartesc_index::ukmax, mod_urban_admin::urban_do, urban_input(), 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::read_sounding().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ land_interporation()

subroutine mod_realinput::land_interporation ( integer, intent(in)  kmax,
integer, intent(in)  imax,
integer, intent(in)  jmax,
integer, intent(in)  oimax,
integer, intent(in)  ojmax,
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,n_rad_dir,n_rad_rgn), intent(out)  albg,
real(rp), dimension(kmax,imax,jmax), intent(inout)  tg_org,
real(rp), dimension(kmax,imax,jmax), intent(inout)  strg_org,
real(rp), dimension(kmax,imax,jmax), intent(inout)  smds_org,
real(rp), dimension(imax,jmax), intent(inout)  lst_org,
real(rp), dimension(imax,jmax,n_rad_dir,n_rad_rgn), intent(inout)  albg_org,
real(rp), dimension(oimax,ojmax), intent(inout)  sst_org,
real(rp), dimension(imax,jmax), intent(in)  lmask_org,
real(rp), dimension(ia,ja), intent(in)  lsmask_nest,
real(rp), dimension(imax,jmax), intent(in)  topo_org,
real(rp), dimension(kmax), intent(in)  lz_org,
real(rp), dimension(imax,jmax), intent(in)  llon_org,
real(rp), dimension(imax,jmax), intent(in)  llat_org,
real(rp), dimension(lkmax), intent(in)  LCZ,
real(rp), dimension(ia), intent(in)  CX,
real(rp), dimension(ja), intent(in)  CY,
real(rp), dimension(ia,ja), intent(in)  LON,
real(rp), dimension(ia,ja), intent(in)  LAT,
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_correction,
integer, intent(in)  intrp_iter_max,
logical, intent(in)  ol_interp 
)

Definition at line 3297 of file mod_realinput.F90.

3297  use scale_prc, only: &
3298  prc_abort
3299  use scale_const, only: &
3300  undef => const_undef, &
3301  eps => const_eps, &
3302  i_sw => const_i_sw, &
3303  i_lw => const_i_lw, &
3304  pi => const_pi, &
3305  laps => const_laps
3306  use scale_interp, only: &
3307  interp_factor2d, &
3308  interp_factor3d, &
3309  interp_interp2d, &
3311  use scale_mapprojection, only: &
3312  mapprojection_lonlat2xy
3313  use scale_comm_cartesc, only: &
3314  comm_vars8, &
3315  comm_wait
3316  use scale_filter, only: &
3317  filter_hyperdiff
3318  use scale_topography, only: &
3320  use scale_landuse, only: &
3321  landuse_pft_nmax, &
3323  use mod_land_vars, only: &
3325  implicit none
3326  integer, intent(in) :: kmax, imax, jmax
3327  integer, intent(in) :: oimax, ojmax
3328  real(RP), intent(out) :: tg(LKMAX,IA,JA)
3329  real(RP), intent(out) :: strg(LKMAX,IA,JA)
3330  real(RP), intent(out) :: lst(IA,JA)
3331  real(RP), intent(out) :: albg(IA,JA,N_RAD_DIR,N_RAD_RGN)
3332  real(RP), intent(inout) :: tg_org(kmax,imax,jmax)
3333  real(RP), intent(inout) :: strg_org(kmax,imax,jmax)
3334  real(RP), intent(inout) :: smds_org(kmax,imax,jmax)
3335  real(RP), intent(inout) :: lst_org(imax,jmax)
3336  real(RP), intent(inout) :: albg_org(imax,jmax,N_RAD_DIR,N_RAD_RGN)
3337  real(RP), intent(inout) :: sst_org(oimax,ojmax)
3338  real(RP), intent(in) :: lmask_org(imax,jmax)
3339  real(RP), intent(in) :: lsmask_nest(IA,JA)
3340  real(RP), intent(in) :: topo_org(imax,jmax)
3341  real(RP), intent(in) :: lz_org(kmax)
3342  real(RP), intent(in) :: llon_org(imax,jmax)
3343  real(RP), intent(in) :: llat_org(imax,jmax)
3344  real(RP), intent(in) :: LCZ(LKMAX)
3345  real(RP), intent(in) :: CX(IA)
3346  real(RP), intent(in) :: CY(JA)
3347  real(RP), intent(in) :: LON(IA,JA)
3348  real(RP), intent(in) :: LAT(IA,JA)
3349  real(RP), intent(in) :: maskval_tg
3350  real(RP), intent(in) :: maskval_strg
3351  real(RP), intent(in) :: init_landwater_ratio
3352 ! real(RP), intent(in) :: init_landwater_ratio_each(LANDUSE_PFT_nmax)
3353  logical, intent(in) :: use_file_landwater
3354  logical, intent(in) :: use_waterratio
3355  logical, intent(in) :: soilwater_ds2vc_flag
3356  logical, intent(in) :: elevation_correction
3357  integer, intent(in) :: intrp_iter_max
3358  logical, intent(in) :: ol_interp
3359 
3360  real(RP) :: lmask(imax,jmax)
3361  real(RP) :: smds(LKMAX,IA,JA)
3362 
3363  ! data for interporation
3364  real(RP) :: hfact_l(imax,jmax,itp_nh_ol)
3365  integer :: igrd_l (imax,jmax,itp_nh_ol)
3366  integer :: jgrd_l (imax,jmax,itp_nh_ol)
3367  real(RP) :: lX_org (imax,jmax)
3368  real(RP) :: lY_org (imax,jmax)
3369  logical :: zonal, pole
3370  integer :: igrd ( IA,JA,itp_nh_l)
3371  integer :: jgrd ( IA,JA,itp_nh_l)
3372  real(RP) :: hfact( IA,JA,itp_nh_l)
3373  integer :: kgrdl (LKMAX,2,IA,JA,itp_nh_l)
3374  real(RP) :: vfactl(LKMAX, IA,JA,itp_nh_l)
3375 
3376 
3377  real(RP) :: sst_land(imax,jmax)
3378  real(RP) :: work (imax,jmax)
3379  real(RP) :: work2(imax,jmax)
3380 
3381  real(RP) :: lz3d_org(kmax,imax,jmax)
3382  real(RP) :: lcz_3D(LKMAX,IA,JA)
3383 
3384  ! elevation correction
3385  real(RP) :: topo(IA,JA)
3386  real(RP) :: tdiff
3387 
3388  real(RP) :: one2d(IA,JA)
3389  real(RP) :: one3d(LKMAX,IA,JA)
3390 
3391  integer :: k, i, j, m, n
3392 
3393 
3394  ! Surface skin temp: interpolate over the ocean
3395  if ( i_intrp_land_sfc_temp .ne. i_intrp_off ) then
3396  select case( i_intrp_land_sfc_temp )
3397  case( i_intrp_mask )
3398  call make_mask( lmask, lst_org, imax, jmax, landdata=.true.)
3399  !$omp parallel do
3400  do j = 1, jmax
3401  do i = 1, imax
3402  if ( lmask_org(i,j) .ne. undef ) lmask(i,j) = lmask_org(i,j)
3403  end do
3404  end do
3405  case( i_intrp_fill )
3406  call make_mask( lmask, lst_org, imax, jmax, landdata=.true.)
3407  case default
3408  log_error("land_interporation",*) 'INTRP_LAND_SFC_TEMP is invalid.'
3409  call prc_abort
3410  end select
3411  call interp_oceanland_data(lst_org, lmask, imax, jmax, .true., intrp_iter_max)
3412  end if
3413 
3414  if ( ol_interp ) then
3415  ! interpolation facter between outer land grid and ocean grid
3416  call interp_factor2d( itp_nh_ol, & ! [IN]
3417  oimax, ojmax, & ! [IN]
3418  imax, jmax, & ! [IN]
3419  olon_org(:,:), & ! [IN]
3420  olat_org(:,:), & ! [IN]
3421  llon_org(:,:), & ! [IN]
3422  llat_org(:,:), & ! [IN]
3423  igrd_l(:,:,:), & ! [OUT]
3424  jgrd_l(:,:,:), & ! [OUT]
3425  hfact_l(:,:,:) ) ! [OUT]
3426 
3427  ! sst on land grid
3428  call interp_interp2d( itp_nh_ol, & ! [IN]
3429  oimax, ojmax, & ! [IN]
3430  imax, jmax, & ! [IN]
3431  igrd_l(:,:,:), & ! [IN]
3432  jgrd_l(:,:,:), & ! [IN]
3433  hfact_l(:,:,:), & ! [IN]
3434  sst_org(:,:), & ! [IN]
3435  sst_land(:,:) ) ! [OUT]
3436  else
3437  !$omp parallel do
3438  do j = 1, jmax
3439  do i = 1, imax
3440  sst_land(i,j) = sst_org(i,j)
3441  end do
3442  end do
3443  end if
3444 
3445  !$omp parallel do
3446  do j = 1, jmax
3447  do i = 1, imax
3448  if ( topo_org(i,j) > undef + eps ) then ! ignore UNDEF value
3449  sst_land(i,j) = sst_land(i,j) - topo_org(i,j) * laps
3450  end if
3451  end do
3452  end do
3453 
3454  call replace_misval_map( lst_org, sst_land, imax, jmax, "SKINT" )
3455 
3456  ! replace missing value
3457  !$omp parallel do
3458  do j = 1, jmax
3459  do i = 1, imax
3460 ! if ( skinw_org(i,j) == UNDEF ) skinw_org(i,j) = 0.0_RP
3461 ! if ( snowq_org(i,j) == UNDEF ) snowq_org(i,j) = 0.0_RP
3462 ! if ( snowt_org(i,j) == UNDEF ) snowt_org(i,j) = TEM00
3463  do m = 1, n_rad_dir
3464  if( albg_org(i,j,m,i_r_ir ) == undef ) albg_org(i,j,m,i_r_ir ) = 0.03_rp ! emissivity of general ground surface : 0.95-0.98
3465  if( albg_org(i,j,m,i_r_nir) == undef ) albg_org(i,j,m,i_r_nir) = 0.22_rp
3466  if( albg_org(i,j,m,i_r_vis) == undef ) albg_org(i,j,m,i_r_vis) = 0.22_rp
3467  end do
3468  end do
3469  end do
3470 
3471  ! Land temp: interpolate over the ocean
3472  if ( i_intrp_land_temp .ne. i_intrp_off ) then
3473  do k = 1, kmax
3474  !$omp parallel do
3475  do j = 1, jmax
3476  do i = 1, imax
3477  work(i,j) = tg_org(k,i,j)
3478  end do
3479  end do
3480  select case( i_intrp_land_temp )
3481  case( i_intrp_mask )
3482  call make_mask( lmask, work, imax, jmax, landdata=.true.)
3483  !$omp parallel do
3484  do j = 1, jmax
3485  do i = 1, imax
3486  if ( lmask_org(i,j) .ne. undef ) lmask(i,j) = lmask_org(i,j)
3487  end do
3488  end do
3489  case( i_intrp_fill )
3490  call make_mask( lmask, work, imax, jmax, landdata=.true.)
3491  end select
3492  call interp_oceanland_data( work, lmask, imax, jmax, .true., intrp_iter_max )
3493  !replace land temp using skin temp
3494  call replace_misval_map( work, lst_org, imax, jmax, "STEMP")
3495  !$omp parallel do
3496  do j = 1, jmax
3497  do i = 1, imax
3498  tg_org(k,i,j) = work(i,j)
3499  end do
3500  end do
3501  end do
3502  end if
3503 
3504 
3505  ! fill grid data
3506  !$omp parallel do collapse(2)
3507  do j = 1, jmax
3508  do i = 1, imax
3509  lz3d_org(:,i,j) = lz_org(:)
3510  end do
3511  end do
3512 
3513  !$omp parallel do collapse(2)
3514  do j = 1, ja
3515  do i = 1, ia
3516  lcz_3d(:,i,j) = lcz(:)
3517  enddo
3518  enddo
3519 
3520  select case( itp_type_l )
3521  case ( i_intrp_linear )
3522 
3523  if ( imax == 1 .or. jmax == 1 ) then
3524  log_error("land_interporation",*) 'LINER interpolation requires nx, ny > 1'
3525  log_error_cont(*) 'Use "DIST-WEIGHT" as INTRP_TYPE of PARAM_MKINIT_REAL_LAND'
3526  call prc_abort
3527  end if
3528 
3529  !$omp parallel do
3530  do j = 1, jmax
3531  do i = 1, imax
3532  work(i,j) = sign( min( abs(llat_org(i,j)), pi * 0.499999_rp ), llat_org(i,j) )
3533  end do
3534  end do
3535 
3536  call mapprojection_lonlat2xy( imax, 1, imax, jmax, 1, jmax, &
3537  llon_org(:,:), work(:,:), & ! [IN]
3538  lx_org(:,:), ly_org(:,:) ) ! [OUT]
3539 
3540  zonal = ( maxval(llon_org) - minval(llon_org) ) > 2.0_rp * pi * 0.9_rp
3541  pole = ( maxval(llat_org) > pi * 0.5_rp * 0.9_rp ) .or. ( minval(llat_org) < - pi * 0.5_rp * 0.9_rp )
3542  call interp_factor3d( kmax, 1, kmax, & ! [IN]
3543  imax, jmax, & ! [IN]
3544  lkmax, lks, lke, & ! [IN]
3545  ia, ja, & ! [IN]
3546  lx_org(:,:), & ! [IN]
3547  ly_org(:,:), & ! [IN]
3548  lz3d_org(:,:,:), & ! [IN]
3549  cx(:), cy(:), & ! [IN]
3550  lcz_3d(:,:,:), & ! [IN]
3551  igrd( :,:,:), & ! [OUT]
3552  jgrd( :,:,:), & ! [OUT]
3553  hfact( :,:,:), & ! [OUT]
3554  kgrdl(:,:,:,:,:), & ! [OUT]
3555  vfactl(:, :,:,:), & ! [OUT]
3556  flag_extrap = .true., & ! [IN, optional]
3557  zonal = zonal, & ! [IN, optional]
3558  pole = pole ) ! [IN, optional]
3559 
3560  case ( i_intrp_dstwgt )
3561 
3562  call interp_factor3d( itp_nh_l, & ! [IN]
3563  kmax, 1, kmax, & ! [IN]
3564  imax, jmax, & ! [IN]
3565  lkmax, lks, lke, & ! [IN]
3566  ia, ja, & ! [IN]
3567  llon_org(:,:), & ! [IN]
3568  llat_org(:,:), & ! [IN]
3569  lz3d_org(:,:,:), & ! [IN]
3570  lon(:,:), lat(:,:), & ! [IN]
3571  lcz_3d(:,:,:), & ! [IN]
3572  igrd( :,:,:), & ! [OUT]
3573  jgrd( :,:,:), & ! [OUT]
3574  hfact( :,:,:), & ! [OUT]
3575  kgrdl(:,:,:,:,:), & ! [OUT]
3576  vfactl(:, :,:,:), & ! [OUT]
3577  flag_extrap = .true. ) ! [IN, optional]
3578 
3579  end select
3580 
3581  call interp_interp2d( itp_nh_l, & ! [IN]
3582  imax, jmax, & ! [IN]
3583  ia, ja, & ! [IN]
3584  igrd(:,:,:), & ! [IN]
3585  jgrd(:,:,:), & ! [IN]
3586  hfact(:,:,:), & ! [IN]
3587  lst_org(:,:), & ! [IN]
3588  lst(:,:) ) ! [OUT]
3589  if ( filter_niter > 0 ) then
3590  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
3591  lst(:,:), filter_order, filter_niter )
3592  call comm_vars8( lst(:,:), 1 )
3593  call comm_wait ( lst(:,:), 1, .false. )
3594  end if
3595 
3596 
3597  if ( filter_niter > 0 ) then
3598  !$omp parallel do
3599  do j = 1, ja
3600  do i = 1, ia
3601  one2d(i,j) = 1.0_rp
3602  end do
3603  end do
3604  end if
3605 
3606  do n = 1, n_rad_rgn
3607  do m = 1, n_rad_dir
3608 
3609  call interp_interp2d( itp_nh_l, & ! [IN]
3610  imax, jmax, & ! [IN]
3611  ia, ja, & ! [IN]
3612  igrd(:,:,:), & ! [IN]
3613  jgrd(:,:,:), & ! [IN]
3614  hfact(:,:,:), & ! [IN]
3615  albg_org(:,:,m,n), & ! [IN]
3616  albg(:,:,m,n) ) ! [OUT]
3617  if ( filter_niter > 0 ) then
3618  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
3619  albg(:,:,m,n), filter_order, filter_niter, &
3620  limiter_sign = one2d(:,:) )
3621  call comm_vars8( albg(:,:,m,n), 1 )
3622  call comm_wait ( albg(:,:,m,n), 1, .false. )
3623  end if
3624  end do
3625  end do
3626 
3627  call interp_interp3d( itp_nh_l, &
3628  kmax, 1, kmax, &
3629  imax, jmax, &
3630  lkmax, lks, lke, &
3631  ia, ja, &
3632  igrd( :,:,:), & ! [IN]
3633  jgrd( :,:,:), & ! [IN]
3634  hfact( :,:,:), & ! [IN]
3635  kgrdl(:,:,:,:,:), & ! [IN]
3636  vfactl(:, :,:,:), & ! [IN]
3637  lz3d_org(:,:,:), & ! [IN]
3638  lcz_3d(:,:,:), & ! [IN]
3639  tg_org(:,:,:), & ! [IN]
3640  tg(:,:,:) ) ! [OUT]
3641 
3642  !$omp parallel do
3643  do j = 1, ja
3644  do i = 1, ia
3645  tg(lkmax,i,j) = tg(lkmax-1,i,j)
3646  enddo
3647  enddo
3648 
3649  ! replace values over the ocean
3650  do k = 1, lkmax
3651  call replace_misval_const( tg(k,:,:), maskval_tg, lsmask_nest )
3652  enddo
3653  if ( filter_niter > 0 ) then
3654  call filter_hyperdiff( lkmax, 1, lkmax, ia, isb, ieb, ja, jsb, jeb, &
3655  tg(:,:,:), filter_order, filter_niter )
3656  call comm_vars8( tg(:,:,:), 1 )
3657  call comm_wait ( tg(:,:,:), 1, .false. )
3658  end if
3659 
3660 
3661  ! elevation correction
3662  if ( elevation_correction ) then
3663  call interp_interp2d( itp_nh_l, & ! [IN]
3664  imax, jmax, & ! [IN]
3665  ia, ja, & ! [IN]
3666  igrd(:,:,:), & ! [IN]
3667  jgrd(:,:,:), & ! [IN]
3668  hfact(:,:,:), & ! [IN]
3669  topo_org(:,:), & ! [IN]
3670  topo(:,:) ) ! [OUT]
3671  if ( filter_niter > 0 ) then
3672  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
3673  topo(:,:), filter_order, filter_niter )
3674  call comm_vars8( topo(:,:), 1 )
3675  call comm_wait ( topo(:,:), 1, .false. )
3676  end if
3677 
3678  !$omp parallel do collapse(2) &
3679  !$omp private(tdiff)
3680  do j = 1, ja
3681  do i = 1, ia
3682  if ( topo(i,j) > undef + eps ) then ! ignore UNDEF value
3683  tdiff = ( topography_zsfc(i,j) - topo(i,j) ) * laps
3684  lst(i,j) = lst(i,j) - tdiff
3685  do k = 1, lkmax
3686  tg(k,i,j) = tg(k,i,j) - tdiff
3687  end do
3688  end if
3689  end do
3690  end do
3691 
3692  end if
3693 
3694 
3695 
3696  ! Land water: interpolate over the ocean
3697  if( use_file_landwater )then
3698 
3699  if ( use_waterratio ) then
3700 
3701  if ( i_intrp_land_water .ne. i_intrp_off ) then
3702  do k = 1, kmax
3703  !$omp parallel do
3704  do j = 1, jmax
3705  do i = 1, imax
3706  work(i,j) = smds_org(k,i,j)
3707  end do
3708  end do
3709  select case( i_intrp_land_water )
3710  case( i_intrp_mask )
3711  call make_mask( lmask, work, imax, jmax, landdata=.true.)
3712  !$omp parallel do
3713  do j = 1, jmax
3714  do i = 1, imax
3715  if ( lmask_org(i,j) .ne. undef ) lmask(i,j) = lmask_org(i,j)
3716  end do
3717  end do
3718  case( i_intrp_fill )
3719  call make_mask( lmask, work, imax, jmax, landdata=.true.)
3720  end select
3721  call interp_oceanland_data(work, lmask, imax, jmax, .true., intrp_iter_max)
3722  !$omp parallel do
3723  do j = 1, jmax
3724  do i = 1, imax
3725 ! work2(i,j) = init_landwater_ratio_each( LANDUSE_index_PFT(i,j,1) )
3726  work2(i,j) = init_landwater_ratio
3727  end do
3728  end do
3729  !replace missing value to init_landwater_ratio_each
3730  call replace_misval_map( work, work2, imax, jmax, "SMOISDS")
3731  !$omp parallel do
3732  do j = 1, jmax
3733  do i = 1, imax
3734  smds_org(k,i,j) = work(i,j)
3735  end do
3736  end do
3737  enddo
3738  end if
3739 
3740  call interp_interp3d( itp_nh_l, &
3741  kmax, 1, kmax, &
3742  imax, jmax, &
3743  lkmax, lks, lke, &
3744  ia, ja, &
3745  igrd( :,:,:), & ! [IN]
3746  jgrd( :,:,:), & ! [IN]
3747  hfact( :,:,:), & ! [IN]
3748  kgrdl(:,:,:,:,:), & ! [IN]
3749  vfactl(:, :,:,:), & ! [IN]
3750  lz3d_org(:,:,:), & ! [IN]
3751  lcz_3d(:,:,:), & ! [IN]
3752  smds_org(:,:,:), & ! [IN]
3753  smds(:,:,:) ) ! [OUT]
3754 
3755  do k = 1, lkmax-1
3756  strg(k,:,:) = convert_ws2vwc( smds(k,:,:), critical=soilwater_ds2vc_flag )
3757  end do
3758 
3759  else
3760 
3761  if ( i_intrp_land_water .ne. i_intrp_off ) then
3762  do k = 1, kmax
3763  !$omp parallel do
3764  do j = 1, jmax
3765  do i = 1, imax
3766  work(i,j) = strg_org(k,i,j)
3767  end do
3768  end do
3769  select case( i_intrp_land_water )
3770  case( i_intrp_mask )
3771  call make_mask( lmask, work, imax, jmax, landdata=.true.)
3772  !$omp parallel do
3773  do j = 1, jmax
3774  do i = 1, imax
3775  if ( lmask_org(i,j) .ne. undef ) lmask(i,j) = lmask_org(i,j)
3776  end do
3777  end do
3778  case( i_intrp_fill )
3779  call make_mask( lmask, work, imax, jmax, landdata=.true.)
3780  end select
3781  call interp_oceanland_data(work, lmask, imax, jmax, .true., intrp_iter_max)
3782  !$omp parallel do
3783  do j = 1, jmax
3784  do i = 1, imax
3785  lmask(i,j) = maskval_strg
3786  end do
3787  end do
3788  !replace missing value to init_landwater_ratio
3789  call replace_misval_map( work, lmask, imax, jmax, "SMOIS")
3790  !$omp parallel do
3791  do j = 1, jmax
3792  do i = 1, imax
3793  strg_org(k,i,j) = work(i,j)
3794  end do
3795  end do
3796  enddo
3797  end if
3798 
3799  call interp_interp3d( itp_nh_l, &
3800  kmax, 1, kmax, &
3801  imax, jmax, &
3802  lkmax, lks, lke, &
3803  ia, ja, &
3804  igrd( :,:,:), & ! [IN]
3805  jgrd( :,:,:), & ! [IN]
3806  hfact( :,:,:), & ! [IN]
3807  kgrdl(:,:,:,:,:), & ! [IN]
3808  vfactl(:, :,:,:), & ! [IN]
3809  lz3d_org(:,:,:), & ! [IN]
3810  lcz_3d(:,:,:), & ! [IN]
3811  strg_org(:,:,:), & ! [IN]
3812  strg(:,:,:) ) ! [OUT]
3813  end if
3814 
3815  ! replace values over the ocean
3816  do k = 1, lkmax-1
3817  call replace_misval_const( strg(k,:,:), maskval_strg, lsmask_nest )
3818  enddo
3819 
3820  !$omp parallel do collapse(2)
3821  do j = 1, ja
3822  do i = 1, ia
3823  do k = 1, lkmax
3824  strg(k,i,j) = max( min( strg(k,i,j), 1.0_rp ), 0.0_rp )
3825  end do
3826  end do
3827  end do
3828 
3829  if ( filter_niter > 0 ) then
3830  !$omp parallel do collapse(2)
3831  do j = 1, ja
3832  do i = 1, ia
3833  do k = 1, lkmax-1
3834  one3d(k,i,j) = 1.0_rp
3835  end do
3836  end do
3837  end do
3838  call filter_hyperdiff( lkmax, 1, lkmax-1, ia, isb, ieb, ja, jsb, jeb, &
3839  strg(:,:,:), filter_order, filter_niter, &
3840  limiter_sign = one3d(:,:,:) )
3841  call comm_vars8( strg(:,:,:), 1 )
3842  call comm_wait ( strg(:,:,:), 1, .false. )
3843  end if
3844 
3845  !$omp parallel do
3846  do j = 1, ja
3847  do i = 1, ia
3848  strg(lkmax,i,j) = strg(lkmax-1,i,j)
3849  enddo
3850  enddo
3851 
3852  else ! not read from boundary file
3853 
3854  do k = 1, lkmax
3855  !$omp parallel do
3856  do j = 1, ja
3857  do i = 1, ia
3858 ! work(i,j) = init_landwater_ratio_each( LANDUSE_index_PFT(i,j,1) )
3859  work(i,j) = init_landwater_ratio
3860  end do
3861  end do
3862  ! conversion from water saturation [fraction] to volumetric water content [m3/m3]
3863  strg(k,:,:) = convert_ws2vwc( work(:,:), critical=soilwater_ds2vc_flag )
3864  end do
3865 
3866  endif ! use_file_waterratio
3867 
3868 
3869  return

References scale_const::const_eps, scale_const::const_i_lw, scale_const::const_i_sw, scale_const::const_laps, scale_const::const_pi, scale_const::const_undef, mod_land_vars::convert_ws2vwc(), scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_atmos_grid_cartesc_index::ieb, scale_interp::interp_interp2d(), scale_interp::interp_interp3d(), scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::jsb, scale_landuse::landuse_index_pft, scale_landuse::landuse_pft_nmax, scale_land_grid_cartesc_index::lke, scale_land_grid_cartesc_index::lks, make_mask(), scale_prc::prc_abort(), replace_misval_const(), replace_misval_map(), and scale_topography::topography_zsfc.

Referenced by realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ocean_interporation()

subroutine mod_realinput::ocean_interporation ( integer, intent(in)  imax,
integer, intent(in)  jmax,
real(rp), dimension (imax,jmax), intent(in)  sst_org,
real(rp), dimension (imax,jmax), intent(in)  tw_org,
real(rp), dimension(imax,jmax,n_rad_dir,n_rad_rgn), intent(inout)  albw_org,
real(rp), dimension (imax,jmax), intent(inout)  z0w_org,
real(rp), dimension(ia), intent(in)  CX,
real(rp), dimension(ja), intent(in)  CY,
logical, intent(in)  elevation_correction_ocean,
real(rp), intent(in)  init_ocean_alb_lw,
real(rp), intent(in)  init_ocean_alb_sw,
real(rp), intent(in)  init_ocean_z0w,
logical, intent(in)  first_surface,
logical, intent(in)  update_coord,
real(rp), dimension (ia,ja), intent(out)  sst,
real(rp), dimension (ia,ja), intent(out)  tw,
real(rp), dimension(ia,ja,n_rad_dir,n_rad_rgn), intent(out)  albw,
real(rp), dimension (ia,ja), intent(out)  z0w 
)

Definition at line 3881 of file mod_realinput.F90.

3881  use scale_const, only: &
3882  undef => const_undef, &
3883  pi => const_pi, &
3884  laps => const_laps
3885  use scale_topography, only: &
3887  use scale_interp, only: &
3888  interp_factor2d, &
3890  use scale_filter, only: &
3891  filter_hyperdiff
3892  use scale_mapprojection, only: &
3893  mapprojection_lonlat2xy
3894  use scale_comm_cartesc, only: &
3895  comm_vars8, &
3896  comm_wait
3897  implicit none
3898  integer, intent(in) :: imax, jmax
3899  real(RP), intent(in) :: sst_org (imax,jmax)
3900  real(RP), intent(in) :: tw_org (imax,jmax)
3901  real(RP), intent(inout) :: albw_org(imax,jmax,N_RAD_DIR,N_RAD_RGN)
3902  real(RP), intent(inout) :: z0w_org (imax,jmax)
3903  real(RP), intent(in) :: CX(IA)
3904  real(RP), intent(in) :: CY(JA)
3905  logical, intent(in) :: elevation_correction_ocean
3906  real(RP), intent(in) :: init_ocean_alb_lw
3907  real(RP), intent(in) :: init_ocean_alb_sw
3908  real(RP), intent(in) :: init_ocean_z0w
3909  logical, intent(in) :: first_surface
3910  logical, intent(in) :: update_coord
3911 
3912  real(RP), intent(out) :: sst (IA,JA)
3913  real(RP), intent(out) :: tw (IA,JA)
3914  real(RP), intent(out) :: albw(IA,JA,N_RAD_DIR,N_RAD_RGN)
3915  real(RP), intent(out) :: z0w (IA,JA)
3916 
3917  ! for interpolation
3918  real(RP) :: oX_org(imax,jmax)
3919  real(RP) :: oY_org(imax,jmax)
3920  logical :: zonal, pole
3921 
3922  real(RP) :: one(IA,JA)
3923  real(RP) :: tdiff
3924 
3925  integer :: i, j, m, n
3926 
3927  !$omp parallel do
3928  do j = 1, jmax
3929  do i = 1, imax
3930  do m = 1, n_rad_dir
3931  if ( albw_org(i,j,m,i_r_ir ) == undef ) albw_org(i,j,m,i_r_ir ) = init_ocean_alb_lw
3932  if ( albw_org(i,j,m,i_r_nir) == undef ) albw_org(i,j,m,i_r_nir) = init_ocean_alb_sw
3933  if ( albw_org(i,j,m,i_r_vis) == undef ) albw_org(i,j,m,i_r_vis) = init_ocean_alb_sw
3934  if ( albw_org(i,j,m,i_r_vis) == undef ) albw_org(i,j,m,i_r_vis) = init_ocean_alb_sw
3935  end do
3936  if ( z0w_org(i,j) == undef ) z0w_org(i,j) = init_ocean_z0w
3937  end do
3938  end do
3939 
3940  if ( first_surface .or. update_coord ) then
3941 
3942  ! interporation for ocean variables
3943 
3944  select case( itp_type_a )
3945  case ( i_intrp_linear )
3946 
3947  if ( imax == 1 .or. jmax == 1 ) then
3948  log_error("ocean_interporation",*) 'LINER interpolation requires nx, ny > 1'
3949  log_error_cont(*) 'Use "DIST-WEIGHT" as INTRP_TYPE of PARAM_MKINIT_REAL_OCEAN'
3950  call prc_abort
3951  end if
3952 
3953  !$omp parallel do
3954  do j = 1, jmax
3955  do i = 1, imax
3956  olat_org(i,j) = sign( min( abs(olat_org(i,j)), pi * 0.499999_rp ), olat_org(i,j) )
3957  end do
3958  end do
3959 
3960  call mapprojection_lonlat2xy( imax, 1, imax, jmax, 1, jmax, &
3961  olon_org(:,:), olat_org(:,:), & ! [IN]
3962  ox_org(:,:), oy_org(:,:) ) ! [OUT]
3963 
3964  zonal = ( maxval(olon_org) - minval(olon_org) ) > 2.0_rp * pi * 0.9_rp
3965  pole = ( maxval(olat_org) > pi * 0.5_rp * 0.9_rp ) .or. ( minval(olat_org) < - pi * 0.5_rp * 0.9_rp )
3966  call interp_factor2d( imax, jmax, & ! [IN]
3967  ia, ja, & ! [IN]
3968  ox_org(:,:), & ! [IN]
3969  oy_org(:,:), & ! [IN]
3970  cx(:), cy(:), & ! [IN]
3971  oigrd(:,:,:), & ! [OUT]
3972  ojgrd(:,:,:), & ! [OUT]
3973  ohfact(:,:,:), & ! [OUT]
3974  zonal = zonal, & ! [IN]
3975  pole = pole ) ! [IN]
3976 
3977  case ( i_intrp_dstwgt )
3978 
3979  call interp_factor2d( itp_nh_o, & ! [IN]
3980  imax, jmax, & ! [IN]
3981  ia, ja, & ! [IN]
3982  olon_org(:,:), & ! [IN]
3983  olat_org(:,:), & ! [IN]
3984  lon(:,:), lat(:,:), & ! [IN]
3985  oigrd(:,:,:), & ! [OUT]
3986  ojgrd(:,:,:), & ! [OUT]
3987  ohfact(:,:,:) ) ! [OUT]
3988 
3989  end select
3990 
3991  end if
3992 
3993  call interp_interp2d( itp_nh_o, & ! [IN]
3994  imax, jmax, & ! [IN]
3995  ia, ja, & ! [IN]
3996  oigrd(:,:,:), & ! [IN]
3997  ojgrd(:,:,:), & ! [IN]
3998  ohfact(:,:,:), & ! [IN]
3999  tw_org(:,:), & ! [IN]
4000  tw(:,:) ) ! [OUT]
4001  if ( filter_niter > 0 ) then
4002  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
4003  tw(:,:), filter_order, filter_niter )
4004  call comm_vars8( tw(:,:), 1 )
4005  call comm_wait ( tw(:,:), 1, .false. )
4006  end if
4007 
4008  call interp_interp2d( itp_nh_o, & ! [IN]
4009  imax, jmax, & ! [IN]
4010  ia, ja, & ! [IN]
4011  oigrd(:,:,:), & ! [IN]
4012  ojgrd(:,:,:), & ! [IN]
4013  ohfact(:,:,:), & ! [IN]
4014  sst_org(:,:), & ! [IN]
4015  sst(:,:) ) ! [OUT]
4016  if ( filter_niter > 0 ) then
4017  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
4018  sst(:,:), filter_order, filter_niter )
4019  call comm_vars8( sst(:,:), 1 )
4020  call comm_wait ( sst(:,:), 1, .false. )
4021  end if
4022 
4023  ! elevation correction
4024  if ( elevation_correction_ocean ) then
4025 
4026  !$omp parallel do &
4027  !$omp private(tdiff)
4028  do j = 1, ja
4029  do i = 1, ia
4030  tdiff = topography_zsfc(i,j) * laps
4031  sst(i,j) = sst(i,j) - tdiff
4032  tw(i,j) = tw(i,j) - tdiff
4033  end do
4034  end do
4035 
4036  end if
4037 
4038 
4039  if ( filter_niter > 0 ) then
4040  !$omp parallel do
4041  do j = 1, ja
4042  do i = 1, ia
4043  one(i,j) = 1.0_rp
4044  end do
4045  end do
4046  end if
4047 
4048  do n = 1, n_rad_rgn
4049  do m = 1, n_rad_dir
4050 
4051  call interp_interp2d( itp_nh_o, & ! [IN]
4052  imax, jmax, & ! [IN]
4053  ia, ja, & ! [IN]
4054  oigrd(:,:,:), & ! [IN]
4055  ojgrd(:,:,:), & ! [IN]
4056  ohfact(:,:,:), & ! [IN]
4057  albw_org(:,:,m,n), & ! [IN]
4058  albw(:,:,m,n) ) ! [OUT]
4059  if ( filter_niter > 0 ) then
4060  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
4061  albw(:,:,m,n), filter_order, filter_niter, &
4062  limiter_sign = one(:,:) )
4063  call comm_vars8( albw(:,:,m,n), 1 )
4064  call comm_wait ( albw(:,:,m,n), 1, .false. )
4065  end if
4066 
4067  end do
4068  end do
4069 
4070  call interp_interp2d( itp_nh_o, & ! [IN]
4071  imax, jmax, & ! [IN]
4072  ia, ja, & ! [IN]
4073  oigrd(:,:,:), & ! [IN]
4074  ojgrd(:,:,:), & ! [IN]
4075  ohfact(:,:,:), & ! [IN]
4076  z0w_org(:,:), & ! [IN]
4077  z0w(:,:) ) ! [OUT]
4078  if ( filter_niter > 0 ) then
4079  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
4080  z0w(:,:), filter_order, filter_niter, &
4081  limiter_sign = one(:,:) )
4082  call comm_vars8( z0w(:,:), 1 )
4083  call comm_wait ( z0w(:,:), 1, .false. )
4084  end if
4085 
4086 
4087  return

References scale_const::const_laps, scale_const::const_pi, scale_const::const_undef, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_atmos_grid_cartesc_index::ieb, scale_interp::interp_interp2d(), scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::jsb, scale_prc::prc_abort(), and scale_topography::topography_zsfc.

Referenced by realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_input()

subroutine mod_realinput::urban_input ( real(rp), dimension (ia,ja), intent(in)  lst,
real(rp), dimension (ia,ja,n_rad_dir,n_rad_rgn), intent(in)  albg,
real(rp), dimension(ia,ja), intent(out)  tc_urb,
real(rp), dimension(ia,ja), intent(out)  qc_urb,
real(rp), dimension(ia,ja), intent(out)  uc_urb,
real(rp), dimension (ia,ja), intent(out)  ust,
real(rp), dimension (ia,ja,n_rad_dir,n_rad_rgn), intent(out)  albu 
)

Definition at line 4095 of file mod_realinput.F90.

4095  use mod_atmos_vars, only: &
4096  dens, &
4097  momx, &
4098  momy, &
4099  rhot, &
4100  qtrc
4101  use scale_atmos_hydrometeor, only: &
4102  i_qv
4103  use scale_atmos_thermodyn, only: &
4104  thermodyn_specific_heat => atmos_thermodyn_specific_heat, &
4105  thermodyn_rhot2temp_pres => atmos_thermodyn_rhot2temp_pres
4106  use scale_comm_cartesc, only: &
4107  comm_vars8, &
4108  comm_wait
4109  implicit none
4110  real(RP), intent(in) :: lst (IA,JA)
4111  real(RP), intent(in) :: albg (IA,JA,N_RAD_DIR,N_RAD_RGN)
4112  real(RP), intent(out) :: tc_urb(IA,JA)
4113  real(RP), intent(out) :: qc_urb(IA,JA)
4114  real(RP), intent(out) :: uc_urb(IA,JA)
4115  real(RP), intent(out) :: ust (IA,JA)
4116  real(RP), intent(out) :: albu (IA,JA,N_RAD_DIR,N_RAD_RGN)
4117 
4118  real(RP) :: temp, pres
4119  real(RP) :: Qdry
4120  real(RP) :: Rtot
4121  real(RP) :: CVtot
4122  real(RP) :: CPtot
4123 
4124  integer :: i, j
4125 
4126  ! urban data
4127 
4128  !$omp parallel do collapse(2) &
4129  !$omp private(Qdry,Rtot,CVtot,CPtot,temp,pres)
4130  do j = 1, ja
4131  do i = 1, ia
4132  call thermodyn_specific_heat( qa, &
4133  qtrc(ks,i,j,:), &
4134  tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), & ! [IN]
4135  qdry, rtot, cvtot, cptot ) ! [OUT]
4136  call thermodyn_rhot2temp_pres( dens(ks,i,j), rhot(ks,i,j), rtot, cvtot, cptot, & ! [IN]
4137  temp, pres ) ! [OUT]
4138 
4139  tc_urb(i,j) = temp
4140  if ( i_qv > 0 ) then
4141  qc_urb(i,j) = qtrc(ks,i,j,i_qv)
4142  else
4143  qc_urb(i,j) = 0.0_rp
4144  end if
4145  enddo
4146  enddo
4147 
4148  !$omp parallel do
4149  do j = 1, ja-1
4150  do i = 1, ia-1
4151  uc_urb(i,j) = max(sqrt( ( momx(ks,i,j) / (dens(ks,i+1, j)+dens(ks,i,j)) * 2.0_rp )**2.0_rp &
4152  + ( momy(ks,i,j) / (dens(ks, i,j+1)+dens(ks,i,j)) * 2.0_rp )**2.0_rp ), &
4153  0.01_rp)
4154  enddo
4155  enddo
4156  !$omp parallel do
4157  do j = 1, ja-1
4158  uc_urb(ia,j) = max(sqrt( ( momx(ks,ia,j) / dens(ks,ia,j ) )**2.0_rp &
4159  + ( momy(ks,ia,j) / (dens(ks,ia,j+1)+dens(ks,ia,j)) * 2.0_rp )**2.0_rp ), &
4160  0.01_rp)
4161  enddo
4162  !$omp parallel do
4163  do i = 1, ia-1
4164  uc_urb(i,ja) = max(sqrt( ( momx(ks,i,ja) / (dens(ks,i+1,ja)+dens(ks,i,ja)) * 2.0_rp )**2.0_rp &
4165  + ( momy(ks,i,ja) / dens(ks,i ,ja) )**2.0_rp ), 0.01_rp)
4166  enddo
4167  uc_urb(ia,ja) = max(sqrt( ( momx(ks,ia,ja) / dens(ks,ia,ja) )**2.0_rp &
4168  + ( momy(ks,ia,ja) / dens(ks,ia,ja) )**2.0_rp ), 0.01_rp)
4169 
4170  call comm_vars8( uc_urb, 1 )
4171  call comm_wait ( uc_urb, 1, .false. )
4172 
4173 
4174 !!$ ! Urban surface temp: interpolate over the ocean
4175 !!$ if ( i_INTRP_URB_SFC_TEMP .ne. i_intrp_off ) then
4176 !!$ select case( i_INTRP_URB_SFC_TEMP )
4177 !!$ case( i_intrp_mask )
4178 !!$ call make_mask( lmask, ust_org, imax, jmax, landdata=.true.)
4179 !!$ !$omp parallel do
4180 !!$ do j = 1, jmax
4181 !!$ do i = 1, imax
4182 !!$ if ( lmask_org(i,j) .ne. UNDEF ) lmask(i,j) = lmask_org(i,j)
4183 !!$ end do
4184 !!$ end do
4185 !!$ case( i_intrp_fill )
4186 !!$ call make_mask( lmask, ust_org, imax, jmax, landdata=.true.)
4187 !!$ case default
4188 !!$ LOG_ERROR("urban_input",*) 'INTRP_URB_SFC_TEMP is invalid.'
4189 !!$ call PRC_abort
4190 !!$ end select
4191 !!$ call interp_OceanLand_data(ust_org, lmask, imax, jmax, .true., intrp_iter_max)
4192 !!$ end if
4193 !!$
4194 !!$ !$omp parallel do
4195 !!$ do j = 1, jmax
4196 !!$ do i = 1, imax
4197 !!$ if ( ust_org(i,j) == UNDEF ) ust_org(i,j) = lst_org(i,j)
4198 !!$ end do
4199 !!$ end do
4200 !!$
4201 !!$ call INTERP_interp2d( itp_nh_l, & ! [IN]
4202 !!$ imax, jmax, & ! [IN]
4203 !!$ IA, JA, & ! [IN]
4204 !!$ igrd (:,:,:), & ! [IN]
4205 !!$ jgrd (:,:,:), & ! [IN]
4206 !!$ hfact (:,:,:), & ! [IN]
4207 !!$ ust_org (:,:), & ! [IN]
4208 !!$ ust (:,:) ) ! [OUT]
4209 !!$ if ( FILTER_NITER > 0 ) then
4210 !!$ call FILTER_hyperdiff( IA, ISB, IEB, JA, JSB, JEB, &
4211 !!$ ust(:,:), FILTER_ORDER, FILTER_NITER )
4212 !!$ call COMM_vars8( ust(:,:), 1 )
4213 !!$ call COMM_wait ( ust(:,:), 1, .false. )
4214 !!$ end if
4215 !!$
4216 !!$ !$omp parallel do
4217 !!$ do j = 1, JA
4218 !!$ do i = 1, IA
4219 !!$ if( abs(lsmask_nest(i,j)-0.0_RP) < EPS ) then ! ocean grid
4220 !!$ ust(i,j) = sst(i,j,nn)
4221 !!$ endif
4222 !!$ enddo
4223 !!$ enddo
4224 !!$
4225 !!$ !$omp parallel do &
4226 !!$ !$omp private(tdiff)
4227 !!$ do j = 1, JA
4228 !!$ do i = 1, IA
4229 !!$ if ( topo(i,j) > 0.0_RP ) then ! ignore UNDEF value
4230 !!$ tdiff = ( TOPOGRAPHY_Zsfc(i,j) - topo(i,j) ) * LAPS
4231 !!$ ust(i,j) = ust(i,j) - tdiff
4232 !!$ end if
4233 !!$ end do
4234 !!$ end do
4235 
4236 
4237  !$omp parallel do
4238  do j = 1, ja
4239  do i = 1, ia
4240  ust(i,j) = lst(i,j)
4241  end do
4242  end do
4243 
4244 
4245  ! copy albedo of land to urban
4246  !$omp parallel do
4247  do j = 1, ja
4248  do i = 1, ia
4249  albu(i,j,:,:) = albg(i,j,:,:)
4250  enddo
4251  enddo
4252 
4253  return

References mod_atmos_vars::dens, scale_atmos_hydrometeor::i_qv, scale_atmos_grid_cartesc_index::ks, mod_atmos_vars::momx, mod_atmos_vars::momy, scale_tracer::qa, mod_atmos_vars::qtrc, mod_atmos_vars::rhot, scale_tracer::tracer_cp, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by realinput_surface().

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 4263 of file mod_realinput.F90.

4263  use scale_const, only: &
4264  eps => const_eps, &
4265  undef => const_undef
4266  implicit none
4267  real(RP), intent(out) :: gmask(:,:)
4268  real(RP), intent(in) :: data(:,:)
4269  integer, intent(in) :: nx
4270  integer, intent(in) :: ny
4271  logical, intent(in) :: landdata ! .true. => land data , .false. => ocean data
4272 
4273  real(RP) :: dd
4274  integer :: i,j
4275 
4276  if( landdata )then
4277  !$omp parallel do
4278  do j = 1, ny
4279  do i = 1, nx
4280  gmask(i,j) = 1.0_rp ! gmask=1 will be skip in "interp_OceanLand_data"
4281  end do
4282  end do
4283  dd = 0.0_rp
4284  else
4285  !$omp parallel do
4286  do j = 1, ny
4287  do i = 1, nx
4288  gmask(i,j) = 0.0_rp ! gmask=0 will be skip in "interp_OceanLand_data"
4289  end do
4290  end do
4291  dd = 1.0_rp
4292  endif
4293 
4294  !$omp parallel do
4295  do j = 1, ny
4296  do i = 1, nx
4297  if( abs(data(i,j) - undef) < sqrt(eps) )then
4298  gmask(i,j) = dd
4299  endif
4300  enddo
4301  enddo
4302 
4303  return

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

Referenced by land_interporation(), and realinput_surface().

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 4438 of file mod_realinput.F90.

4438  use scale_const, only: &
4439  eps => const_eps
4440  implicit none
4441  real(RP), intent(inout) :: data(:,:)
4442  real(RP), intent(in) :: maskval
4443  real(RP), intent(in) :: frac_land(:,:)
4444  integer :: i, j
4445 
4446  !$omp parallel do
4447  do j = 1, ja
4448  do i = 1, ia
4449  if( abs(frac_land(i,j)-0.0_rp) < eps )then ! ocean grid
4450  data(i,j) = maskval
4451  endif
4452  enddo
4453  enddo
4454 

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

Referenced by land_interporation().

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 4459 of file mod_realinput.F90.

4459  use scale_const, only: &
4460  eps => const_eps, &
4461  undef => const_undef
4462  implicit none
4463 
4464  real(RP), intent(inout) :: data(:,:)
4465  real(RP), intent(in) :: maskval(:,:)
4466  integer, intent(in) :: nx, ny
4467  character(len=*), intent(in) :: elem
4468 
4469  integer :: i, j
4470  logical :: error
4471 
4472  error = .false.
4473  !$omp parallel do
4474  do j = 1, ny
4475  if ( error ) cycle
4476  do i = 1, nx
4477  if( abs(data(i,j) - undef) < sqrt(eps) )then
4478  if( abs(maskval(i,j) - undef) < sqrt(eps) )then
4479  log_error("replace_misval_map",*) "data for mask of "//trim(elem)//"(",i,",",j,") includes missing value."
4480  error = .true.
4481  exit
4482  else
4483  data(i,j) = maskval(i,j)
4484  endif
4485  endif
4486  enddo
4487  enddo
4488 
4489  if ( error ) then
4490  log_error_cont(*) "Please check input data of SKINTEMP or SST. "
4491  call prc_abort
4492  end if
4493 
4494  return

References scale_const::const_eps, scale_const::const_undef, and scale_prc::prc_abort().

Referenced by land_interporation(), and realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ iscale

integer, parameter, public mod_realinput::iscale = 1

Definition at line 72 of file mod_realinput.F90.

72  integer, public, parameter :: iSCALE = 1

Referenced by realinput_surface().

◆ iwrfarw

integer, parameter, public mod_realinput::iwrfarw = 2

Definition at line 73 of file mod_realinput.F90.

73  integer, public, parameter :: iWRFARW = 2

Referenced by realinput_surface().

◆ igrads

integer, parameter, public mod_realinput::igrads = 4

Definition at line 75 of file mod_realinput.F90.

75  integer, public, parameter :: iGrADS = 4

Referenced by realinput_surface().

mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
Definition: mod_atmos_phy_sf_vars.F90:66
scale_interp::interp_interp2d
subroutine, public interp_interp2d(npoints, IA_ref, JA_ref, IA, JA, idx_i, idx_j, hfact, val_ref, val, threshold_undef, wsum, val2)
Definition: scale_interp.F90:1232
mod_ocean_vars::ocean_ocn_z0m
real(rp), dimension(:,:), allocatable, public ocean_ocn_z0m
surface roughness length for momentum, open ocean [m]
Definition: mod_ocean_vars.F90:65
scale_landuse::landuse_fact_ocean
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
Definition: scale_landuse.F90:44
mod_land_vars::land_temp
real(rp), dimension(:,:,:), allocatable, public land_temp
temperature of each soil layer [K]
Definition: mod_land_vars.F90:61
mod_urban_admin::urban_do
logical, public urban_do
Definition: mod_urban_admin.F90:32
mod_urban_vars::urban_sfc_temp
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
Definition: mod_urban_vars.F90:74
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
mod_ocean_vars::ocean_ice_temp
real(rp), dimension(:,:), allocatable, public ocean_ice_temp
sea ice temperature [K]
Definition: mod_ocean_vars.F90:73
mod_urban_vars::urban_tb
real(rp), dimension(:,:), allocatable, public urban_tb
Definition: mod_urban_vars.F90:64
mod_ocean_vars::ocean_sfc_z0e
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0e
ocean surface roughness length for vapor [m]
Definition: mod_ocean_vars.F90:71
scale_ocean_phy_ice_simple
module ocean / physics / ice / simple
Definition: scale_ocean_phy_ice_simple.F90:12
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0e
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
Definition: mod_atmos_phy_sf_vars.F90:68
mod_land_vars::land_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public land_sfc_albedo
land surface albedo (direct/diffuse,IR/near-IR/VIS) (0-1)
Definition: mod_land_vars.F90:65
mod_ocean_vars::ice_flag
logical, public ice_flag
Definition: mod_ocean_vars.F90:141
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
mod_urban_vars::urban_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public urban_sfc_albedo
Definition: mod_urban_vars.F90:75
scale_atmos_phy_rd_common::i_lw
integer, parameter, public i_lw
Definition: scale_atmos_phy_rd_common.F90:37
scale_interp
module INTERPOLATION
Definition: scale_interp.F90:12
scale_land_grid_cartesc_index::lkmax
integer, public lkmax
Definition: scale_land_grid_cartesC_index.F90:32
mod_land_vars::land_water
real(rp), dimension(:,:,:), allocatable, public land_water
moisture of each soil layer [m3/m3]
Definition: mod_land_vars.F90:62
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
mod_urban_vars::urban_trl
real(rp), dimension(:,:,:), allocatable, public urban_trl
Definition: mod_urban_vars.F90:60
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
mod_atmos_phy_sf_vars
module ATMOSPHERIC Surface Variables
Definition: mod_atmos_phy_sf_vars.F90:12
scale_ocean_grid_cartesc_index::okmax
integer, public okmax
Definition: scale_ocean_grid_cartesC_index.F90:32
mod_urban_vars::urban_raing
real(rp), dimension(:,:), allocatable, public urban_raing
Definition: mod_urban_vars.F90:71
scale_topography
module TOPOGRAPHY
Definition: scale_topography.F90:11
mod_urban_vars::urban_qc
real(rp), dimension(:,:), allocatable, public urban_qc
Definition: mod_urban_vars.F90:67
mod_ocean_vars::ocean_sfc_z0m
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0m
ocean surface roughness length for momentum [m]
Definition: mod_ocean_vars.F90:69
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:79
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
Definition: mod_atmos_phy_sf_vars.F90:67
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_temp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
Definition: mod_atmos_phy_sf_vars.F90:64
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:80
scale_const::const_i_sw
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:97
scale_urban_grid_cartesc_index::ukmax
integer, public ukmax
Definition: scale_urban_grid_cartesC_index.F90:32
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const::const_pi
real(rp), public const_pi
pi
Definition: scale_const.F90:31
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:75
mod_ocean_vars::ocean_sfc_z0h
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0h
ocean surface roughness length for heat [m]
Definition: mod_ocean_vars.F90:70
mod_ocean_admin
module Ocean admin
Definition: mod_ocean_admin.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
mod_ocean_vars::ocean_uvel
real(rp), dimension(:,:,:), allocatable, public ocean_uvel
ocean zonal velocity [m/s]
Definition: mod_ocean_vars.F90:62
mod_land_admin
module Land admin
Definition: mod_land_admin.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
Definition: mod_atmos_phy_sf_vars.F90:65
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:76
mod_urban_vars::urban_tgl
real(rp), dimension(:,:,:), allocatable, public urban_tgl
Definition: mod_urban_vars.F90:62
scale_landuse::landuse_index_pft
integer, dimension(:,:,:), allocatable, public landuse_index_pft
index of PFT for each mosaic
Definition: scale_landuse.F90:67
mod_land_vars
module LAND Variables
Definition: mod_land_vars.F90:11
scale_landuse::landuse_pft_nmax
integer, public landuse_pft_nmax
number of plant functional type(PFT)
Definition: scale_landuse.F90:63
scale_ocean_phy_ice_simple::ocean_phy_ice_freezetemp
real(rp), public ocean_phy_ice_freezetemp
Definition: scale_ocean_phy_ice_simple.F90:44
mod_atmos_vars::momx
real(rp), dimension(:,:,:), allocatable, target, public momx
Definition: mod_atmos_vars.F90:77
mod_ocean_vars
module OCEAN Variables
Definition: mod_ocean_vars.F90:12
mod_atmos_vars::momy
real(rp), dimension(:,:,:), allocatable, target, public momy
Definition: mod_atmos_vars.F90:78
scale_time
module TIME
Definition: scale_time.F90:11
scale_interp::interp_interp3d
subroutine, public interp_interp3d(npoints, KA_ref, KS_ref, KE_ref, IA_ref, JA_ref, KA, KS, KE, IA, JA, idx_i, idx_j, hfact, idx_k, vfact, hgt_ref, hgt, val_ref, val, spline, logwgt, threshold_undef, wsum, val2)
Definition: scale_interp.F90:1322
scale_const::const_i_lw
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:96
mod_urban_vars::urban_tc
real(rp), dimension(:,:), allocatable, public urban_tc
Definition: mod_urban_vars.F90:66
mod_urban_vars::urban_tg
real(rp), dimension(:,:), allocatable, public urban_tg
Definition: mod_urban_vars.F90:65
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:77
mod_land_vars::convert_ws2vwc
real(rp) function, dimension(lia, lja), public convert_ws2vwc(WS, critical)
conversion from water saturation [fraction] to volumetric water content [m3/m3]
Definition: mod_land_vars.F90:1253
scale_ocean_grid_cartesc_index::oimax
integer, public oimax
Definition: scale_ocean_grid_cartesC_index.F90:33
mod_urban_vars::urban_uc
real(rp), dimension(:,:), allocatable, public urban_uc
Definition: mod_urban_vars.F90:68
scale_atmos_phy_rd_common::i_sw
integer, parameter, public i_sw
Definition: scale_atmos_phy_rd_common.F90:38
scale_landuse::landuse_fact_land
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
Definition: scale_landuse.F90:45
scale_topography::topography_zsfc
real(rp), dimension(:,:), allocatable, public topography_zsfc
absolute ground height [m]
Definition: scale_topography.F90:38
mod_land_vars::land_sfc_temp
real(rp), dimension(:,:), allocatable, public land_sfc_temp
land surface skin temperature [K]
Definition: mod_land_vars.F90:64
mod_ocean_admin::ocean_do
logical, public ocean_do
Definition: mod_ocean_admin.F90:32
mod_ocean_vars::ocean_ice_mass
real(rp), dimension(:,:), allocatable, public ocean_ice_mass
sea ice mass [kg]
Definition: mod_ocean_vars.F90:74
scale_const::const_tem00
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
mod_ocean_vars::ocean_vvel
real(rp), dimension(:,:,:), allocatable, public ocean_vvel
ocean meridional velocity [m/s]
Definition: mod_ocean_vars.F90:63
scale_time::time_gettimelabel
subroutine, public time_gettimelabel(timelabel)
generate time label
Definition: scale_time.F90:91
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
mod_land_admin::land_do
logical, public land_do
Definition: mod_land_admin.F90:41
mod_urban_vars::urban_rainb
real(rp), dimension(:,:), allocatable, public urban_rainb
Definition: mod_urban_vars.F90:70
mod_ocean_vars::ocean_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public ocean_sfc_albedo
ocean surface albedo (direct/diffuse,IR/near-IR/VIS) (0-1)
Definition: mod_ocean_vars.F90:68
mod_urban_vars::urban_tr
real(rp), dimension(:,:), allocatable, public urban_tr
Definition: mod_urban_vars.F90:63
mod_ocean_vars::ocean_sfc_temp
real(rp), dimension(:,:), allocatable, public ocean_sfc_temp
ocean surface skin temperature [K]
Definition: mod_ocean_vars.F90:67
scale_filter
module FILTER
Definition: scale_filter.F90:11
mod_atmos_admin::atmos_phy_mp_type
character(len=h_short), public atmos_phy_mp_type
Definition: mod_atmos_admin.F90:36
scale_atmos_thermodyn
module atmosphere / thermodyn
Definition: scale_atmos_thermodyn.F90:11
mod_urban_vars::urban_roff
real(rp), dimension(:,:), allocatable, public urban_roff
Definition: mod_urban_vars.F90:91
scale_const::const_laps
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:58
scale_landuse::landuse_fact_urban
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
Definition: scale_landuse.F90:46
scale_land_grid_cartesc_index::lke
integer, public lke
Definition: scale_land_grid_cartesC_index.F90:41
mod_ocean_vars::ocean_salt
real(rp), dimension(:,:,:), allocatable, public ocean_salt
ocean salinity [PSU]
Definition: mod_ocean_vars.F90:61
mod_urban_vars::urban_tbl
real(rp), dimension(:,:,:), allocatable, public urban_tbl
Definition: mod_urban_vars.F90:61
scale_landuse
module LANDUSE
Definition: scale_landuse.F90:19
mod_urban_admin
module Urban admin
Definition: mod_urban_admin.F90:11
mod_land_vars::land_ice
real(rp), dimension(:,:,:), allocatable, public land_ice
ice of each soil layer [m3/m3]
Definition: mod_land_vars.F90:63
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_ocean_grid_cartesc_index::ojmax
integer, public ojmax
Definition: scale_ocean_grid_cartesC_index.F90:34
mod_urban_vars::urban_rainr
real(rp), dimension(:,:), allocatable, public urban_rainr
Definition: mod_urban_vars.F90:69
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
scale_land_grid_cartesc_index::lks
integer, public lks
Definition: scale_land_grid_cartesC_index.F90:40
mod_urban_vars
module URBAN Variables
Definition: mod_urban_vars.F90:12
mod_ocean_vars::ocean_temp
real(rp), dimension(:,:,:), allocatable, public ocean_temp
ocean temperature [K]
Definition: mod_ocean_vars.F90:60
scale_ocean_grid_cartesc_index::oks
integer, public oks
Definition: scale_ocean_grid_cartesC_index.F90:37